Actual source code: dmimpl.h
petsc-3.9.0 2018-04-07
3: #if !defined(_DMIMPL_H)
4: #define _DMIMPL_H
6: #include <petscdm.h>
7: #include <petsc/private/petscimpl.h>
8: #include <petsc/private/petscdsimpl.h>
10: PETSC_EXTERN PetscBool DMRegisterAllCalled;
11: PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
12: typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt field, MatNullSpace *nullSpace);
14: typedef struct _DMOps *DMOps;
15: struct _DMOps {
16: PetscErrorCode (*view)(DM,PetscViewer);
17: PetscErrorCode (*load)(DM,PetscViewer);
18: PetscErrorCode (*clone)(DM,DM*);
19: PetscErrorCode (*setfromoptions)(PetscOptionItems*,DM);
20: PetscErrorCode (*setup)(DM);
21: PetscErrorCode (*createdefaultsection)(DM);
22: PetscErrorCode (*createdefaultconstraints)(DM);
23: PetscErrorCode (*createglobalvector)(DM,Vec*);
24: PetscErrorCode (*createlocalvector)(DM,Vec*);
25: PetscErrorCode (*getlocaltoglobalmapping)(DM);
26: PetscErrorCode (*createfieldis)(DM,PetscInt*,char***,IS**);
27: PetscErrorCode (*createcoordinatedm)(DM,DM*);
29: PetscErrorCode (*getcoloring)(DM,ISColoringType,ISColoring*);
30: PetscErrorCode (*creatematrix)(DM, Mat*);
31: PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*);
32: PetscErrorCode (*createrestriction)(DM,DM,Mat*);
33: PetscErrorCode (*createmassmatrix)(DM,DM,Mat*);
34: PetscErrorCode (*getaggregates)(DM,DM,Mat*);
35: PetscErrorCode (*hascreateinjection)(DM,PetscBool*);
36: PetscErrorCode (*getinjection)(DM,DM,Mat*);
38: PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
39: PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*);
40: PetscErrorCode (*refinehierarchy)(DM,PetscInt,DM*);
41: PetscErrorCode (*coarsenhierarchy)(DM,PetscInt,DM*);
42: PetscErrorCode (*adaptlabel)(DM,DMLabel,DM*);
43: PetscErrorCode (*adaptmetric)(DM,Vec,DMLabel,DM*);
45: PetscErrorCode (*globaltolocalbegin)(DM,Vec,InsertMode,Vec);
46: PetscErrorCode (*globaltolocalend)(DM,Vec,InsertMode,Vec);
47: PetscErrorCode (*localtoglobalbegin)(DM,Vec,InsertMode,Vec);
48: PetscErrorCode (*localtoglobalend)(DM,Vec,InsertMode,Vec);
49: PetscErrorCode (*localtolocalbegin)(DM,Vec,InsertMode,Vec);
50: PetscErrorCode (*localtolocalend)(DM,Vec,InsertMode,Vec);
52: PetscErrorCode (*destroy)(DM);
54: PetscErrorCode (*computevariablebounds)(DM,Vec,Vec);
56: PetscErrorCode (*createsubdm)(DM,PetscInt,const PetscInt*,IS*,DM*);
57: PetscErrorCode (*createsuperdm)(DM*,PetscInt,IS**,DM*);
58: PetscErrorCode (*createfielddecomposition)(DM,PetscInt*,char***,IS**,DM**);
59: PetscErrorCode (*createdomaindecomposition)(DM,PetscInt*,char***,IS**,IS**,DM**);
60: PetscErrorCode (*createddscatters)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
62: PetscErrorCode (*getdimpoints)(DM,PetscInt,PetscInt*,PetscInt*);
63: PetscErrorCode (*locatepoints)(DM,Vec,DMPointLocationType,PetscSF);
64: PetscErrorCode (*getneighbors)(DM,PetscInt*,const PetscMPIInt**);
65: PetscErrorCode (*getboundingbox)(DM,PetscReal*,PetscReal*);
66: PetscErrorCode (*getlocalboundingbox)(DM,PetscReal*,PetscReal*);
67: PetscErrorCode (*locatepointssubdomain)(DM,Vec,PetscMPIInt**);
69: PetscErrorCode (*projectfunctionlocal)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
70: PetscErrorCode (*projectfunctionlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
71: PetscErrorCode (*projectfieldlocal)(DM,PetscReal,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[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
72: PetscErrorCode (*projectfieldlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],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[],PetscInt,const PetscScalar[],PetscScalar[]),InsertMode,Vec);
73: PetscErrorCode (*computel2diff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
74: PetscErrorCode (*computel2gradientdiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal [],const PetscReal[],PetscInt, PetscScalar *,void *),void **,Vec,const PetscReal[],PetscReal *);
75: PetscErrorCode (*computel2fielddiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
76: };
78: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
79: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
80: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
82: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
83: struct _DMCoarsenHookLink {
84: PetscErrorCode (*coarsenhook)(DM,DM,void*); /* Run once, when coarse DM is created */
85: PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*); /* Run each time a new problem is restricted to a coarse grid */
86: void *ctx;
87: DMCoarsenHookLink next;
88: };
90: typedef struct _DMRefineHookLink *DMRefineHookLink;
91: struct _DMRefineHookLink {
92: PetscErrorCode (*refinehook)(DM,DM,void*); /* Run once, when a fine DM is created */
93: PetscErrorCode (*interphook)(DM,Mat,DM,void*); /* Run each time a new problem is interpolated to a fine grid */
94: void *ctx;
95: DMRefineHookLink next;
96: };
98: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
99: struct _DMSubDomainHookLink {
100: PetscErrorCode (*ddhook)(DM,DM,void*);
101: PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*);
102: void *ctx;
103: DMSubDomainHookLink next;
104: };
106: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
107: struct _DMGlobalToLocalHookLink {
108: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
109: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
110: void *ctx;
111: DMGlobalToLocalHookLink next;
112: };
114: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
115: struct _DMLocalToGlobalHookLink {
116: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
117: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
118: void *ctx;
119: DMLocalToGlobalHookLink next;
120: };
122: typedef enum {DMVEC_STATUS_IN,DMVEC_STATUS_OUT} DMVecStatus;
123: typedef struct _DMNamedVecLink *DMNamedVecLink;
124: struct _DMNamedVecLink {
125: Vec X;
126: char *name;
127: DMVecStatus status;
128: DMNamedVecLink next;
129: };
131: typedef struct _DMWorkLink *DMWorkLink;
132: struct _DMWorkLink {
133: size_t bytes;
134: void *mem;
135: DMWorkLink next;
136: };
138: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
140: struct _n_DMLabelLink {
141: DMLabel label;
142: PetscBool output;
143: struct _n_DMLabelLink *next;
144: };
145: typedef struct _n_DMLabelLink *DMLabelLink;
147: struct _n_DMLabelLinkList {
148: PetscInt refct;
149: DMLabelLink next;
150: };
151: typedef struct _n_DMLabelLinkList *DMLabelLinkList;
153: typedef struct _n_Boundary *DMBoundary;
155: struct _n_Boundary {
156: DSBoundary dsboundary;
157: DMLabel label;
158: DMBoundary next;
159: };
161: PETSC_EXTERN PetscErrorCode DMDestroyLabelLinkList(DM);
163: struct _p_DM {
164: PETSCHEADER(struct _DMOps);
165: Vec localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
166: Vec globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
167: DMNamedVecLink namedglobal;
168: DMNamedVecLink namedlocal;
169: DMWorkLink workin,workout;
170: DMLabelLinkList labels; /* Linked list of labels */
171: DMLabel depthLabel; /* Optimized access to depth label */
172: void *ctx; /* a user context */
173: PetscErrorCode (*ctxdestroy)(void**);
174: Vec x; /* location at which the functions/Jacobian are computed */
175: ISColoringType coloringtype;
176: MatFDColoring fd;
177: VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
178: MatType mattype; /* type of matrix created with DMCreateMatrix() */
179: PetscInt bs;
180: ISLocalToGlobalMapping ltogmap;
181: PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
182: PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
183: 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 */
184: 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 */
185: void *data;
186: /* Hierarchy / Submeshes */
187: DM coarseMesh;
188: DM fineMesh;
189: DMCoarsenHookLink coarsenhook; /* For transfering auxiliary problem data to coarser grids */
190: DMRefineHookLink refinehook;
191: DMSubDomainHookLink subdomainhook;
192: DMGlobalToLocalHookLink gtolhook;
193: DMLocalToGlobalHookLink ltoghook;
194: /* Topology */
195: PetscInt dim; /* The topological dimension */
196: /* Flexible communication */
197: PetscSF sfMigration; /* SF for point distribution created during distribution */
198: PetscSF sf; /* SF for parallel point overlap */
199: PetscSF defaultSF; /* SF for parallel dof overlap using default section */
200: PetscSF sfNatural; /* SF mapping to the "natural" ordering */
201: PetscBool useNatural; /* Create the natural SF */
202: /* Allows a non-standard data layout */
203: PetscSection defaultSection; /* Layout for local vectors */
204: PetscSection defaultGlobalSection; /* Layout for global vectors */
205: PetscLayout map;
206: /* Constraints */
207: PetscSection defaultConstraintSection;
208: Mat defaultConstraintMat;
209: /* Coordinates */
210: PetscInt dimEmbed; /* The dimension of the embedding space */
211: DM coordinateDM; /* Layout for coordinates (default section) */
212: Vec coordinates; /* Coordinate values in global vector */
213: Vec coordinatesLocal; /* Coordinate values in local vector */
214: PetscBool periodic; /* Is the DM periodic? */
215: PetscReal *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
216: DMBoundaryType *bdtype; /* Indicates type of topological boundary */
217: /* Null spaces -- of course I should make this have a variable number of fields */
218: /* I now believe this might not be the right way: see below */
219: NullSpaceFunc nullspaceConstructors[10];
220: /* Fields are represented by objects */
221: PetscDS prob;
222: DMBoundary boundary; /* List of boundary conditions */
223: /* Output structures */
224: DM dmBC; /* The DM with boundary conditions in the global DM */
225: PetscInt outputSequenceNum; /* The current sequence number for output */
226: PetscReal outputSequenceVal; /* The current sequence value for output */
228: PetscObject dmksp,dmsnes,dmts;
229: };
231: PETSC_EXTERN PetscLogEvent DM_Convert, DM_GlobalToLocal, DM_LocalToGlobal, DM_LocatePoints, DM_Coarsen, DM_Refine, DM_CreateInterpolation, DM_CreateRestriction;
233: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
234: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);
235: PETSC_EXTERN PetscErrorCode DMCreateSubDM_Section_Private(DM,PetscInt,const PetscInt[],IS*,DM*);
236: PETSC_EXTERN PetscErrorCode DMCreateSuperDM_Section_Private(DM[],PetscInt,IS**,DM*);
238: /*
240: Composite Vectors
242: Single global representation
243: Individual global representations
244: Single local representation
245: Individual local representations
247: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
249: DM da_u, da_v, da_p
251: DM dm_velocities
253: DM dm
255: DMDACreate(,&da_u);
256: DMDACreate(,&da_v);
257: DMCompositeCreate(,&dm_velocities);
258: DMCompositeAddDM(dm_velocities,(DM)du);
259: DMCompositeAddDM(dm_velocities,(DM)dv);
261: DMDACreate(,&da_p);
262: DMCompositeCreate(,&dm_velocities);
263: DMCompositeAddDM(dm,(DM)dm_velocities);
264: DMCompositeAddDM(dm,(DM)dm_p);
267: Access parts of composite vectors (Composite only)
268: ---------------------------------
269: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
270: ADD for local vector -
272: Element access
273: --------------
274: From global vectors
275: -DAVecGetArray - for DMDA
276: -VecGetArray - for DMSliced
277: ADD for DMComposite??? maybe
279: From individual vector
280: -DAVecGetArray - for DMDA
281: -VecGetArray -for sliced
282: ADD for DMComposite??? maybe
284: From single local vector
285: ADD * single local vector as arrays?
287: Communication
288: -------------
289: DMGlobalToLocal - global vector to single local vector
291: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
293: Obtaining vectors
294: -----------------
295: DMCreateGlobal/Local
296: DMGetGlobal/Local
297: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
300: ????? individual global vectors ????
302: */
304: #if defined(PETSC_HAVE_HDF5)
305: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
306: #endif
308: #endif