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

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

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

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

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

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

 35: typedef struct _DMOps *DMOps;
 36: struct _DMOps {
 37:   PetscErrorCode (*view)(DM, PetscViewer);
 38:   PetscErrorCode (*load)(DM, PetscViewer);
 39:   PetscErrorCode (*clone)(DM, DM *);
 40:   PetscErrorCode (*setfromoptions)(DM, PetscOptionItems *);
 41:   PetscErrorCode (*setup)(DM);
 42:   PetscErrorCode (*createlocalsection)(DM);
 43:   PetscErrorCode (*createsectionpermutation)(DM, IS *, PetscBT *);
 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 *, 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[], PetscInt[], PetscInt[]);
 89:   PetscErrorCode (*locatepointssubdomain)(DM, Vec, PetscMPIInt **);
 90:   PetscErrorCode (*snaptogeommodel)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);

 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: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*NullSpaceFn)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace);

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

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

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

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

349: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM, Vec *);
350: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM, Vec *);

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

354: /*

356:           Composite Vectors

358:       Single global representation
359:       Individual global representations
360:       Single local representation
361:       Individual local representations

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

365:        DM da_u, da_v, da_p

367:        DM dm_velocities

369:        DM dm

371:        DMDACreate(,&da_u);
372:        DMDACreate(,&da_v);
373:        DMCompositeCreate(,&dm_velocities);
374:        DMCompositeAddDM(dm_velocities,(DM)du);
375:        DMCompositeAddDM(dm_velocities,(DM)dv);

377:        DMDACreate(,&da_p);
378:        DMCompositeCreate(,&dm_velocities);
379:        DMCompositeAddDM(dm,(DM)dm_velocities);
380:        DMCompositeAddDM(dm,(DM)dm_p);

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

387:     Element access
388:     --------------
389:       From global vectors
390:          -DAVecGetArray   - for DMDA
391:          -VecGetArray - for DMSliced
392:          ADD for DMComposite???  maybe

394:       From individual vector
395:           -DAVecGetArray - for DMDA
396:           -VecGetArray -for sliced
397:          ADD for DMComposite??? maybe

399:       From single local vector
400:           ADD         * single local vector as arrays?

402:    Communication
403:    -------------
404:       DMGlobalToLocal - global vector to single local vector

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

408:    Obtaining vectors
409:    -----------------
410:       DMCreateGlobal/Local
411:       DMGetGlobal/Local
412:       DMCompositeGetLocalVectors   - gives individual local work vectors and arrays

414:     ?????   individual global vectors   ????

416: */

418: #if defined(PETSC_HAVE_HDF5)
419: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char[], PetscInt, PetscScalar *, PetscViewer);
420: PETSC_EXTERN PetscErrorCode DMSequenceGetLength_HDF5_Internal(DM, const char[], PetscInt *, PetscViewer);
421: #endif

423: static inline PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
424: {
425:   PetscFunctionBeginHot;
426: #if defined(PETSC_USE_DEBUG)
427:   {
428:     PetscInt dof;

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

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

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

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

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

538: PETSC_EXTERN PetscErrorCode DMGetCoordinateDegree_Internal(DM, PetscInt *);
539: PETSC_INTERN PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM, PetscReal[], PetscReal[], PetscInt[], PetscInt[]);
540: PETSC_INTERN PetscErrorCode DMGetIsoperiodicPointSF_Internal(DM dm, PetscSF *sf);

542: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
543: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
544: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);

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

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

552: PETSC_INTERN PetscErrorCode DMCompleteBCLabels_Internal(DM dm);
553: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
554: PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
555: PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
556: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
557: PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);
558: PETSC_INTERN PetscInt       PetscGCD(PetscInt a, PetscInt b);