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: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
110: struct _DMCoarsenHookLink {
111: PetscErrorCode (*coarsenhook)(DM, DM, void *); /* Run once, when coarse DM is created */
112: PetscErrorCode (*restricthook)(DM, Mat, Vec, Mat, DM, void *); /* Run each time a new problem is restricted to a coarse grid */
113: void *ctx;
114: DMCoarsenHookLink next;
115: };
117: typedef struct _DMRefineHookLink *DMRefineHookLink;
118: struct _DMRefineHookLink {
119: PetscErrorCode (*refinehook)(DM, DM, void *); /* Run once, when a fine DM is created */
120: PetscErrorCode (*interphook)(DM, Mat, DM, void *); /* Run each time a new problem is interpolated to a fine grid */
121: void *ctx;
122: DMRefineHookLink next;
123: };
125: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
126: struct _DMSubDomainHookLink {
127: PetscErrorCode (*ddhook)(DM, DM, void *);
128: PetscErrorCode (*restricthook)(DM, VecScatter, VecScatter, DM, void *);
129: void *ctx;
130: DMSubDomainHookLink next;
131: };
133: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
134: struct _DMGlobalToLocalHookLink {
135: PetscErrorCode (*beginhook)(DM, Vec, InsertMode, Vec, void *);
136: PetscErrorCode (*endhook)(DM, Vec, InsertMode, Vec, void *);
137: void *ctx;
138: DMGlobalToLocalHookLink next;
139: };
141: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
142: struct _DMLocalToGlobalHookLink {
143: PetscErrorCode (*beginhook)(DM, Vec, InsertMode, Vec, void *);
144: PetscErrorCode (*endhook)(DM, Vec, InsertMode, Vec, void *);
145: void *ctx;
146: DMLocalToGlobalHookLink next;
147: };
149: typedef enum {
150: DMVEC_STATUS_IN,
151: DMVEC_STATUS_OUT
152: } DMVecStatus;
153: typedef struct _DMNamedVecLink *DMNamedVecLink;
154: struct _DMNamedVecLink {
155: Vec X;
156: char *name;
157: DMVecStatus status;
158: DMNamedVecLink next;
159: };
161: typedef struct _DMWorkLink *DMWorkLink;
162: struct _DMWorkLink {
163: size_t bytes;
164: void *mem;
165: DMWorkLink next;
166: };
168: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
170: struct _n_DMLabelLink {
171: DMLabel label;
172: PetscBool output;
173: struct _n_DMLabelLink *next;
174: };
175: typedef struct _n_DMLabelLink *DMLabelLink;
177: typedef struct _n_Boundary *DMBoundary;
179: struct _n_Boundary {
180: DSBoundary dsboundary;
181: DMLabel label;
182: DMBoundary next;
183: };
185: typedef struct _n_Field {
186: PetscObject disc; /* Field discretization, or a PetscContainer with the field name */
187: DMLabel label; /* Label defining the domain of definition of the field */
188: PetscBool adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
189: PetscBool avoidTensor; /* Flag to avoid defining field over tensor cells */
190: } RegionField;
192: typedef struct _n_Space {
193: PetscDS ds; /* Approximation space in this domain */
194: DMLabel label; /* Label defining the domain of definition of the discretization */
195: IS fields; /* Map from DS field numbers to original field numbers in the DM */
196: } DMSpace;
198: struct _p_UniversalLabel {
199: DMLabel label; /* The universal label */
200: PetscInt Nl; /* Number of labels encoded */
201: char **names; /* The label names */
202: PetscInt *indices; /* The original indices in the input DM */
203: PetscInt Nv; /* Total number of values in all the labels */
204: PetscInt *bits; /* Starting bit for values of each label */
205: PetscInt *masks; /* Masks to pull out label value bits for each label */
206: PetscInt *offsets; /* Starting offset for label values for each label */
207: PetscInt *values; /* Original label values before renumbering */
208: };
210: PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);
212: #define MAXDMMONITORS 5
214: typedef struct {
215: PetscInt dim; /* The dimension of the embedding space */
216: DM dm; /* Layout for coordinates */
217: Vec x; /* Coordinate values in global vector */
218: Vec xl; /* Coordinate values in local vector */
219: DMField field; /* Coordinates as an abstract field */
220: } DMCoordinates;
222: struct _p_DM {
223: PETSCHEADER(struct _DMOps);
224: Vec localin[DM_MAX_WORK_VECTORS], localout[DM_MAX_WORK_VECTORS];
225: Vec globalin[DM_MAX_WORK_VECTORS], globalout[DM_MAX_WORK_VECTORS];
226: DMNamedVecLink namedglobal;
227: DMNamedVecLink namedlocal;
228: DMWorkLink workin, workout;
229: DMLabelLink labels; /* Linked list of labels */
230: DMLabel depthLabel; /* Optimized access to depth label */
231: DMLabel celltypeLabel; /* Optimized access to celltype label */
232: void *ctx; /* a user context */
233: PetscErrorCode (*ctxdestroy)(void **);
234: ISColoringType coloringtype;
235: MatFDColoring fd;
236: VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
237: MatType mattype; /* type of matrix created with DMCreateMatrix() */
238: PetscInt bind_below; /* Local size threshold (in entries/rows) below which Vec/Mat objects are bound to CPU */
239: PetscInt bs;
240: ISLocalToGlobalMapping ltogmap;
241: PetscBool prealloc_skip; // Flag indicating the DMCreateMatrix() should not preallocate (only set sizes and local-to-global)
242: PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
243: PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
244: 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 */
245: 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 */
246: PetscBool setfromoptionscalled;
247: void *data;
248: /* Hierarchy / Submeshes */
249: DM coarseMesh;
250: DM fineMesh;
251: DMCoarsenHookLink coarsenhook; /* For transfering auxiliary problem data to coarser grids */
252: DMRefineHookLink refinehook;
253: DMSubDomainHookLink subdomainhook;
254: DMGlobalToLocalHookLink gtolhook;
255: DMLocalToGlobalHookLink ltoghook;
256: /* Topology */
257: PetscInt dim; /* The topological dimension */
258: /* Auxiliary data */
259: PetscHMapAux auxData; /* Auxiliary DM and Vec for region denoted by the key */
260: /* Flexible communication */
261: PetscSF sfMigration; /* SF for point distribution created during distribution */
262: PetscSF sf; /* SF for parallel point overlap */
263: PetscSF sectionSF; /* SF for parallel dof overlap using section */
264: PetscSF sfNatural; /* SF mapping to the "natural" ordering */
265: PetscBool useNatural; /* Create the natural SF */
266: /* Allows a non-standard data layout */
267: PetscBool adjacency[2]; /* [use cone() or support() first, use the transitive closure] */
268: PetscSection localSection; /* Layout for local vectors */
269: PetscSection globalSection; /* Layout for global vectors */
270: PetscLayout map;
271: /* Constraints */
272: struct {
273: PetscSection section;
274: Mat mat;
275: Vec bias;
276: } defaultConstraint;
277: /* Basis transformation */
278: DM transformDM; /* Layout for basis transformation */
279: Vec transform; /* Basis transformation matrices */
280: void *transformCtx; /* Basis transformation context */
281: PetscErrorCode (*transformSetUp)(DM, void *);
282: PetscErrorCode (*transformDestroy)(DM, void *);
283: PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
284: /* Coordinates */
285: DMCoordinates coordinates[2]; /* Coordinates, 0 is default, 1 is possible DG coordinate field for periodicity */
286: /* Periodicity */
287: PetscReal *Lstart, *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
288: PetscBool sparseLocalize; /* Localize coordinates only for cells near periodic boundary */
289: /* Null spaces -- of course I should make this have a variable number of fields */
290: NullSpaceFunc nullspaceConstructors[10];
291: NullSpaceFunc nearnullspaceConstructors[10];
292: /* Fields are represented by objects */
293: PetscInt Nf; /* Number of fields defined on the total domain */
294: RegionField *fields; /* Array of discretization fields with regions of validity */
295: DMBoundary boundary; /* List of boundary conditions */
296: PetscInt Nds; /* Number of discrete systems defined on the total domain */
297: DMSpace *probs; /* Array of discrete systems */
298: /* Output structures */
299: DM dmBC; /* The DM with boundary conditions in the global DM */
300: PetscInt outputSequenceNum; /* The current sequence number for output */
301: PetscReal outputSequenceVal; /* The current sequence value for output */
302: PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
303: PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
304: void *monitorcontext[MAXDMMONITORS];
305: PetscInt numbermonitors;
307: PetscObject dmksp, dmsnes, dmts;
308: #ifdef PETSC_HAVE_LIBCEED
309: Ceed ceed; /* LibCEED context */
310: CeedElemRestriction ceedERestrict; /* Map from the local vector (Lvector) to the cells (Evector) */
311: #endif
312: };
314: PETSC_EXTERN PetscLogEvent DM_Convert;
315: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
316: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
317: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
318: PETSC_EXTERN PetscLogEvent DM_Coarsen;
319: PETSC_EXTERN PetscLogEvent DM_Refine;
320: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
321: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
322: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
323: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
324: PETSC_EXTERN PetscLogEvent DM_CreateMassMatrix;
325: PETSC_EXTERN PetscLogEvent DM_Load;
326: PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
328: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM, Vec *);
329: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM, Vec *);
331: PETSC_EXTERN PetscErrorCode DMView_GLVis(DM, PetscViewer, PetscErrorCode (*)(DM, PetscViewer));
333: /*
335: Composite Vectors
337: Single global representation
338: Individual global representations
339: Single local representation
340: Individual local representations
342: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
344: DM da_u, da_v, da_p
346: DM dm_velocities
348: DM dm
350: DMDACreate(,&da_u);
351: DMDACreate(,&da_v);
352: DMCompositeCreate(,&dm_velocities);
353: DMCompositeAddDM(dm_velocities,(DM)du);
354: DMCompositeAddDM(dm_velocities,(DM)dv);
356: DMDACreate(,&da_p);
357: DMCompositeCreate(,&dm_velocities);
358: DMCompositeAddDM(dm,(DM)dm_velocities);
359: DMCompositeAddDM(dm,(DM)dm_p);
361: Access parts of composite vectors (Composite only)
362: ---------------------------------
363: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
364: ADD for local vector -
366: Element access
367: --------------
368: From global vectors
369: -DAVecGetArray - for DMDA
370: -VecGetArray - for DMSliced
371: ADD for DMComposite??? maybe
373: From individual vector
374: -DAVecGetArray - for DMDA
375: -VecGetArray -for sliced
376: ADD for DMComposite??? maybe
378: From single local vector
379: ADD * single local vector as arrays?
381: Communication
382: -------------
383: DMGlobalToLocal - global vector to single local vector
385: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
387: Obtaining vectors
388: -----------------
389: DMCreateGlobal/Local
390: DMGetGlobal/Local
391: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
393: ????? individual global vectors ????
395: */
397: #if defined(PETSC_HAVE_HDF5)
398: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
399: #endif
401: static inline PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
402: {
404: #if defined(PETSC_USE_DEBUG)
405: {
406: PetscInt dof;
408: *start = *end = 0; /* Silence overzealous compiler warning */
410: PetscSectionGetOffset(dm->localSection, point, start);
411: PetscSectionGetDof(dm->localSection, point, &dof);
412: *end = *start + dof;
413: }
414: #else
415: {
416: const PetscSection s = dm->localSection;
417: *start = s->atlasOff[point - s->pStart];
418: *end = *start + s->atlasDof[point - s->pStart];
419: }
420: #endif
421: return 0;
422: }
424: static inline PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
425: {
426: #if defined(PETSC_USE_DEBUG)
427: {
428: PetscInt dof;
429: *start = *end = 0; /* Silence overzealous compiler warning */
431: PetscSectionGetFieldOffset(dm->localSection, point, field, start);
432: PetscSectionGetFieldDof(dm->localSection, point, field, &dof);
433: *end = *start + dof;
434: }
435: #else
436: {
437: const PetscSection s = dm->localSection->field[field];
438: *start = s->atlasOff[point - s->pStart];
439: *end = *start + s->atlasDof[point - s->pStart];
440: }
441: #endif
442: return 0;
443: }
445: static inline PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
446: {
447: #if defined(PETSC_USE_DEBUG)
448: {
449: PetscInt dof, cdof;
450: *start = *end = 0; /* Silence overzealous compiler warning */
453: PetscSectionGetOffset(dm->globalSection, point, start);
454: PetscSectionGetDof(dm->globalSection, point, &dof);
455: PetscSectionGetConstraintDof(dm->globalSection, point, &cdof);
456: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
457: }
458: #else
459: {
460: const PetscSection s = dm->globalSection;
461: const PetscInt dof = s->atlasDof[point - s->pStart];
462: const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
463: *start = s->atlasOff[point - s->pStart];
464: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
465: }
466: #endif
467: return 0;
468: }
470: static inline PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
471: {
472: #if defined(PETSC_USE_DEBUG)
473: {
474: PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
475: *start = *end = 0; /* Silence overzealous compiler warning */
478: PetscSectionGetOffset(dm->globalSection, point, start);
479: PetscSectionGetOffset(dm->localSection, point, &loff);
480: PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff);
481: PetscSectionGetFieldDof(dm->localSection, point, field, &fdof);
482: PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof);
483: *start = *start < 0 ? *start - (lfoff - loff) : *start + lfoff - loff;
484: for (f = 0; f < field; ++f) {
485: PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof);
486: *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
487: }
488: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
489: }
490: #else
491: {
492: const PetscSection s = dm->localSection;
493: const PetscSection fs = dm->localSection->field[field];
494: const PetscSection gs = dm->globalSection;
495: const PetscInt loff = s->atlasOff[point - s->pStart];
496: const PetscInt goff = gs->atlasOff[point - s->pStart];
497: const PetscInt lfoff = fs->atlasOff[point - s->pStart];
498: const PetscInt fdof = fs->atlasDof[point - s->pStart];
499: const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
500: PetscInt ffcdof = 0, f;
502: for (f = 0; f < field; ++f) {
503: const PetscSection ffs = dm->localSection->field[f];
504: ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
505: }
506: *start = goff + (goff < 0 ? loff - lfoff + ffcdof : lfoff - loff - ffcdof);
507: *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
508: }
509: #endif
510: return 0;
511: }
513: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
514: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
515: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);
517: PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
518: PETSC_INTERN PetscErrorCode DMSetField_Internal(DM, PetscInt, DMLabel, PetscObject);
520: PETSC_INTERN PetscErrorCode DMSetLabelValue_Fast(DM, DMLabel *, const char[], PetscInt, PetscInt);
522: PETSC_INTERN PetscErrorCode DMCompleteBCLabels_Internal(DM dm);
523: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
524: PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
525: PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
526: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
527: PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);
528: PETSC_INTERN PetscInt PetscGCD(PetscInt a, PetscInt b);
530: #endif