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, PETSC_NULLPTR)
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: 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;
275: // Affine transform applied in DMGlobalToLocal
276: struct {
277: VecScatter affine_to_local;
278: Vec affine;
279: PetscErrorCode (*setup)(DM);
280: } periodic;
281: /* Constraints */
282: struct {
283: PetscSection section;
284: Mat mat;
285: Vec bias;
286: } defaultConstraint;
287: /* Basis transformation */
288: DM transformDM; /* Layout for basis transformation */
289: Vec transform; /* Basis transformation matrices */
290: void *transformCtx; /* Basis transformation context */
291: PetscErrorCode (*transformSetUp)(DM, void *);
292: PetscErrorCode (*transformDestroy)(DM, void *);
293: PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
294: /* Coordinates */
295: DMCoordinates coordinates[2]; /* Coordinates, 0 is default, 1 is possible DG coordinate field for periodicity */
296: /* Periodicity */
297: PetscReal *Lstart, *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
298: PetscBool sparseLocalize; /* Localize coordinates only for cells near periodic boundary */
299: /* Null spaces -- of course I should make this have a variable number of fields */
300: NullSpaceFunc nullspaceConstructors[10];
301: NullSpaceFunc nearnullspaceConstructors[10];
302: /* Fields are represented by objects */
303: PetscInt Nf; /* Number of fields defined on the total domain */
304: RegionField *fields; /* Array of discretization fields with regions of validity */
305: DMBoundary boundary; /* List of boundary conditions */
306: PetscInt Nds; /* Number of discrete systems defined on the total domain */
307: DMSpace *probs; /* Array of discrete systems */
308: /* Output structures */
309: DM dmBC; /* The DM with boundary conditions in the global DM */
310: PetscInt outputSequenceNum; /* The current sequence number for output */
311: PetscReal outputSequenceVal; /* The current sequence value for output */
312: PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
313: PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
314: void *monitorcontext[MAXDMMONITORS];
315: PetscInt numbermonitors;
316: /* Configuration */
317: 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 */
319: PetscObject dmksp, dmsnes, dmts;
320: #ifdef PETSC_HAVE_LIBCEED
321: Ceed ceed; /* LibCEED context */
322: CeedElemRestriction ceedERestrict; /* Map from the local vector (Lvector) to the cells (Evector) */
323: #endif
324: };
326: PETSC_EXTERN PetscLogEvent DM_Convert;
327: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
328: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
329: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
330: PETSC_EXTERN PetscLogEvent DM_Coarsen;
331: PETSC_EXTERN PetscLogEvent DM_Refine;
332: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
333: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
334: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
335: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
336: PETSC_EXTERN PetscLogEvent DM_CreateMassMatrix;
337: PETSC_EXTERN PetscLogEvent DM_Load;
338: PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
340: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM, Vec *);
341: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM, Vec *);
343: PETSC_EXTERN PetscErrorCode DMView_GLVis(DM, PetscViewer, PetscErrorCode (*)(DM, PetscViewer));
345: /*
347: Composite Vectors
349: Single global representation
350: Individual global representations
351: Single local representation
352: Individual local representations
354: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
356: DM da_u, da_v, da_p
358: DM dm_velocities
360: DM dm
362: DMDACreate(,&da_u);
363: DMDACreate(,&da_v);
364: DMCompositeCreate(,&dm_velocities);
365: DMCompositeAddDM(dm_velocities,(DM)du);
366: DMCompositeAddDM(dm_velocities,(DM)dv);
368: DMDACreate(,&da_p);
369: DMCompositeCreate(,&dm_velocities);
370: DMCompositeAddDM(dm,(DM)dm_velocities);
371: DMCompositeAddDM(dm,(DM)dm_p);
373: Access parts of composite vectors (Composite only)
374: ---------------------------------
375: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
376: ADD for local vector -
378: Element access
379: --------------
380: From global vectors
381: -DAVecGetArray - for DMDA
382: -VecGetArray - for DMSliced
383: ADD for DMComposite??? maybe
385: From individual vector
386: -DAVecGetArray - for DMDA
387: -VecGetArray -for sliced
388: ADD for DMComposite??? maybe
390: From single local vector
391: ADD * single local vector as arrays?
393: Communication
394: -------------
395: DMGlobalToLocal - global vector to single local vector
397: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
399: Obtaining vectors
400: -----------------
401: DMCreateGlobal/Local
402: DMGetGlobal/Local
403: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
405: ????? individual global vectors ????
407: */
409: #if defined(PETSC_HAVE_HDF5)
410: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
411: #endif
413: static inline PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
414: {
415: PetscFunctionBeginHot;
416: #if defined(PETSC_USE_DEBUG)
417: {
418: PetscInt dof;
420: *start = *end = 0; /* Silence overzealous compiler warning */
421: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
422: PetscCall(PetscSectionGetOffset(dm->localSection, point, start));
423: PetscCall(PetscSectionGetDof(dm->localSection, point, &dof));
424: *end = *start + dof;
425: }
426: #else
427: {
428: const PetscSection s = dm->localSection;
429: *start = s->atlasOff[point - s->pStart];
430: *end = *start + s->atlasDof[point - s->pStart];
431: }
432: #endif
433: PetscFunctionReturn(PETSC_SUCCESS);
434: }
436: static inline PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
437: {
438: PetscFunctionBegin;
439: #if defined(PETSC_USE_DEBUG)
440: {
441: PetscInt dof;
442: *start = *end = 0; /* Silence overzealous compiler warning */
443: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
444: PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, start));
445: PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &dof));
446: *end = *start + dof;
447: }
448: #else
449: {
450: const PetscSection s = dm->localSection->field[field];
451: *start = s->atlasOff[point - s->pStart];
452: *end = *start + s->atlasDof[point - s->pStart];
453: }
454: #endif
455: PetscFunctionReturn(PETSC_SUCCESS);
456: }
458: static inline PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
459: {
460: PetscFunctionBegin;
461: #if defined(PETSC_USE_DEBUG)
462: {
463: PetscInt dof, cdof;
464: *start = *end = 0; /* Silence overzealous compiler warning */
465: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
466: PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
467: PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
468: PetscCall(PetscSectionGetDof(dm->globalSection, point, &dof));
469: PetscCall(PetscSectionGetConstraintDof(dm->globalSection, point, &cdof));
470: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
471: }
472: #else
473: {
474: const PetscSection s = dm->globalSection;
475: const PetscInt dof = s->atlasDof[point - s->pStart];
476: const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
477: *start = s->atlasOff[point - s->pStart];
478: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
479: }
480: #endif
481: PetscFunctionReturn(PETSC_SUCCESS);
482: }
484: static inline PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
485: {
486: PetscFunctionBegin;
487: #if defined(PETSC_USE_DEBUG)
488: {
489: PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
490: *start = *end = 0; /* Silence overzealous compiler warning */
491: PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
492: PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
493: PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
494: PetscCall(PetscSectionGetOffset(dm->localSection, point, &loff));
495: PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff));
496: PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &fdof));
497: PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof));
498: *start = *start < 0 ? *start - (lfoff - loff) : *start + lfoff - loff;
499: for (f = 0; f < field; ++f) {
500: PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof));
501: *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
502: }
503: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
504: }
505: #else
506: {
507: const PetscSection s = dm->localSection;
508: const PetscSection fs = dm->localSection->field[field];
509: const PetscSection gs = dm->globalSection;
510: const PetscInt loff = s->atlasOff[point - s->pStart];
511: const PetscInt goff = gs->atlasOff[point - s->pStart];
512: const PetscInt lfoff = fs->atlasOff[point - s->pStart];
513: const PetscInt fdof = fs->atlasDof[point - s->pStart];
514: const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
515: PetscInt ffcdof = 0, f;
517: for (f = 0; f < field; ++f) {
518: const PetscSection ffs = dm->localSection->field[f];
519: ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
520: }
521: *start = goff + (goff < 0 ? loff - lfoff + ffcdof : lfoff - loff - ffcdof);
522: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
523: }
524: #endif
525: PetscFunctionReturn(PETSC_SUCCESS);
526: }
528: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
529: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
530: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);
532: PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
533: PETSC_INTERN PetscErrorCode DMSetField_Internal(DM, PetscInt, DMLabel, PetscObject);
535: PETSC_INTERN PetscErrorCode DMSetLabelValue_Fast(DM, DMLabel *, const char[], PetscInt, PetscInt);
537: PETSC_INTERN PetscErrorCode DMCompleteBCLabels_Internal(DM dm);
538: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
539: PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
540: PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
541: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
542: PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);
543: PETSC_INTERN PetscInt PetscGCD(PetscInt a, PetscInt b);
545: #endif