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