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 (*createdefaultconstraints)(DM);
 45:   PetscErrorCode (*createglobalvector)(DM, Vec *);
 46:   PetscErrorCode (*createlocalvector)(DM, Vec *);
 47:   PetscErrorCode (*getlocaltoglobalmapping)(DM);
 48:   PetscErrorCode (*createfieldis)(DM, PetscInt *, char ***, IS **);
 49:   PetscErrorCode (*createcoordinatedm)(DM, DM *);
 50:   PetscErrorCode (*createcoordinatefield)(DM, DMField *);

 52:   PetscErrorCode (*getcoloring)(DM, ISColoringType, ISColoring *);
 53:   PetscErrorCode (*creatematrix)(DM, Mat *);
 54:   PetscErrorCode (*createinterpolation)(DM, DM, Mat *, Vec *);
 55:   PetscErrorCode (*createrestriction)(DM, DM, Mat *);
 56:   PetscErrorCode (*createmassmatrix)(DM, DM, Mat *);
 57:   PetscErrorCode (*createmassmatrixlumped)(DM, Vec *);
 58:   PetscErrorCode (*hascreateinjection)(DM, PetscBool *);
 59:   PetscErrorCode (*createinjection)(DM, DM, Mat *);

 61:   PetscErrorCode (*refine)(DM, MPI_Comm, DM *);
 62:   PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *);
 63:   PetscErrorCode (*refinehierarchy)(DM, PetscInt, DM *);
 64:   PetscErrorCode (*coarsenhierarchy)(DM, PetscInt, DM *);
 65:   PetscErrorCode (*extrude)(DM, PetscInt, DM *);

 67:   PetscErrorCode (*globaltolocalbegin)(DM, Vec, InsertMode, Vec);
 68:   PetscErrorCode (*globaltolocalend)(DM, Vec, InsertMode, Vec);
 69:   PetscErrorCode (*localtoglobalbegin)(DM, Vec, InsertMode, Vec);
 70:   PetscErrorCode (*localtoglobalend)(DM, Vec, InsertMode, Vec);
 71:   PetscErrorCode (*localtolocalbegin)(DM, Vec, InsertMode, Vec);
 72:   PetscErrorCode (*localtolocalend)(DM, Vec, InsertMode, Vec);

 74:   PetscErrorCode (*destroy)(DM);

 76:   PetscErrorCode (*computevariablebounds)(DM, Vec, Vec);

 78:   PetscErrorCode (*createsubdm)(DM, PetscInt, const PetscInt *, IS *, DM *);
 79:   PetscErrorCode (*createsuperdm)(DM *, PetscInt, IS **, DM *);
 80:   PetscErrorCode (*createfielddecomposition)(DM, PetscInt *, char ***, IS **, DM **);
 81:   PetscErrorCode (*createdomaindecomposition)(DM, PetscInt *, char ***, IS **, IS **, DM **);
 82:   PetscErrorCode (*createddscatters)(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);

 84:   PetscErrorCode (*getdimpoints)(DM, PetscInt, PetscInt *, PetscInt *);
 85:   PetscErrorCode (*locatepoints)(DM, Vec, DMPointLocationType, PetscSF);
 86:   PetscErrorCode (*getneighbors)(DM, PetscInt *, const PetscMPIInt **);
 87:   PetscErrorCode (*getboundingbox)(DM, PetscReal *, PetscReal *);
 88:   PetscErrorCode (*getlocalboundingbox)(DM, PetscReal *, PetscReal *);
 89:   PetscErrorCode (*locatepointssubdomain)(DM, Vec, PetscMPIInt **);

 91:   PetscErrorCode (*projectfunctionlocal)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
 92:   PetscErrorCode (*projectfunctionlabellocal)(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
 93:   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);
 94:   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);
 95:   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);
 96:   PetscErrorCode (*computel2diff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
 97:   PetscErrorCode (*computel2gradientdiff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
 98:   PetscErrorCode (*computel2fielddiff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);

100:   PetscErrorCode (*getcompatibility)(DM, DM, PetscBool *, PetscBool *);
101: };

103: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
104: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
105: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);

107: PETSC_INTERN PetscErrorCode DMCountNonCyclicReferences(PetscObject, PetscInt *);

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:   PetscDS dsIn;   /* Approximation space for input to this domain (now only used for cohesive cells) */
195:   DMLabel label;  /* Label defining the domain of definition of the discretization */
196:   IS      fields; /* Map from DS field numbers to original field numbers in the DM */
197: } DMSpace;

199: struct _p_UniversalLabel {
200:   DMLabel   label;   /* The universal label */
201:   PetscInt  Nl;      /* Number of labels encoded */
202:   char    **names;   /* The label names */
203:   PetscInt *indices; /* The original indices in the input DM */
204:   PetscInt  Nv;      /* Total number of values in all the labels */
205:   PetscInt *bits;    /* Starting bit for values of each label */
206:   PetscInt *masks;   /* Masks to pull out label value bits for each label */
207:   PetscInt *offsets; /* Starting offset for label values for each label */
208:   PetscInt *values;  /* Original label values before renumbering */
209: };

211: PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);

213: #define MAXDMMONITORS 5

215: typedef struct {
216:   PetscInt dim;   /* The dimension of the embedding space */
217:   DM       dm;    /* Layout for coordinates */
218:   Vec      x;     /* Coordinate values in global vector */
219:   Vec      xl;    /* Coordinate values in local  vector */
220:   DMField  field; /* Coordinates as an abstract field */
221: } DMCoordinates;

223: struct _p_DM {
224:   PETSCHEADER(struct _DMOps);
225:   Vec            localin[DM_MAX_WORK_VECTORS], localout[DM_MAX_WORK_VECTORS];
226:   Vec            globalin[DM_MAX_WORK_VECTORS], globalout[DM_MAX_WORK_VECTORS];
227:   DMNamedVecLink namedglobal;
228:   DMNamedVecLink namedlocal;
229:   DMWorkLink     workin, workout;
230:   DMLabelLink    labels;        /* Linked list of labels */
231:   DMLabel        depthLabel;    /* Optimized access to depth label */
232:   DMLabel        celltypeLabel; /* Optimized access to celltype label */
233:   void          *ctx;           /* a user context */
234:   PetscErrorCode (*ctxdestroy)(void **);
235:   ISColoringType         coloringtype;
236:   MatFDColoring          fd;
237:   VecType                vectype;    /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
238:   MatType                mattype;    /* type of matrix created with DMCreateMatrix() */
239:   PetscInt               bind_below; /* Local size threshold (in entries/rows) below which Vec/Mat objects are bound to CPU */
240:   PetscInt               bs;
241:   DMBlockingType         blocking_type;
242:   ISLocalToGlobalMapping ltogmap;
243:   PetscBool              prealloc_skip;      // Flag indicating the DMCreateMatrix() should not preallocate (only set sizes and local-to-global)
244:   PetscBool              prealloc_only;      /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
245:   PetscBool              structure_only;     /* Flag indicating the DMCreateMatrix() create matrix structure without values */
246:   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 */
247:   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 */
248:   PetscBool              setfromoptionscalled;
249:   void                  *data;
250:   /* Hierarchy / Submeshes */
251:   DM                      coarseMesh;
252:   DM                      fineMesh;
253:   DMCoarsenHookLink       coarsenhook; /* For transferring auxiliary problem data to coarser grids */
254:   DMRefineHookLink        refinehook;
255:   DMSubDomainHookLink     subdomainhook;
256:   DMGlobalToLocalHookLink gtolhook;
257:   DMLocalToGlobalHookLink ltoghook;
258:   /* Topology */
259:   PetscInt dim; /* The topological dimension */
260:   /* Auxiliary data */
261:   PetscHMapAux auxData; /* Auxiliary DM and Vec for region denoted by the key */
262:   /* Flexible communication */
263:   PetscSF   sfMigration; /* SF for point distribution created during distribution */
264:   PetscSF   sf;          /* SF for parallel point overlap */
265:   PetscSF   sectionSF;   /* SF for parallel dof overlap using section */
266:   PetscSF   sfNatural;   /* SF mapping to the "natural" ordering */
267:   PetscBool useNatural;  /* Create the natural SF */
268:   /* Allows a non-standard data layout */
269:   PetscBool    adjacency[2];  /* [use cone() or support() first, use the transitive closure] */
270:   PetscSection localSection;  /* Layout for local vectors */
271:   PetscSection globalSection; /* Layout for global vectors */
272:   PetscLayout  map;
273:   // Affine transform applied in DMGlobalToLocal
274:   struct {
275:     VecScatter affine_to_local;
276:     Vec        affine;
277:     PetscErrorCode (*setup)(DM);
278:   } periodic;
279:   /* Constraints */
280:   struct {
281:     PetscSection section;
282:     Mat          mat;
283:     Vec          bias;
284:   } defaultConstraint;
285:   /* Basis transformation */
286:   DM    transformDM;  /* Layout for basis transformation */
287:   Vec   transform;    /* Basis transformation matrices */
288:   void *transformCtx; /* Basis transformation context */
289:   PetscErrorCode (*transformSetUp)(DM, void *);
290:   PetscErrorCode (*transformDestroy)(DM, void *);
291:   PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
292:   /* Coordinates */
293:   DMCoordinates coordinates[2]; /* Coordinates, 0 is default, 1 is possible DG coordinate field for periodicity */
294:   /* Periodicity */
295:   PetscReal *Lstart, *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
296:   PetscBool  sparseLocalize;       /* Localize coordinates only for cells near periodic boundary */
297:   /* Null spaces -- of course I should make this have a variable number of fields */
298:   NullSpaceFunc nullspaceConstructors[10];
299:   NullSpaceFunc nearnullspaceConstructors[10];
300:   /* Fields are represented by objects */
301:   PetscInt     Nf;       /* Number of fields defined on the total domain */
302:   RegionField *fields;   /* Array of discretization fields with regions of validity */
303:   DMBoundary   boundary; /* List of boundary conditions */
304:   PetscInt     Nds;      /* Number of discrete systems defined on the total domain */
305:   DMSpace     *probs;    /* Array of discrete systems */
306:   /* Output structures */
307:   DM        dmBC;              /* The DM with boundary conditions in the global DM */
308:   PetscInt  outputSequenceNum; /* The current sequence number for output */
309:   PetscReal outputSequenceVal; /* The current sequence value for output */
310:   PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
311:   PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
312:   void    *monitorcontext[MAXDMMONITORS];
313:   PetscInt numbermonitors;
314:   /* Configuration */
315:   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 */

317:   PetscObject dmksp, dmsnes, dmts;
318: #ifdef PETSC_HAVE_LIBCEED
319:   Ceed                ceed;          // LibCEED context
320:   CeedElemRestriction ceedERestrict; // Map from the local vector (Lvector) to the cells (Evector)
321: #endif
322:   DMCeed dmceed; // CEED operator and data for this problem
323: };

325: PETSC_EXTERN PetscLogEvent DM_Convert;
326: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
327: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
328: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
329: PETSC_EXTERN PetscLogEvent DM_Coarsen;
330: PETSC_EXTERN PetscLogEvent DM_Refine;
331: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
332: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
333: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
334: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
335: PETSC_EXTERN PetscLogEvent DM_CreateMassMatrix;
336: PETSC_EXTERN PetscLogEvent DM_Load;
337: PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
338: PETSC_EXTERN PetscLogEvent DM_ProjectFunction;

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);