Actual source code: dmimpl.h
1: #pragma once
3: #include <petscdm.h>
4: #ifdef PETSC_HAVE_LIBCEED
5: #include <petscdmceed.h>
6: #endif
7: #include <petsc/private/petscimpl.h>
8: #include <petsc/private/petscdsimpl.h>
9: #include <petsc/private/sectionimpl.h>
11: PETSC_EXTERN PetscBool DMRegisterAllCalled;
12: PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
13: typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace);
15: typedef struct _PetscHashAuxKey {
16: DMLabel label;
17: PetscInt value;
18: PetscInt part;
19: } PetscHashAuxKey;
21: #define PetscHashAuxKeyHash(key) PetscHashCombine(PetscHashCombine(PetscHashPointer((key).label), PetscHashInt((key).value)), PetscHashInt((key).part))
23: #define PetscHashAuxKeyEqual(k1, k2) (((k1).label == (k2).label) ? (((k1).value == (k2).value) ? ((k1).part == (k2).part) : 0) : 0)
25: PETSC_HASH_MAP(HMapAux, PetscHashAuxKey, Vec, PetscHashAuxKeyHash, PetscHashAuxKeyEqual, PETSC_NULLPTR)
27: struct _n_DMGeneratorFunctionList {
28: PetscErrorCode (*generate)(DM, PetscBool, DM *);
29: PetscErrorCode (*refine)(DM, PetscReal *, DM *);
30: PetscErrorCode (*adapt)(DM, Vec, DMLabel, DMLabel, DM *);
31: char *name;
32: PetscInt dim;
33: DMGeneratorFunctionList next;
34: };
36: typedef struct _DMOps *DMOps;
37: struct _DMOps {
38: PetscErrorCode (*view)(DM, PetscViewer);
39: PetscErrorCode (*load)(DM, PetscViewer);
40: PetscErrorCode (*clone)(DM, DM *);
41: PetscErrorCode (*setfromoptions)(DM, PetscOptionItems *);
42: PetscErrorCode (*setup)(DM);
43: PetscErrorCode (*createlocalsection)(DM);
44: PetscErrorCode (*createsectionpermutation)(DM, IS *, PetscBT *);
45: PetscErrorCode (*createdefaultconstraints)(DM);
46: PetscErrorCode (*createglobalvector)(DM, Vec *);
47: PetscErrorCode (*createlocalvector)(DM, Vec *);
48: PetscErrorCode (*getlocaltoglobalmapping)(DM);
49: PetscErrorCode (*createfieldis)(DM, PetscInt *, char ***, IS **);
50: PetscErrorCode (*createcoordinatedm)(DM, DM *);
51: PetscErrorCode (*createcoordinatefield)(DM, DMField *);
53: PetscErrorCode (*getcoloring)(DM, ISColoringType, ISColoring *);
54: PetscErrorCode (*creatematrix)(DM, Mat *);
55: PetscErrorCode (*createinterpolation)(DM, DM, Mat *, Vec *);
56: PetscErrorCode (*createrestriction)(DM, DM, Mat *);
57: PetscErrorCode (*createmassmatrix)(DM, DM, Mat *);
58: PetscErrorCode (*createmassmatrixlumped)(DM, Vec *, Vec *);
59: PetscErrorCode (*hascreateinjection)(DM, PetscBool *);
60: PetscErrorCode (*createinjection)(DM, DM, Mat *);
62: PetscErrorCode (*refine)(DM, MPI_Comm, DM *);
63: PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *);
64: PetscErrorCode (*refinehierarchy)(DM, PetscInt, DM *);
65: PetscErrorCode (*coarsenhierarchy)(DM, PetscInt, DM *);
66: PetscErrorCode (*extrude)(DM, PetscInt, DM *);
68: PetscErrorCode (*globaltolocalbegin)(DM, Vec, InsertMode, Vec);
69: PetscErrorCode (*globaltolocalend)(DM, Vec, InsertMode, Vec);
70: PetscErrorCode (*localtoglobalbegin)(DM, Vec, InsertMode, Vec);
71: PetscErrorCode (*localtoglobalend)(DM, Vec, InsertMode, Vec);
72: PetscErrorCode (*localtolocalbegin)(DM, Vec, InsertMode, Vec);
73: PetscErrorCode (*localtolocalend)(DM, Vec, InsertMode, Vec);
75: PetscErrorCode (*destroy)(DM);
77: PetscErrorCode (*computevariablebounds)(DM, Vec, Vec);
79: PetscErrorCode (*createsubdm)(DM, PetscInt, const PetscInt *, IS *, DM *);
80: PetscErrorCode (*createsuperdm)(DM *, PetscInt, IS **, DM *);
81: PetscErrorCode (*createfielddecomposition)(DM, PetscInt *, char ***, IS **, DM **);
82: PetscErrorCode (*createdomaindecomposition)(DM, PetscInt *, char ***, IS **, IS **, DM **);
83: PetscErrorCode (*createddscatters)(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
85: PetscErrorCode (*getdimpoints)(DM, PetscInt, PetscInt *, PetscInt *);
86: PetscErrorCode (*locatepoints)(DM, Vec, DMPointLocationType, PetscSF);
87: PetscErrorCode (*getneighbors)(DM, PetscInt *, const PetscMPIInt **);
88: PetscErrorCode (*getboundingbox)(DM, PetscReal *, PetscReal *);
89: PetscErrorCode (*getlocalboundingbox)(DM, PetscReal[], PetscReal[], PetscInt[], PetscInt[]);
90: PetscErrorCode (*locatepointssubdomain)(DM, Vec, PetscMPIInt **);
91: PetscErrorCode (*snaptogeommodel)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
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: PetscDS dsIn; /* Approximation space for input to this domain (now only used for cohesive cells) */
197: DMLabel label; /* Label defining the domain of definition of the discretization */
198: IS fields; /* Map from DS field numbers to original field numbers in the DM */
199: } DMSpace;
201: struct _p_UniversalLabel {
202: DMLabel label; /* The universal label */
203: PetscInt Nl; /* Number of labels encoded */
204: char **names; /* The label names */
205: PetscInt *indices; /* The original indices in the input DM */
206: PetscInt Nv; /* Total number of values in all the labels */
207: PetscInt *bits; /* Starting bit for values of each label */
208: PetscInt *masks; /* Masks to pull out label value bits for each label */
209: PetscInt *offsets; /* Starting offset for label values for each label */
210: PetscInt *values; /* Original label values before renumbering */
211: };
213: PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);
215: #define MAXDMMONITORS 5
217: typedef struct {
218: PetscInt dim; /* The dimension of the embedding space */
219: DM dm; /* Layout for coordinates */
220: Vec x; /* Coordinate values in global vector */
221: Vec xl; /* Coordinate values in local vector */
222: DMField field; /* Coordinates as an abstract field */
223: } DMCoordinates;
225: struct _p_DM {
226: PETSCHEADER(struct _DMOps);
227: Vec localin[DM_MAX_WORK_VECTORS], localout[DM_MAX_WORK_VECTORS];
228: Vec globalin[DM_MAX_WORK_VECTORS], globalout[DM_MAX_WORK_VECTORS];
229: DMNamedVecLink namedglobal;
230: DMNamedVecLink namedlocal;
231: DMWorkLink workin, workout;
232: DMLabelLink labels; /* Linked list of labels */
233: DMLabel depthLabel; /* Optimized access to depth label */
234: DMLabel celltypeLabel; /* Optimized access to celltype label */
235: void *ctx; /* a user context */
236: PetscErrorCode (*ctxdestroy)(void **);
237: ISColoringType coloringtype;
238: MatFDColoring fd;
239: VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
240: MatType mattype; /* type of matrix created with DMCreateMatrix() */
241: PetscInt bind_below; /* Local size threshold (in entries/rows) below which Vec/Mat objects are bound to CPU */
242: PetscInt bs;
243: DMBlockingType blocking_type;
244: ISLocalToGlobalMapping ltogmap;
245: PetscBool prealloc_skip; // Flag indicating the DMCreateMatrix() should not preallocate (only set sizes and local-to-global)
246: PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
247: PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
248: 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 */
249: 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 */
250: PetscBool setfromoptionscalled;
251: void *data;
252: /* Hierarchy / Submeshes */
253: DM coarseMesh;
254: DM fineMesh;
255: DMCoarsenHookLink coarsenhook; /* For transferring auxiliary problem data to coarser grids */
256: DMRefineHookLink refinehook;
257: DMSubDomainHookLink subdomainhook;
258: DMGlobalToLocalHookLink gtolhook;
259: DMLocalToGlobalHookLink ltoghook;
260: /* Topology */
261: PetscInt dim; /* The topological dimension */
262: /* Auxiliary data */
263: PetscHMapAux auxData; /* Auxiliary DM and Vec for region denoted by the key */
264: /* Flexible communication */
265: PetscSF sfMigration; /* SF for point distribution created during distribution */
266: PetscSF sf; /* SF for parallel point overlap */
267: PetscSF sectionSF; /* SF for parallel dof overlap using section */
268: PetscSF sfNatural; /* SF mapping to the "natural" ordering */
269: PetscBool useNatural; /* Create the natural SF */
270: /* Allows a non-standard data layout */
271: PetscBool adjacency[2]; /* [use cone() or support() first, use the transitive closure] */
272: PetscSection localSection; /* Layout for local vectors */
273: PetscSection globalSection; /* Layout for global vectors */
274: PetscLayout map; /* Parallel division of unknowns across processes */
275: DMReorderDefaultFlag reorderSection; /* Reorder the local section by default */
276: MatOrderingType reorderSectionType; /* The type of reordering */
278: // Affine transform applied in DMGlobalToLocal
279: struct {
280: PetscInt num_affines;
281: VecScatter *affine_to_local;
282: Vec *affine;
283: PetscErrorCode (*setup)(DM);
284: } periodic;
285: /* Constraints */
286: struct {
287: PetscSection section;
288: Mat mat;
289: Vec bias;
290: } defaultConstraint;
291: /* Basis transformation */
292: DM transformDM; /* Layout for basis transformation */
293: Vec transform; /* Basis transformation matrices */
294: void *transformCtx; /* Basis transformation context */
295: PetscErrorCode (*transformSetUp)(DM, void *);
296: PetscErrorCode (*transformDestroy)(DM, void *);
297: PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
298: /* Coordinates */
299: DMCoordinates coordinates[2]; /* Coordinates, 0 is default, 1 is possible DG coordinate field for periodicity */
300: /* Periodicity */
301: PetscReal *Lstart, *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
302: PetscBool sparseLocalize; /* Localize coordinates only for cells near periodic boundary */
303: /* Null spaces -- of course I should make this have a variable number of fields */
304: NullSpaceFunc nullspaceConstructors[10];
305: NullSpaceFunc nearnullspaceConstructors[10];
306: /* Fields are represented by objects */
307: PetscInt Nf; /* Number of fields defined on the total domain */
308: RegionField *fields; /* Array of discretization fields with regions of validity */
309: DMBoundary boundary; /* List of boundary conditions */
310: PetscInt Nds; /* Number of discrete systems defined on the total domain */
311: DMSpace *probs; /* Array of discrete systems */
312: /* Output structures */
313: DM dmBC; /* The DM with boundary conditions in the global DM */
314: PetscBool ignorePermOutput; /* Ignore the local section permutation on output */
315: PetscInt outputSequenceNum; /* The current sequence number for output */
316: PetscReal outputSequenceVal; /* The current sequence value for output */
317: PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
318: PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
319: void *monitorcontext[MAXDMMONITORS];
320: PetscInt numbermonitors;
321: /* Configuration */
322: 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 */
324: PetscObject dmksp, dmsnes, dmts;
325: #ifdef PETSC_HAVE_LIBCEED
326: Ceed ceed; // LibCEED context
327: CeedElemRestriction ceedERestrict; // Map from the local vector (Lvector) to the cells (Evector)
328: #endif
329: DMCeed dmceed; // CEED operator and data for this problem
330: };
332: PETSC_EXTERN PetscLogEvent DM_Convert;
333: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
334: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
335: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
336: PETSC_EXTERN PetscLogEvent DM_Coarsen;
337: PETSC_EXTERN PetscLogEvent DM_Refine;
338: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
339: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
340: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
341: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
342: PETSC_EXTERN PetscLogEvent DM_CreateMassMatrix;
343: PETSC_EXTERN PetscLogEvent DM_Load;
344: PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
345: PETSC_EXTERN PetscLogEvent DM_ProjectFunction;
347: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM, Vec *);
348: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM, Vec *);
350: PETSC_EXTERN PetscErrorCode DMView_GLVis(DM, PetscViewer, PetscErrorCode (*)(DM, PetscViewer));
352: /*
354: Composite Vectors
356: Single global representation
357: Individual global representations
358: Single local representation
359: Individual local representations
361: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
363: DM da_u, da_v, da_p
365: DM dm_velocities
367: DM dm
369: DMDACreate(,&da_u);
370: DMDACreate(,&da_v);
371: DMCompositeCreate(,&dm_velocities);
372: DMCompositeAddDM(dm_velocities,(DM)du);
373: DMCompositeAddDM(dm_velocities,(DM)dv);
375: DMDACreate(,&da_p);
376: DMCompositeCreate(,&dm_velocities);
377: DMCompositeAddDM(dm,(DM)dm_velocities);
378: DMCompositeAddDM(dm,(DM)dm_p);
380: Access parts of composite vectors (Composite only)
381: ---------------------------------
382: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
383: ADD for local vector -
385: Element access
386: --------------
387: From global vectors
388: -DAVecGetArray - for DMDA
389: -VecGetArray - for DMSliced
390: ADD for DMComposite??? maybe
392: From individual vector
393: -DAVecGetArray - for DMDA
394: -VecGetArray -for sliced
395: ADD for DMComposite??? maybe
397: From single local vector
398: ADD * single local vector as arrays?
400: Communication
401: -------------
402: DMGlobalToLocal - global vector to single local vector
404: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
406: Obtaining vectors
407: -----------------
408: DMCreateGlobal/Local
409: DMGetGlobal/Local
410: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
412: ????? individual global vectors ????
414: */
416: #if defined(PETSC_HAVE_HDF5)
417: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char[], PetscInt, PetscScalar *, PetscViewer);
418: PETSC_EXTERN PetscErrorCode DMSequenceGetLength_HDF5_Internal(DM, const char[], PetscInt *, PetscViewer);
419: #endif
421: static inline PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
422: {
423: PetscFunctionBeginHot;
424: #if defined(PETSC_USE_DEBUG)
425: {
426: PetscInt dof;
428: *start = *end = 0; /* Silence overzealous compiler warning */
429: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
430: PetscCall(PetscSectionGetOffset(dm->localSection, point, start));
431: PetscCall(PetscSectionGetDof(dm->localSection, point, &dof));
432: *end = *start + dof;
433: }
434: #else
435: {
436: const PetscSection s = dm->localSection;
437: *start = s->atlasOff[point - s->pStart];
438: *end = *start + s->atlasDof[point - s->pStart];
439: }
440: #endif
441: PetscFunctionReturn(PETSC_SUCCESS);
442: }
444: static inline PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
445: {
446: PetscFunctionBegin;
447: #if defined(PETSC_USE_DEBUG)
448: {
449: PetscInt dof;
450: *start = *end = 0; /* Silence overzealous compiler warning */
451: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
452: PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, start));
453: PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &dof));
454: *end = *start + dof;
455: }
456: #else
457: {
458: const PetscSection s = dm->localSection->field[field];
459: *start = s->atlasOff[point - s->pStart];
460: *end = *start + s->atlasDof[point - s->pStart];
461: }
462: #endif
463: PetscFunctionReturn(PETSC_SUCCESS);
464: }
466: static inline PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
467: {
468: PetscFunctionBegin;
469: #if defined(PETSC_USE_DEBUG)
470: {
471: PetscInt dof, cdof;
472: *start = *end = 0; /* Silence overzealous compiler warning */
473: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
474: PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
475: PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
476: PetscCall(PetscSectionGetDof(dm->globalSection, point, &dof));
477: PetscCall(PetscSectionGetConstraintDof(dm->globalSection, point, &cdof));
478: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
479: }
480: #else
481: {
482: const PetscSection s = dm->globalSection;
483: const PetscInt dof = s->atlasDof[point - s->pStart];
484: const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
485: *start = s->atlasOff[point - s->pStart];
486: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
487: }
488: #endif
489: PetscFunctionReturn(PETSC_SUCCESS);
490: }
492: static inline PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
493: {
494: PetscFunctionBegin;
495: #if defined(PETSC_USE_DEBUG)
496: {
497: PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
498: *start = *end = 0; /* Silence overzealous compiler warning */
499: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
500: PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
501: PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
502: PetscCall(PetscSectionGetOffset(dm->localSection, point, &loff));
503: PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff));
504: PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &fdof));
505: PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof));
506: *start = *start < 0 ? *start - (lfoff - loff) : *start + lfoff - loff;
507: for (f = 0; f < field; ++f) {
508: PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof));
509: *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
510: }
511: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
512: }
513: #else
514: {
515: const PetscSection s = dm->localSection;
516: const PetscSection fs = dm->localSection->field[field];
517: const PetscSection gs = dm->globalSection;
518: const PetscInt loff = s->atlasOff[point - s->pStart];
519: const PetscInt goff = gs->atlasOff[point - s->pStart];
520: const PetscInt lfoff = fs->atlasOff[point - s->pStart];
521: const PetscInt fdof = fs->atlasDof[point - s->pStart];
522: const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
523: PetscInt ffcdof = 0, f;
525: for (f = 0; f < field; ++f) {
526: const PetscSection ffs = dm->localSection->field[f];
527: ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
528: }
529: *start = goff + (goff < 0 ? loff - lfoff + ffcdof : lfoff - loff - ffcdof);
530: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
531: }
532: #endif
533: PetscFunctionReturn(PETSC_SUCCESS);
534: }
536: PETSC_EXTERN PetscErrorCode DMGetCoordinateDegree_Internal(DM, PetscInt *);
537: PETSC_INTERN PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM, PetscReal[], PetscReal[], PetscInt[], PetscInt[]);
539: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
540: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
541: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);
543: PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
544: PETSC_INTERN PetscErrorCode DMSetField_Internal(DM, PetscInt, DMLabel, PetscObject);
546: PETSC_INTERN PetscErrorCode DMSetLabelValue_Fast(DM, DMLabel *, const char[], PetscInt, PetscInt);
547: PETSC_INTERN PetscErrorCode DMGetPoints_Internal(DM, DMLabel, PetscInt, PetscInt, IS *);
549: PETSC_INTERN PetscErrorCode DMCompleteBCLabels_Internal(DM dm);
550: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
551: PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
552: PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
553: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
554: PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);
555: PETSC_INTERN PetscInt PetscGCD(PetscInt a, PetscInt b);