Actual source code: dmimpl.h
2: #ifndef _DMIMPL_H
3: #define _DMIMPL_H
5: #include <petscdm.h>
6: #ifdef PETSC_HAVE_LIBCEED
7: #include <petscdmceed.h>
8: #endif
9: #include <petsc/private/petscimpl.h>
10: #include <petsc/private/petscdsimpl.h>
11: #include <petsc/private/sectionimpl.h>
13: PETSC_EXTERN PetscBool DMRegisterAllCalled;
14: PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
15: typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace);
17: typedef struct _PetscHashAuxKey {
18: DMLabel label;
19: PetscInt value;
20: PetscInt part;
21: } PetscHashAuxKey;
23: #define PetscHashAuxKeyHash(key) PetscHashCombine(PetscHashCombine(PetscHashPointer((key).label), PetscHashInt((key).value)), PetscHashInt((key).part))
25: #define PetscHashAuxKeyEqual(k1, k2) (((k1).label == (k2).label) ? (((k1).value == (k2).value) ? ((k1).part == (k2).part) : 0) : 0)
27: PETSC_HASH_MAP(HMapAux, PetscHashAuxKey, Vec, PetscHashAuxKeyHash, PetscHashAuxKeyEqual, NULL)
29: struct _n_DMGeneratorFunctionList {
30: PetscErrorCode (*generate)(DM, PetscBool, DM *);
31: PetscErrorCode (*refine)(DM, PetscReal *, DM *);
32: PetscErrorCode (*adapt)(DM, Vec, DMLabel, DMLabel, DM *);
33: char *name;
34: PetscInt dim;
35: DMGeneratorFunctionList next;
36: };
38: typedef struct _DMOps *DMOps;
39: struct _DMOps {
40: PetscErrorCode (*view)(DM, PetscViewer);
41: PetscErrorCode (*load)(DM, PetscViewer);
42: PetscErrorCode (*clone)(DM, DM *);
43: PetscErrorCode (*setfromoptions)(DM, PetscOptionItems *);
44: PetscErrorCode (*setup)(DM);
45: PetscErrorCode (*createlocalsection)(DM);
46: PetscErrorCode (*createdefaultconstraints)(DM);
47: PetscErrorCode (*createglobalvector)(DM, Vec *);
48: PetscErrorCode (*createlocalvector)(DM, Vec *);
49: PetscErrorCode (*getlocaltoglobalmapping)(DM);
50: PetscErrorCode (*createfieldis)(DM, PetscInt *, char ***, IS **);
51: PetscErrorCode (*createcoordinatedm)(DM, DM *);
52: PetscErrorCode (*createcoordinatefield)(DM, DMField *);
54: PetscErrorCode (*getcoloring)(DM, ISColoringType, ISColoring *);
55: PetscErrorCode (*creatematrix)(DM, Mat *);
56: PetscErrorCode (*createinterpolation)(DM, DM, Mat *, Vec *);
57: PetscErrorCode (*createrestriction)(DM, DM, Mat *);
58: PetscErrorCode (*createmassmatrix)(DM, DM, Mat *);
59: PetscErrorCode (*createmassmatrixlumped)(DM, Vec *);
60: PetscErrorCode (*hascreateinjection)(DM, PetscBool *);
61: PetscErrorCode (*createinjection)(DM, DM, Mat *);
63: PetscErrorCode (*refine)(DM, MPI_Comm, DM *);
64: PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *);
65: PetscErrorCode (*refinehierarchy)(DM, PetscInt, DM *);
66: PetscErrorCode (*coarsenhierarchy)(DM, PetscInt, DM *);
67: PetscErrorCode (*extrude)(DM, PetscInt, DM *);
69: PetscErrorCode (*globaltolocalbegin)(DM, Vec, InsertMode, Vec);
70: PetscErrorCode (*globaltolocalend)(DM, Vec, InsertMode, Vec);
71: PetscErrorCode (*localtoglobalbegin)(DM, Vec, InsertMode, Vec);
72: PetscErrorCode (*localtoglobalend)(DM, Vec, InsertMode, Vec);
73: PetscErrorCode (*localtolocalbegin)(DM, Vec, InsertMode, Vec);
74: PetscErrorCode (*localtolocalend)(DM, Vec, InsertMode, Vec);
76: PetscErrorCode (*destroy)(DM);
78: PetscErrorCode (*computevariablebounds)(DM, Vec, Vec);
80: PetscErrorCode (*createsubdm)(DM, PetscInt, const PetscInt *, IS *, DM *);
81: PetscErrorCode (*createsuperdm)(DM *, PetscInt, IS **, DM *);
82: PetscErrorCode (*createfielddecomposition)(DM, PetscInt *, char ***, IS **, DM **);
83: PetscErrorCode (*createdomaindecomposition)(DM, PetscInt *, char ***, IS **, IS **, DM **);
84: PetscErrorCode (*createddscatters)(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
86: PetscErrorCode (*getdimpoints)(DM, PetscInt, PetscInt *, PetscInt *);
87: PetscErrorCode (*locatepoints)(DM, Vec, DMPointLocationType, PetscSF);
88: PetscErrorCode (*getneighbors)(DM, PetscInt *, const PetscMPIInt **);
89: PetscErrorCode (*getboundingbox)(DM, PetscReal *, PetscReal *);
90: PetscErrorCode (*getlocalboundingbox)(DM, PetscReal *, PetscReal *);
91: PetscErrorCode (*locatepointssubdomain)(DM, Vec, PetscMPIInt **);
93: PetscErrorCode (*projectfunctionlocal)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
94: PetscErrorCode (*projectfunctionlabellocal)(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
95: 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);
96: 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);
97: PetscErrorCode (*projectbdfieldlabellocal)(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[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);
98: PetscErrorCode (*computel2diff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
99: PetscErrorCode (*computel2gradientdiff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
100: PetscErrorCode (*computel2fielddiff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
102: PetscErrorCode (*getcompatibility)(DM, DM, PetscBool *, PetscBool *);
103: };
105: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
106: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
107: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
109: PETSC_INTERN PetscErrorCode DMCountNonCyclicReferences(PetscObject, PetscInt *);
111: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
112: struct _DMCoarsenHookLink {
113: PetscErrorCode (*coarsenhook)(DM, DM, void *); /* Run once, when coarse DM is created */
114: PetscErrorCode (*restricthook)(DM, Mat, Vec, Mat, DM, void *); /* Run each time a new problem is restricted to a coarse grid */
115: void *ctx;
116: DMCoarsenHookLink next;
117: };
119: typedef struct _DMRefineHookLink *DMRefineHookLink;
120: struct _DMRefineHookLink {
121: PetscErrorCode (*refinehook)(DM, DM, void *); /* Run once, when a fine DM is created */
122: PetscErrorCode (*interphook)(DM, Mat, DM, void *); /* Run each time a new problem is interpolated to a fine grid */
123: void *ctx;
124: DMRefineHookLink next;
125: };
127: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
128: struct _DMSubDomainHookLink {
129: PetscErrorCode (*ddhook)(DM, DM, void *);
130: PetscErrorCode (*restricthook)(DM, VecScatter, VecScatter, DM, void *);
131: void *ctx;
132: DMSubDomainHookLink next;
133: };
135: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
136: struct _DMGlobalToLocalHookLink {
137: PetscErrorCode (*beginhook)(DM, Vec, InsertMode, Vec, void *);
138: PetscErrorCode (*endhook)(DM, Vec, InsertMode, Vec, void *);
139: void *ctx;
140: DMGlobalToLocalHookLink next;
141: };
143: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
144: struct _DMLocalToGlobalHookLink {
145: PetscErrorCode (*beginhook)(DM, Vec, InsertMode, Vec, void *);
146: PetscErrorCode (*endhook)(DM, Vec, InsertMode, Vec, void *);
147: void *ctx;
148: DMLocalToGlobalHookLink next;
149: };
151: typedef enum {
152: DMVEC_STATUS_IN,
153: DMVEC_STATUS_OUT
154: } DMVecStatus;
155: typedef struct _DMNamedVecLink *DMNamedVecLink;
156: struct _DMNamedVecLink {
157: Vec X;
158: char *name;
159: DMVecStatus status;
160: DMNamedVecLink next;
161: };
163: typedef struct _DMWorkLink *DMWorkLink;
164: struct _DMWorkLink {
165: size_t bytes;
166: void *mem;
167: DMWorkLink next;
168: };
170: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
172: struct _n_DMLabelLink {
173: DMLabel label;
174: PetscBool output;
175: struct _n_DMLabelLink *next;
176: };
177: typedef struct _n_DMLabelLink *DMLabelLink;
179: typedef struct _n_Boundary *DMBoundary;
181: struct _n_Boundary {
182: DSBoundary dsboundary;
183: DMLabel label;
184: DMBoundary next;
185: };
187: typedef struct _n_Field {
188: PetscObject disc; /* Field discretization, or a PetscContainer with the field name */
189: DMLabel label; /* Label defining the domain of definition of the field */
190: PetscBool adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
191: PetscBool avoidTensor; /* Flag to avoid defining field over tensor cells */
192: } RegionField;
194: typedef struct _n_Space {
195: PetscDS ds; /* Approximation space in this domain */
196: DMLabel label; /* Label defining the domain of definition of the discretization */
197: IS fields; /* Map from DS field numbers to original field numbers in the DM */
198: } DMSpace;
200: struct _p_UniversalLabel {
201: DMLabel label; /* The universal label */
202: PetscInt Nl; /* Number of labels encoded */
203: char **names; /* The label names */
204: PetscInt *indices; /* The original indices in the input DM */
205: PetscInt Nv; /* Total number of values in all the labels */
206: PetscInt *bits; /* Starting bit for values of each label */
207: PetscInt *masks; /* Masks to pull out label value bits for each label */
208: PetscInt *offsets; /* Starting offset for label values for each label */
209: PetscInt *values; /* Original label values before renumbering */
210: };
212: PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);
214: #define MAXDMMONITORS 5
216: typedef struct {
217: PetscInt dim; /* The dimension of the embedding space */
218: DM dm; /* Layout for coordinates */
219: Vec x; /* Coordinate values in global vector */
220: Vec xl; /* Coordinate values in local vector */
221: DMField field; /* Coordinates as an abstract field */
222: } DMCoordinates;
224: struct _p_DM {
225: PETSCHEADER(struct _DMOps);
226: Vec localin[DM_MAX_WORK_VECTORS], localout[DM_MAX_WORK_VECTORS];
227: Vec globalin[DM_MAX_WORK_VECTORS], globalout[DM_MAX_WORK_VECTORS];
228: DMNamedVecLink namedglobal;
229: DMNamedVecLink namedlocal;
230: DMWorkLink workin, workout;
231: DMLabelLink labels; /* Linked list of labels */
232: DMLabel depthLabel; /* Optimized access to depth label */
233: DMLabel celltypeLabel; /* Optimized access to celltype label */
234: void *ctx; /* a user context */
235: PetscErrorCode (*ctxdestroy)(void **);
236: ISColoringType coloringtype;
237: MatFDColoring fd;
238: VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
239: MatType mattype; /* type of matrix created with DMCreateMatrix() */
240: PetscInt bind_below; /* Local size threshold (in entries/rows) below which Vec/Mat objects are bound to CPU */
241: PetscInt bs;
242: DMBlockingType blocking_type;
243: ISLocalToGlobalMapping ltogmap;
244: PetscBool prealloc_skip; // Flag indicating the DMCreateMatrix() should not preallocate (only set sizes and local-to-global)
245: PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
246: PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
247: 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 */
248: 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 */
249: PetscBool setfromoptionscalled;
250: void *data;
251: /* Hierarchy / Submeshes */
252: DM coarseMesh;
253: DM fineMesh;
254: DMCoarsenHookLink coarsenhook; /* For transferring auxiliary problem data to coarser grids */
255: DMRefineHookLink refinehook;
256: DMSubDomainHookLink subdomainhook;
257: DMGlobalToLocalHookLink gtolhook;
258: DMLocalToGlobalHookLink ltoghook;
259: /* Topology */
260: PetscInt dim; /* The topological dimension */
261: /* Auxiliary data */
262: PetscHMapAux auxData; /* Auxiliary DM and Vec for region denoted by the key */
263: /* Flexible communication */
264: PetscSF sfMigration; /* SF for point distribution created during distribution */
265: PetscSF sf; /* SF for parallel point overlap */
266: PetscSF sectionSF; /* SF for parallel dof overlap using section */
267: PetscSF sfNatural; /* SF mapping to the "natural" ordering */
268: PetscBool useNatural; /* Create the natural SF */
269: /* Allows a non-standard data layout */
270: PetscBool adjacency[2]; /* [use cone() or support() first, use the transitive closure] */
271: PetscSection localSection; /* Layout for local vectors */
272: PetscSection globalSection; /* Layout for global vectors */
273: PetscLayout map;
274: // Affine transform applied in DMGlobalToLocal
275: struct {
276: VecScatter affine_to_local;
277: Vec affine;
278: PetscErrorCode (*setup)(DM);
279: } periodic;
280: /* Constraints */
281: struct {
282: PetscSection section;
283: Mat mat;
284: Vec bias;
285: } defaultConstraint;
286: /* Basis transformation */
287: DM transformDM; /* Layout for basis transformation */
288: Vec transform; /* Basis transformation matrices */
289: void *transformCtx; /* Basis transformation context */
290: PetscErrorCode (*transformSetUp)(DM, void *);
291: PetscErrorCode (*transformDestroy)(DM, void *);
292: PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
293: /* Coordinates */
294: DMCoordinates coordinates[2]; /* Coordinates, 0 is default, 1 is possible DG coordinate field for periodicity */
295: /* Periodicity */
296: PetscReal *Lstart, *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
297: PetscBool sparseLocalize; /* Localize coordinates only for cells near periodic boundary */
298: /* Null spaces -- of course I should make this have a variable number of fields */
299: NullSpaceFunc nullspaceConstructors[10];
300: NullSpaceFunc nearnullspaceConstructors[10];
301: /* Fields are represented by objects */
302: PetscInt Nf; /* Number of fields defined on the total domain */
303: RegionField *fields; /* Array of discretization fields with regions of validity */
304: DMBoundary boundary; /* List of boundary conditions */
305: PetscInt Nds; /* Number of discrete systems defined on the total domain */
306: DMSpace *probs; /* Array of discrete systems */
307: /* Output structures */
308: DM dmBC; /* The DM with boundary conditions in the global DM */
309: PetscInt outputSequenceNum; /* The current sequence number for output */
310: PetscReal outputSequenceVal; /* The current sequence value for output */
311: PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
312: PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
313: void *monitorcontext[MAXDMMONITORS];
314: PetscInt numbermonitors;
315: /* Configuration */
316: PetscBool cloneOpts; /* Flag indicating that this is a linked clone and should not respond to some options. This is currently used to prevent transformations from also affecting the coordinate DM */
318: PetscObject dmksp, dmsnes, dmts;
319: #ifdef PETSC_HAVE_LIBCEED
320: Ceed ceed; /* LibCEED context */
321: CeedElemRestriction ceedERestrict; /* Map from the local vector (Lvector) to the cells (Evector) */
322: #endif
323: };
325: PETSC_EXTERN PetscLogEvent DM_Convert;
326: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
327: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
328: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
329: PETSC_EXTERN PetscLogEvent DM_Coarsen;
330: PETSC_EXTERN PetscLogEvent DM_Refine;
331: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
332: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
333: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
334: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
335: PETSC_EXTERN PetscLogEvent DM_CreateMassMatrix;
336: PETSC_EXTERN PetscLogEvent DM_Load;
337: PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
339: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM, Vec *);
340: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM, Vec *);
342: PETSC_EXTERN PetscErrorCode DMView_GLVis(DM, PetscViewer, PetscErrorCode (*)(DM, PetscViewer));
344: /*
346: Composite Vectors
348: Single global representation
349: Individual global representations
350: Single local representation
351: Individual local representations
353: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
355: DM da_u, da_v, da_p
357: DM dm_velocities
359: DM dm
361: DMDACreate(,&da_u);
362: DMDACreate(,&da_v);
363: DMCompositeCreate(,&dm_velocities);
364: DMCompositeAddDM(dm_velocities,(DM)du);
365: DMCompositeAddDM(dm_velocities,(DM)dv);
367: DMDACreate(,&da_p);
368: DMCompositeCreate(,&dm_velocities);
369: DMCompositeAddDM(dm,(DM)dm_velocities);
370: DMCompositeAddDM(dm,(DM)dm_p);
372: Access parts of composite vectors (Composite only)
373: ---------------------------------
374: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
375: ADD for local vector -
377: Element access
378: --------------
379: From global vectors
380: -DAVecGetArray - for DMDA
381: -VecGetArray - for DMSliced
382: ADD for DMComposite??? maybe
384: From individual vector
385: -DAVecGetArray - for DMDA
386: -VecGetArray -for sliced
387: ADD for DMComposite??? maybe
389: From single local vector
390: ADD * single local vector as arrays?
392: Communication
393: -------------
394: DMGlobalToLocal - global vector to single local vector
396: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
398: Obtaining vectors
399: -----------------
400: DMCreateGlobal/Local
401: DMGetGlobal/Local
402: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
404: ????? individual global vectors ????
406: */
408: #if defined(PETSC_HAVE_HDF5)
409: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
410: #endif
412: static inline PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
413: {
414: PetscFunctionBeginHot;
415: #if defined(PETSC_USE_DEBUG)
416: {
417: PetscInt dof;
419: *start = *end = 0; /* Silence overzealous compiler warning */
420: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
421: PetscCall(PetscSectionGetOffset(dm->localSection, point, start));
422: PetscCall(PetscSectionGetDof(dm->localSection, point, &dof));
423: *end = *start + dof;
424: }
425: #else
426: {
427: const PetscSection s = dm->localSection;
428: *start = s->atlasOff[point - s->pStart];
429: *end = *start + s->atlasDof[point - s->pStart];
430: }
431: #endif
432: PetscFunctionReturn(PETSC_SUCCESS);
433: }
435: static inline PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
436: {
437: PetscFunctionBegin;
438: #if defined(PETSC_USE_DEBUG)
439: {
440: PetscInt dof;
441: *start = *end = 0; /* Silence overzealous compiler warning */
442: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
443: PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, start));
444: PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &dof));
445: *end = *start + dof;
446: }
447: #else
448: {
449: const PetscSection s = dm->localSection->field[field];
450: *start = s->atlasOff[point - s->pStart];
451: *end = *start + s->atlasDof[point - s->pStart];
452: }
453: #endif
454: PetscFunctionReturn(PETSC_SUCCESS);
455: }
457: static inline PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
458: {
459: PetscFunctionBegin;
460: #if defined(PETSC_USE_DEBUG)
461: {
462: PetscInt dof, cdof;
463: *start = *end = 0; /* Silence overzealous compiler warning */
464: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
465: PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
466: PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
467: PetscCall(PetscSectionGetDof(dm->globalSection, point, &dof));
468: PetscCall(PetscSectionGetConstraintDof(dm->globalSection, point, &cdof));
469: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
470: }
471: #else
472: {
473: const PetscSection s = dm->globalSection;
474: const PetscInt dof = s->atlasDof[point - s->pStart];
475: const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
476: *start = s->atlasOff[point - s->pStart];
477: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
478: }
479: #endif
480: PetscFunctionReturn(PETSC_SUCCESS);
481: }
483: static inline PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
484: {
485: PetscFunctionBegin;
486: #if defined(PETSC_USE_DEBUG)
487: {
488: PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
489: *start = *end = 0; /* Silence overzealous compiler warning */
490: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
491: PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
492: PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
493: PetscCall(PetscSectionGetOffset(dm->localSection, point, &loff));
494: PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff));
495: PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &fdof));
496: PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof));
497: *start = *start < 0 ? *start - (lfoff - loff) : *start + lfoff - loff;
498: for (f = 0; f < field; ++f) {
499: PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof));
500: *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
501: }
502: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
503: }
504: #else
505: {
506: const PetscSection s = dm->localSection;
507: const PetscSection fs = dm->localSection->field[field];
508: const PetscSection gs = dm->globalSection;
509: const PetscInt loff = s->atlasOff[point - s->pStart];
510: const PetscInt goff = gs->atlasOff[point - s->pStart];
511: const PetscInt lfoff = fs->atlasOff[point - s->pStart];
512: const PetscInt fdof = fs->atlasDof[point - s->pStart];
513: const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
514: PetscInt ffcdof = 0, f;
516: for (f = 0; f < field; ++f) {
517: const PetscSection ffs = dm->localSection->field[f];
518: ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
519: }
520: *start = goff + (goff < 0 ? loff - lfoff + ffcdof : lfoff - loff - ffcdof);
521: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
522: }
523: #endif
524: PetscFunctionReturn(PETSC_SUCCESS);
525: }
527: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
528: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
529: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);
531: PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
532: PETSC_INTERN PetscErrorCode DMSetField_Internal(DM, PetscInt, DMLabel, PetscObject);
534: PETSC_INTERN PetscErrorCode DMSetLabelValue_Fast(DM, DMLabel *, const char[], PetscInt, PetscInt);
536: PETSC_INTERN PetscErrorCode DMCompleteBCLabels_Internal(DM dm);
537: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
538: PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
539: PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
540: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
541: PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);
542: PETSC_INTERN PetscInt PetscGCD(PetscInt a, PetscInt b);
544: #endif