Actual source code: dmimpl.h
petsc-3.7.6 2017-04-24
3: #if !defined(_DMIMPL_H)
4: #define _DMIMPL_H
6: #include <petscdm.h>
7: #include <petsc/private/petscimpl.h>
9: PETSC_EXTERN PetscBool DMRegisterAllCalled;
10: PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
11: typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt field, MatNullSpace *nullSpace);
13: typedef struct _DMOps *DMOps;
14: struct _DMOps {
15: PetscErrorCode (*view)(DM,PetscViewer);
16: PetscErrorCode (*load)(DM,PetscViewer);
17: PetscErrorCode (*clone)(DM,DM*);
18: PetscErrorCode (*setfromoptions)(PetscOptionItems*,DM);
19: PetscErrorCode (*setup)(DM);
20: PetscErrorCode (*createdefaultsection)(DM);
21: PetscErrorCode (*createdefaultconstraints)(DM);
22: PetscErrorCode (*createglobalvector)(DM,Vec*);
23: PetscErrorCode (*createlocalvector)(DM,Vec*);
24: PetscErrorCode (*getlocaltoglobalmapping)(DM);
25: PetscErrorCode (*createfieldis)(DM,PetscInt*,char***,IS**);
26: PetscErrorCode (*createcoordinatedm)(DM,DM*);
28: PetscErrorCode (*getcoloring)(DM,ISColoringType,ISColoring*);
29: PetscErrorCode (*creatematrix)(DM, Mat*);
30: PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*);
31: PetscErrorCode (*createrestriction)(DM,DM,Mat*);
32: PetscErrorCode (*getaggregates)(DM,DM,Mat*);
33: PetscErrorCode (*getinjection)(DM,DM,Mat*);
35: PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
36: PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*);
37: PetscErrorCode (*refinehierarchy)(DM,PetscInt,DM*);
38: PetscErrorCode (*coarsenhierarchy)(DM,PetscInt,DM*);
40: PetscErrorCode (*globaltolocalbegin)(DM,Vec,InsertMode,Vec);
41: PetscErrorCode (*globaltolocalend)(DM,Vec,InsertMode,Vec);
42: PetscErrorCode (*localtoglobalbegin)(DM,Vec,InsertMode,Vec);
43: PetscErrorCode (*localtoglobalend)(DM,Vec,InsertMode,Vec);
44: PetscErrorCode (*localtolocalbegin)(DM,Vec,InsertMode,Vec);
45: PetscErrorCode (*localtolocalend)(DM,Vec,InsertMode,Vec);
47: PetscErrorCode (*destroy)(DM);
49: PetscErrorCode (*computevariablebounds)(DM,Vec,Vec);
51: PetscErrorCode (*createsubdm)(DM,PetscInt,PetscInt*,IS*,DM*);
52: PetscErrorCode (*createfielddecomposition)(DM,PetscInt*,char***,IS**,DM**);
53: PetscErrorCode (*createdomaindecomposition)(DM,PetscInt*,char***,IS**,IS**,DM**);
54: PetscErrorCode (*createddscatters)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
56: PetscErrorCode (*getdimpoints)(DM,PetscInt,PetscInt*,PetscInt*);
57: PetscErrorCode (*locatepoints)(DM,Vec,PetscSF);
59: PetscErrorCode (*projectfunctionlocal)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
60: PetscErrorCode (*projectfunctionlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
61: PetscErrorCode (*projectfieldlocal)(DM,Vec,void(**)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscScalar[]),InsertMode,Vec);
62: PetscErrorCode (*computel2diff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
63: PetscErrorCode (*computel2gradientdiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal [],const PetscReal[],PetscInt, PetscScalar *,void *),void **,Vec,const PetscReal[],PetscReal *);
64: PetscErrorCode (*computel2fielddiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
65: };
67: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
68: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
69: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
71: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
72: struct _DMCoarsenHookLink {
73: PetscErrorCode (*coarsenhook)(DM,DM,void*); /* Run once, when coarse DM is created */
74: PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*); /* Run each time a new problem is restricted to a coarse grid */
75: void *ctx;
76: DMCoarsenHookLink next;
77: };
79: typedef struct _DMRefineHookLink *DMRefineHookLink;
80: struct _DMRefineHookLink {
81: PetscErrorCode (*refinehook)(DM,DM,void*); /* Run once, when a fine DM is created */
82: PetscErrorCode (*interphook)(DM,Mat,DM,void*); /* Run each time a new problem is interpolated to a fine grid */
83: void *ctx;
84: DMRefineHookLink next;
85: };
87: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
88: struct _DMSubDomainHookLink {
89: PetscErrorCode (*ddhook)(DM,DM,void*);
90: PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*);
91: void *ctx;
92: DMSubDomainHookLink next;
93: };
95: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
96: struct _DMGlobalToLocalHookLink {
97: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
98: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
99: void *ctx;
100: DMGlobalToLocalHookLink next;
101: };
103: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
104: struct _DMLocalToGlobalHookLink {
105: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
106: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
107: void *ctx;
108: DMLocalToGlobalHookLink next;
109: };
111: typedef enum {DMVEC_STATUS_IN,DMVEC_STATUS_OUT} DMVecStatus;
112: typedef struct _DMNamedVecLink *DMNamedVecLink;
113: struct _DMNamedVecLink {
114: Vec X;
115: char *name;
116: DMVecStatus status;
117: DMNamedVecLink next;
118: };
120: typedef struct _DMWorkLink *DMWorkLink;
121: struct _DMWorkLink {
122: size_t bytes;
123: void *mem;
124: DMWorkLink next;
125: };
127: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
129: struct _n_DMLabelLink {
130: DMLabel label;
131: PetscBool output;
132: struct _n_DMLabelLink *next;
133: };
134: typedef struct _n_DMLabelLink *DMLabelLink;
136: struct _n_DMLabelLinkList {
137: PetscInt refct;
138: DMLabelLink next;
139: };
140: typedef struct _n_DMLabelLinkList *DMLabelLinkList;
142: typedef struct _n_Boundary *DMBoundary;
144: struct _n_Boundary {
145: const char *name;
146: const char *labelname;
147: DMLabel label;
148: PetscBool essential;
149: PetscInt field;
150: PetscInt numcomps;
151: PetscInt *comps;
152: void (*func)();
153: PetscInt numids;
154: PetscInt *ids;
155: void *ctx;
156: DMBoundary next;
157: };
159: struct _n_DMBoundaryLinkList {
160: PetscInt refct;
161: DMBoundary next;
162: };
164: typedef struct _n_DMBoundaryLinkList *DMBoundaryLinkList;
166: PETSC_EXTERN PetscErrorCode DMDestroyLabelLinkList(DM);
167: PETSC_EXTERN PetscErrorCode DMBoundaryDestroy(DMBoundaryLinkList*);
168: PETSC_EXTERN PetscErrorCode DMBoundaryDuplicate(DMBoundaryLinkList,DMBoundaryLinkList*);
170: struct _p_DM {
171: PETSCHEADER(struct _DMOps);
172: Vec localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
173: Vec globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
174: DMNamedVecLink namedglobal;
175: DMNamedVecLink namedlocal;
176: DMWorkLink workin,workout;
177: DMLabelLinkList labels; /* Linked list of labels */
178: DMLabel depthLabel; /* Optimized access to depth label */
179: void *ctx; /* a user context */
180: PetscErrorCode (*ctxdestroy)(void**);
181: Vec x; /* location at which the functions/Jacobian are computed */
182: ISColoringType coloringtype;
183: MatFDColoring fd;
184: VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
185: MatType mattype; /* type of matrix created with DMCreateMatrix() */
186: PetscInt bs;
187: ISLocalToGlobalMapping ltogmap;
188: PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
189: PetscInt levelup,leveldown; /* if the DM has been obtained by refining (or coarsening) this indicates how many times that process has been used to generate this DM */
190: PetscBool setupcalled; /* Indicates that the DM has been set up, methods that modify a DM such that a fresh setup is required should reset this flag */
191: void *data;
192: /* Hierarchy / Submeshes */
193: DM coarseMesh;
194: DM fineMesh;
195: DMCoarsenHookLink coarsenhook; /* For transfering auxiliary problem data to coarser grids */
196: DMRefineHookLink refinehook;
197: DMSubDomainHookLink subdomainhook;
198: DMGlobalToLocalHookLink gtolhook;
199: DMLocalToGlobalHookLink ltoghook;
200: /* Topology */
201: PetscInt dim; /* The topological dimension */
202: /* Flexible communication */
203: PetscSF sf; /* SF for parallel point overlap */
204: PetscSF defaultSF; /* SF for parallel dof overlap using default section */
205: PetscSF sfNatural; /* SF mapping to the "natural" ordering */
206: PetscBool useNatural; /* Create the natural SF */
207: /* Allows a non-standard data layout */
208: PetscSection defaultSection; /* Layout for local vectors */
209: PetscSection defaultGlobalSection; /* Layout for global vectors */
210: PetscLayout map;
211: /* Constraints */
212: PetscSection defaultConstraintSection;
213: Mat defaultConstraintMat;
214: /* Coordinates */
215: PetscInt dimEmbed; /* The dimension of the embedding space */
216: DM coordinateDM; /* Layout for coordinates (default section) */
217: Vec coordinates; /* Coordinate values in global vector */
218: Vec coordinatesLocal; /* Coordinate values in local vector */
219: PetscReal *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
220: DMBoundaryType *bdtype; /* Indicates type of topological boundary */
221: /* Null spaces -- of course I should make this have a variable number of fields */
222: /* I now believe this might not be the right way: see below */
223: NullSpaceFunc nullspaceConstructors[10];
224: /* Fields are represented by objects */
225: PetscDS prob;
226: DMBoundaryLinkList boundary; /* List of boundary conditions */
227: /* Output structures */
228: DM dmBC; /* The DM with boundary conditions in the global DM */
229: PetscInt outputSequenceNum; /* The current sequence number for output */
230: PetscReal outputSequenceVal; /* The current sequence value for output */
232: PetscObject dmksp,dmsnes,dmts;
233: };
235: PETSC_EXTERN PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocatePoints, DM_Coarsen, DM_CreateInterpolation, DM_CreateRestriction;
237: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
238: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);
239: PETSC_EXTERN PetscErrorCode DMCreateSubDM_Section_Private(DM,PetscInt,PetscInt[],IS*,DM*);
241: /*
243: Composite Vectors
245: Single global representation
246: Individual global representations
247: Single local representation
248: Individual local representations
250: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
252: DM da_u, da_v, da_p
254: DM dm_velocities
256: DM dm
258: DMDACreate(,&da_u);
259: DMDACreate(,&da_v);
260: DMCompositeCreate(,&dm_velocities);
261: DMCompositeAddDM(dm_velocities,(DM)du);
262: DMCompositeAddDM(dm_velocities,(DM)dv);
264: DMDACreate(,&da_p);
265: DMCompositeCreate(,&dm_velocities);
266: DMCompositeAddDM(dm,(DM)dm_velocities);
267: DMCompositeAddDM(dm,(DM)dm_p);
270: Access parts of composite vectors (Composite only)
271: ---------------------------------
272: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
273: ADD for local vector -
275: Element access
276: --------------
277: From global vectors
278: -DAVecGetArray - for DMDA
279: -VecGetArray - for DMSliced
280: ADD for DMComposite??? maybe
282: From individual vector
283: -DAVecGetArray - for DMDA
284: -VecGetArray -for sliced
285: ADD for DMComposite??? maybe
287: From single local vector
288: ADD * single local vector as arrays?
290: Communication
291: -------------
292: DMGlobalToLocal - global vector to single local vector
294: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
296: Obtaining vectors
297: -----------------
298: DMCreateGlobal/Local
299: DMGetGlobal/Local
300: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
303: ????? individual global vectors ????
305: */
307: #if defined(PETSC_HAVE_HDF5)
308: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
309: #endif
311: #endif