Actual source code: dmimpl.h

  1: #pragma once

  3: #include <petscdm.h>
  4: #ifdef PETSC_HAVE_LIBCEED
  5: #include <petscdmceed.h>
  6: #endif
  7: #include <petsc/private/petscimpl.h>
  8: #include <petsc/private/petscdsimpl.h>
  9: #include <petsc/private/sectionimpl.h>

 11: PETSC_EXTERN PetscBool      DMRegisterAllCalled;
 12: PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
 13: typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace);

 15: typedef struct _PetscHashAuxKey {
 16:   DMLabel  label;
 17:   PetscInt value;
 18:   PetscInt part;
 19: } PetscHashAuxKey;

 21: #define PetscHashAuxKeyHash(key) PetscHashCombine(PetscHashCombine(PetscHashPointer((key).label), PetscHashInt((key).value)), PetscHashInt((key).part))

 23: #define PetscHashAuxKeyEqual(k1, k2) (((k1).label == (k2).label) ? (((k1).value == (k2).value) ? ((k1).part == (k2).part) : 0) : 0)

 25: PETSC_HASH_MAP(HMapAux, PetscHashAuxKey, Vec, PetscHashAuxKeyHash, PetscHashAuxKeyEqual, PETSC_NULLPTR)

 27: struct _n_DMGeneratorFunctionList {
 28:   PetscErrorCode (*generate)(DM, PetscBool, DM *);
 29:   PetscErrorCode (*refine)(DM, PetscReal *, DM *);
 30:   PetscErrorCode (*adapt)(DM, Vec, DMLabel, DMLabel, DM *);
 31:   char                   *name;
 32:   PetscInt                dim;
 33:   DMGeneratorFunctionList next;
 34: };

 36: typedef struct _DMOps *DMOps;
 37: struct _DMOps {
 38:   PetscErrorCode (*view)(DM, PetscViewer);
 39:   PetscErrorCode (*load)(DM, PetscViewer);
 40:   PetscErrorCode (*clone)(DM, DM *);
 41:   PetscErrorCode (*setfromoptions)(DM, PetscOptionItems *);
 42:   PetscErrorCode (*setup)(DM);
 43:   PetscErrorCode (*createlocalsection)(DM);
 44:   PetscErrorCode (*createsectionpermutation)(DM, IS *, PetscBT *);
 45:   PetscErrorCode (*createdefaultconstraints)(DM);
 46:   PetscErrorCode (*createglobalvector)(DM, Vec *);
 47:   PetscErrorCode (*createlocalvector)(DM, Vec *);
 48:   PetscErrorCode (*getlocaltoglobalmapping)(DM);
 49:   PetscErrorCode (*createfieldis)(DM, PetscInt *, char ***, IS **);
 50:   PetscErrorCode (*createcoordinatedm)(DM, DM *);
 51:   PetscErrorCode (*createcoordinatefield)(DM, DMField *);

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

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

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

 75:   PetscErrorCode (*destroy)(DM);

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

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

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

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

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

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

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

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 {
151:   DMVEC_STATUS_IN,
152:   DMVEC_STATUS_OUT
153: } DMVecStatus;
154: typedef struct _DMNamedVecLink *DMNamedVecLink;
155: struct _DMNamedVecLink {
156:   Vec            X;
157:   char          *name;
158:   DMVecStatus    status;
159:   DMNamedVecLink next;
160: };

162: typedef struct _DMWorkLink *DMWorkLink;
163: struct _DMWorkLink {
164:   size_t     bytes;
165:   void      *mem;
166:   DMWorkLink next;
167: };

169: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users  via DMGetGlobalVector(), DMGetLocalVector() */

171: struct _n_DMLabelLink {
172:   DMLabel                label;
173:   PetscBool              output;
174:   struct _n_DMLabelLink *next;
175: };
176: typedef struct _n_DMLabelLink *DMLabelLink;

178: typedef struct _n_Boundary *DMBoundary;

180: struct _n_Boundary {
181:   DSBoundary dsboundary;
182:   DMLabel    label;
183:   DMBoundary next;
184: };

186: typedef struct _n_Field {
187:   PetscObject disc;         /* Field discretization, or a PetscContainer with the field name */
188:   DMLabel     label;        /* Label defining the domain of definition of the field */
189:   PetscBool   adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
190:   PetscBool   avoidTensor;  /* Flag to avoid defining field over tensor cells */
191: } RegionField;

193: typedef struct _n_Space {
194:   PetscDS ds;     /* Approximation space in this domain */
195:   PetscDS dsIn;   /* Approximation space for input to this domain (now only used for cohesive cells) */
196:   DMLabel label;  /* Label defining the domain of definition of the discretization */
197:   IS      fields; /* Map from DS field numbers to original field numbers in the DM */
198: } DMSpace;

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

212: PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);

214: #define MAXDMMONITORS 5

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

224: struct _p_DM {
225:   PETSCHEADER(struct _DMOps);
226:   Vec            localin[DM_MAX_WORK_VECTORS], localout[DM_MAX_WORK_VECTORS];
227:   Vec            globalin[DM_MAX_WORK_VECTORS], globalout[DM_MAX_WORK_VECTORS];
228:   DMNamedVecLink namedglobal;
229:   DMNamedVecLink namedlocal;
230:   DMWorkLink     workin, workout;
231:   DMLabelLink    labels;        /* Linked list of labels */
232:   DMLabel        depthLabel;    /* Optimized access to depth label */
233:   DMLabel        celltypeLabel; /* Optimized access to celltype label */
234:   void          *ctx;           /* a user context */
235:   PetscErrorCode (*ctxdestroy)(void **);
236:   ISColoringType         coloringtype;
237:   MatFDColoring          fd;
238:   VecType                vectype;    /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
239:   MatType                mattype;    /* type of matrix created with DMCreateMatrix() */
240:   PetscInt               bind_below; /* Local size threshold (in entries/rows) below which Vec/Mat objects are bound to CPU */
241:   PetscInt               bs;
242:   DMBlockingType         blocking_type;
243:   ISLocalToGlobalMapping ltogmap;
244:   PetscBool              prealloc_skip;      // Flag indicating the DMCreateMatrix() should not preallocate (only set sizes and local-to-global)
245:   PetscBool              prealloc_only;      /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
246:   PetscBool              structure_only;     /* Flag indicating the DMCreateMatrix() create matrix structure without values */
247:   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 */
248:   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 */
249:   PetscBool              setfromoptionscalled;
250:   void                  *data;
251:   /* Hierarchy / Submeshes */
252:   DM                      coarseMesh;
253:   DM                      fineMesh;
254:   DMCoarsenHookLink       coarsenhook; /* For transferring auxiliary problem data to coarser grids */
255:   DMRefineHookLink        refinehook;
256:   DMSubDomainHookLink     subdomainhook;
257:   DMGlobalToLocalHookLink gtolhook;
258:   DMLocalToGlobalHookLink ltoghook;
259:   /* Topology */
260:   PetscInt dim; /* The topological dimension */
261:   /* Auxiliary data */
262:   PetscHMapAux auxData; /* Auxiliary DM and Vec for region denoted by the key */
263:   /* Flexible communication */
264:   PetscSF   sfMigration; /* SF for point distribution created during distribution */
265:   PetscSF   sf;          /* SF for parallel point overlap */
266:   PetscSF   sectionSF;   /* SF for parallel dof overlap using section */
267:   PetscSF   sfNatural;   /* SF mapping to the "natural" ordering */
268:   PetscBool useNatural;  /* Create the natural SF */
269:   /* Allows a non-standard data layout */
270:   PetscBool            adjacency[2];       /* [use cone() or support() first, use the transitive closure] */
271:   PetscSection         localSection;       /* Layout for local vectors */
272:   PetscSection         globalSection;      /* Layout for global vectors */
273:   PetscLayout          map;                /* Parallel division of unknowns across processes */
274:   DMReorderDefaultFlag reorderSection;     /* Reorder the local section by default */
275:   MatOrderingType      reorderSectionType; /* The type of reordering */

277:   // Affine transform applied in DMGlobalToLocal
278:   struct {
279:     PetscInt    num_affines;
280:     VecScatter *affine_to_local;
281:     Vec        *affine;
282:     PetscErrorCode (*setup)(DM);
283:   } periodic;
284:   /* Constraints */
285:   struct {
286:     PetscSection section;
287:     Mat          mat;
288:     Vec          bias;
289:   } defaultConstraint;
290:   /* Basis transformation */
291:   DM    transformDM;  /* Layout for basis transformation */
292:   Vec   transform;    /* Basis transformation matrices */
293:   void *transformCtx; /* Basis transformation context */
294:   PetscErrorCode (*transformSetUp)(DM, void *);
295:   PetscErrorCode (*transformDestroy)(DM, void *);
296:   PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
297:   /* Coordinates */
298:   DMCoordinates coordinates[2]; /* Coordinates, 0 is default, 1 is possible DG coordinate field for periodicity */
299:   /* Periodicity */
300:   PetscReal *Lstart, *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
301:   PetscBool  sparseLocalize;       /* Localize coordinates only for cells near periodic boundary */
302:   /* Null spaces -- of course I should make this have a variable number of fields */
303:   NullSpaceFunc nullspaceConstructors[10];
304:   NullSpaceFunc nearnullspaceConstructors[10];
305:   /* Fields are represented by objects */
306:   PetscInt     Nf;       /* Number of fields defined on the total domain */
307:   RegionField *fields;   /* Array of discretization fields with regions of validity */
308:   DMBoundary   boundary; /* List of boundary conditions */
309:   PetscInt     Nds;      /* Number of discrete systems defined on the total domain */
310:   DMSpace     *probs;    /* Array of discrete systems */
311:   /* Output structures */
312:   DM        dmBC;              /* The DM with boundary conditions in the global DM */
313:   PetscBool ignorePermOutput;  /* Ignore the local section permutation on output */
314:   PetscInt  outputSequenceNum; /* The current sequence number for output */
315:   PetscReal outputSequenceVal; /* The current sequence value for output */
316:   PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
317:   PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
318:   void    *monitorcontext[MAXDMMONITORS];
319:   PetscInt numbermonitors;
320:   /* Configuration */
321:   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 */

323:   PetscObject dmksp, dmsnes, dmts;
324: #ifdef PETSC_HAVE_LIBCEED
325:   Ceed                ceed;          // LibCEED context
326:   CeedElemRestriction ceedERestrict; // Map from the local vector (Lvector) to the cells (Evector)
327: #endif
328:   DMCeed dmceed; // CEED operator and data for this problem
329: };

331: PETSC_EXTERN PetscLogEvent DM_Convert;
332: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
333: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
334: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
335: PETSC_EXTERN PetscLogEvent DM_Coarsen;
336: PETSC_EXTERN PetscLogEvent DM_Refine;
337: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
338: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
339: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
340: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
341: PETSC_EXTERN PetscLogEvent DM_CreateMassMatrix;
342: PETSC_EXTERN PetscLogEvent DM_Load;
343: PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
344: PETSC_EXTERN PetscLogEvent DM_ProjectFunction;

346: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM, Vec *);
347: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM, Vec *);

349: PETSC_EXTERN PetscErrorCode DMView_GLVis(DM, PetscViewer, PetscErrorCode (*)(DM, PetscViewer));

351: /*

353:           Composite Vectors

355:       Single global representation
356:       Individual global representations
357:       Single local representation
358:       Individual local representations

360:       Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????

362:        DM da_u, da_v, da_p

364:        DM dm_velocities

366:        DM dm

368:        DMDACreate(,&da_u);
369:        DMDACreate(,&da_v);
370:        DMCompositeCreate(,&dm_velocities);
371:        DMCompositeAddDM(dm_velocities,(DM)du);
372:        DMCompositeAddDM(dm_velocities,(DM)dv);

374:        DMDACreate(,&da_p);
375:        DMCompositeCreate(,&dm_velocities);
376:        DMCompositeAddDM(dm,(DM)dm_velocities);
377:        DMCompositeAddDM(dm,(DM)dm_p);

379:     Access parts of composite vectors (Composite only)
380:     ---------------------------------
381:       DMCompositeGetAccess  - access the global vector as subvectors and array (for redundant arrays)
382:       ADD for local vector -

384:     Element access
385:     --------------
386:       From global vectors
387:          -DAVecGetArray   - for DMDA
388:          -VecGetArray - for DMSliced
389:          ADD for DMComposite???  maybe

391:       From individual vector
392:           -DAVecGetArray - for DMDA
393:           -VecGetArray -for sliced
394:          ADD for DMComposite??? maybe

396:       From single local vector
397:           ADD         * single local vector as arrays?

399:    Communication
400:    -------------
401:       DMGlobalToLocal - global vector to single local vector

403:       DMCompositeScatter/Gather - direct to individual local vectors and arrays   CHANGE name closer to GlobalToLocal?

405:    Obtaining vectors
406:    -----------------
407:       DMCreateGlobal/Local
408:       DMGetGlobal/Local
409:       DMCompositeGetLocalVectors   - gives individual local work vectors and arrays

411:     ?????   individual global vectors   ????

413: */

415: #if defined(PETSC_HAVE_HDF5)
416: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
417: #endif

419: static inline PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
420: {
421:   PetscFunctionBeginHot;
422: #if defined(PETSC_USE_DEBUG)
423:   {
424:     PetscInt dof;

426:     *start = *end = 0; /* Silence overzealous compiler warning */
427:     PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
428:     PetscCall(PetscSectionGetOffset(dm->localSection, point, start));
429:     PetscCall(PetscSectionGetDof(dm->localSection, point, &dof));
430:     *end = *start + dof;
431:   }
432: #else
433:   {
434:     const PetscSection s = dm->localSection;
435:     *start               = s->atlasOff[point - s->pStart];
436:     *end                 = *start + s->atlasDof[point - s->pStart];
437:   }
438: #endif
439:   PetscFunctionReturn(PETSC_SUCCESS);
440: }

442: static inline PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
443: {
444:   PetscFunctionBegin;
445: #if defined(PETSC_USE_DEBUG)
446:   {
447:     PetscInt dof;
448:     *start = *end = 0; /* Silence overzealous compiler warning */
449:     PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
450:     PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, start));
451:     PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &dof));
452:     *end = *start + dof;
453:   }
454: #else
455:   {
456:     const PetscSection s = dm->localSection->field[field];
457:     *start               = s->atlasOff[point - s->pStart];
458:     *end                 = *start + s->atlasDof[point - s->pStart];
459:   }
460: #endif
461:   PetscFunctionReturn(PETSC_SUCCESS);
462: }

464: static inline PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
465: {
466:   PetscFunctionBegin;
467: #if defined(PETSC_USE_DEBUG)
468:   {
469:     PetscInt dof, cdof;
470:     *start = *end = 0; /* Silence overzealous compiler warning */
471:     PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
472:     PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
473:     PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
474:     PetscCall(PetscSectionGetDof(dm->globalSection, point, &dof));
475:     PetscCall(PetscSectionGetConstraintDof(dm->globalSection, point, &cdof));
476:     *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
477:   }
478: #else
479:   {
480:     const PetscSection s    = dm->globalSection;
481:     const PetscInt     dof  = s->atlasDof[point - s->pStart];
482:     const PetscInt     cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
483:     *start                  = s->atlasOff[point - s->pStart];
484:     *end                    = *start + dof - cdof + (dof < 0 ? 1 : 0);
485:   }
486: #endif
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: static inline PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
491: {
492:   PetscFunctionBegin;
493: #if defined(PETSC_USE_DEBUG)
494:   {
495:     PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
496:     *start = *end = 0; /* Silence overzealous compiler warning */
497:     PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
498:     PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
499:     PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
500:     PetscCall(PetscSectionGetOffset(dm->localSection, point, &loff));
501:     PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff));
502:     PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &fdof));
503:     PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof));
504:     *start = *start < 0 ? *start - (lfoff - loff) : *start + lfoff - loff;
505:     for (f = 0; f < field; ++f) {
506:       PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof));
507:       *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
508:     }
509:     *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
510:   }
511: #else
512:   {
513:     const PetscSection s      = dm->localSection;
514:     const PetscSection fs     = dm->localSection->field[field];
515:     const PetscSection gs     = dm->globalSection;
516:     const PetscInt     loff   = s->atlasOff[point - s->pStart];
517:     const PetscInt     goff   = gs->atlasOff[point - s->pStart];
518:     const PetscInt     lfoff  = fs->atlasOff[point - s->pStart];
519:     const PetscInt     fdof   = fs->atlasDof[point - s->pStart];
520:     const PetscInt     fcdof  = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
521:     PetscInt           ffcdof = 0, f;

523:     for (f = 0; f < field; ++f) {
524:       const PetscSection ffs = dm->localSection->field[f];
525:       ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
526:     }
527:     *start = goff + (goff < 0 ? loff - lfoff + ffcdof : lfoff - loff - ffcdof);
528:     *end   = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
529:   }
530: #endif
531:   PetscFunctionReturn(PETSC_SUCCESS);
532: }

534: PETSC_EXTERN PetscErrorCode DMGetCoordinateDegree_Internal(DM, PetscInt *);
535: PETSC_INTERN PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM, PetscReal[], PetscReal[], PetscInt[], PetscInt[]);

537: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
538: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
539: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);

541: PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
542: PETSC_INTERN PetscErrorCode DMSetField_Internal(DM, PetscInt, DMLabel, PetscObject);

544: PETSC_INTERN PetscErrorCode DMSetLabelValue_Fast(DM, DMLabel *, const char[], PetscInt, PetscInt);
545: PETSC_INTERN PetscErrorCode DMGetPoints_Internal(DM, DMLabel, PetscInt, PetscInt, IS *);

547: PETSC_INTERN PetscErrorCode DMCompleteBCLabels_Internal(DM dm);
548: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
549: PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
550: PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
551: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
552: PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);
553: PETSC_INTERN PetscInt       PetscGCD(PetscInt a, PetscInt b);