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: } PetscHashAuxKey;

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

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

 27: PETSC_HASH_MAP(HMapAux, PetscHashAuxKey, Vec, PetscHashAuxKeyHash, PetscHashAuxKeyEqual, NULL)

 29: typedef struct _DMOps *DMOps;
 30: struct _DMOps {
 31:   PetscErrorCode (*view)(DM,PetscViewer);
 32:   PetscErrorCode (*load)(DM,PetscViewer);
 33:   PetscErrorCode (*clone)(DM,DM*);
 34:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,DM);
 35:   PetscErrorCode (*setup)(DM);
 36:   PetscErrorCode (*createlocalsection)(DM);
 37:   PetscErrorCode (*createdefaultconstraints)(DM);
 38:   PetscErrorCode (*createglobalvector)(DM,Vec*);
 39:   PetscErrorCode (*createlocalvector)(DM,Vec*);
 40:   PetscErrorCode (*getlocaltoglobalmapping)(DM);
 41:   PetscErrorCode (*createfieldis)(DM,PetscInt*,char***,IS**);
 42:   PetscErrorCode (*createcoordinatedm)(DM,DM*);
 43:   PetscErrorCode (*createcoordinatefield)(DM,DMField*);

 45:   PetscErrorCode (*getcoloring)(DM,ISColoringType,ISColoring*);
 46:   PetscErrorCode (*creatematrix)(DM, Mat*);
 47:   PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*);
 48:   PetscErrorCode (*createrestriction)(DM,DM,Mat*);
 49:   PetscErrorCode (*createmassmatrix)(DM,DM,Mat*);
 50:   PetscErrorCode (*hascreateinjection)(DM,PetscBool*);
 51:   PetscErrorCode (*createinjection)(DM,DM,Mat*);

 53:   PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
 54:   PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*);
 55:   PetscErrorCode (*refinehierarchy)(DM,PetscInt,DM*);
 56:   PetscErrorCode (*coarsenhierarchy)(DM,PetscInt,DM*);
 57:   PetscErrorCode (*adaptlabel)(DM,DMLabel,DM*);
 58:   PetscErrorCode (*adaptmetric)(DM,Vec,DMLabel,DM*);

 60:   PetscErrorCode (*globaltolocalbegin)(DM,Vec,InsertMode,Vec);
 61:   PetscErrorCode (*globaltolocalend)(DM,Vec,InsertMode,Vec);
 62:   PetscErrorCode (*localtoglobalbegin)(DM,Vec,InsertMode,Vec);
 63:   PetscErrorCode (*localtoglobalend)(DM,Vec,InsertMode,Vec);
 64:   PetscErrorCode (*localtolocalbegin)(DM,Vec,InsertMode,Vec);
 65:   PetscErrorCode (*localtolocalend)(DM,Vec,InsertMode,Vec);

 67:   PetscErrorCode (*destroy)(DM);

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

 71:   PetscErrorCode (*createsubdm)(DM,PetscInt,const PetscInt*,IS*,DM*);
 72:   PetscErrorCode (*createsuperdm)(DM*,PetscInt,IS**,DM*);
 73:   PetscErrorCode (*createfielddecomposition)(DM,PetscInt*,char***,IS**,DM**);
 74:   PetscErrorCode (*createdomaindecomposition)(DM,PetscInt*,char***,IS**,IS**,DM**);
 75:   PetscErrorCode (*createddscatters)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

 77:   PetscErrorCode (*getdimpoints)(DM,PetscInt,PetscInt*,PetscInt*);
 78:   PetscErrorCode (*locatepoints)(DM,Vec,DMPointLocationType,PetscSF);
 79:   PetscErrorCode (*getneighbors)(DM,PetscInt*,const PetscMPIInt**);
 80:   PetscErrorCode (*getboundingbox)(DM,PetscReal*,PetscReal*);
 81:   PetscErrorCode (*getlocalboundingbox)(DM,PetscReal*,PetscReal*);
 82:   PetscErrorCode (*locatepointssubdomain)(DM,Vec,PetscMPIInt**);

 84:   PetscErrorCode (*projectfunctionlocal)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
 85:   PetscErrorCode (*projectfunctionlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
 86:   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);
 87:   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);
 88:   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);
 89:   PetscErrorCode (*computel2diff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
 90:   PetscErrorCode (*computel2gradientdiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal [],const PetscReal[],PetscInt, PetscScalar *,void *),void **,Vec,const PetscReal[],PetscReal *);
 91:   PetscErrorCode (*computel2fielddiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);

 93:   PetscErrorCode (*getcompatibility)(DM,DM,PetscBool*,PetscBool*);
 94: };

 96: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
 97: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
 98: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);

100: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
101: struct _DMCoarsenHookLink {
102:   PetscErrorCode (*coarsenhook)(DM,DM,void*);              /* Run once, when coarse DM is created */
103:   PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*); /* Run each time a new problem is restricted to a coarse grid */
104:   void *ctx;
105:   DMCoarsenHookLink next;
106: };

108: typedef struct _DMRefineHookLink *DMRefineHookLink;
109: struct _DMRefineHookLink {
110:   PetscErrorCode (*refinehook)(DM,DM,void*);     /* Run once, when a fine DM is created */
111:   PetscErrorCode (*interphook)(DM,Mat,DM,void*); /* Run each time a new problem is interpolated to a fine grid */
112:   void *ctx;
113:   DMRefineHookLink next;
114: };

116: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
117: struct _DMSubDomainHookLink {
118:   PetscErrorCode (*ddhook)(DM,DM,void*);
119:   PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*);
120:   void *ctx;
121:   DMSubDomainHookLink next;
122: };

124: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
125: struct _DMGlobalToLocalHookLink {
126:   PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
127:   PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
128:   void *ctx;
129:   DMGlobalToLocalHookLink next;
130: };

132: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
133: struct _DMLocalToGlobalHookLink {
134:   PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
135:   PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
136:   void *ctx;
137:   DMLocalToGlobalHookLink next;
138: };

140: typedef enum {DMVEC_STATUS_IN,DMVEC_STATUS_OUT} DMVecStatus;
141: typedef struct _DMNamedVecLink *DMNamedVecLink;
142: struct _DMNamedVecLink {
143:   Vec X;
144:   char *name;
145:   DMVecStatus status;
146:   DMNamedVecLink next;
147: };

149: typedef struct _DMWorkLink *DMWorkLink;
150: struct _DMWorkLink {
151:   size_t     bytes;
152:   void       *mem;
153:   DMWorkLink next;
154: };

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

158: struct _n_DMLabelLink {
159:   DMLabel              label;
160:   PetscBool            output;
161:   struct _n_DMLabelLink *next;
162: };
163: typedef struct _n_DMLabelLink *DMLabelLink;

165: typedef struct _n_Boundary *DMBoundary;

167: struct _n_Boundary {
168:   DSBoundary  dsboundary;
169:   DMLabel     label;
170:   DMBoundary  next;
171: };

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

180: typedef struct _n_Space {
181:   PetscDS ds;     /* Approximation space in this domain */
182:   DMLabel label;  /* Label defining the domain of definition of the discretization */
183:   IS      fields; /* Map from DS field numbers to original field numbers in the DM */
184: } DMSpace;

186: struct _p_UniversalLabel {
187:   DMLabel    label;   /* The universal label */
188:   PetscInt   Nl;      /* Number of labels encoded */
189:   char     **names;   /* The label names */
190:   PetscInt  *indices; /* The original indices in the input DM */
191:   PetscInt   Nv;      /* Total number of values in all the labels */
192:   PetscInt  *bits;    /* Starting bit for values of each label */
193:   PetscInt  *masks;   /* Masks to pull out label value bits for each label */
194:   PetscInt  *offsets; /* Starting offset for label values for each label */
195:   PetscInt  *values;  /* Original label values before renumbering */
196: };

198: PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);

200: #define MAXDMMONITORS 5

202: struct _p_DM {
203:   PETSCHEADER(struct _DMOps);
204:   Vec                     localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
205:   Vec                     globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
206:   DMNamedVecLink          namedglobal;
207:   DMNamedVecLink          namedlocal;
208:   DMWorkLink              workin,workout;
209:   DMLabelLink             labels;            /* Linked list of labels */
210:   DMLabel                 depthLabel;        /* Optimized access to depth label */
211:   DMLabel                 celltypeLabel;     /* Optimized access to celltype label */
212:   void                    *ctx;    /* a user context */
213:   PetscErrorCode          (*ctxdestroy)(void**);
214:   ISColoringType          coloringtype;
215:   MatFDColoring           fd;
216:   VecType                 vectype;  /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
217:   MatType                 mattype;  /* type of matrix created with DMCreateMatrix() */
218:   PetscInt                bs;
219:   ISLocalToGlobalMapping  ltogmap;
220:   PetscBool               prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
221:   PetscBool               structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
222:   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 */
223:   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 */
224:   PetscBool               setfromoptionscalled;
225:   void                    *data;
226:   /* Hierarchy / Submeshes */
227:   DM                      coarseMesh;
228:   DM                      fineMesh;
229:   DMCoarsenHookLink       coarsenhook; /* For transfering auxiliary problem data to coarser grids */
230:   DMRefineHookLink        refinehook;
231:   DMSubDomainHookLink     subdomainhook;
232:   DMGlobalToLocalHookLink gtolhook;
233:   DMLocalToGlobalHookLink ltoghook;
234:   /* Topology */
235:   PetscInt                dim;                  /* The topological dimension */
236:   /* Auxiliary data */
237:   PetscHMapAux            auxData;              /* Auxiliary DM and Vec for region denoted by the key */
238:   /* Flexible communication */
239:   PetscSF                 sfMigration;          /* SF for point distribution created during distribution */
240:   PetscSF                 sf;                   /* SF for parallel point overlap */
241:   PetscSF                 sectionSF;            /* SF for parallel dof overlap using section */
242:   PetscSF                 sfNatural;            /* SF mapping to the "natural" ordering */
243:   PetscBool               useNatural;           /* Create the natural SF */
244:   /* Allows a non-standard data layout */
245:   PetscBool               adjacency[2];         /* [use cone() or support() first, use the transitive closure] */
246:   PetscSection            localSection;         /* Layout for local vectors */
247:   PetscSection            globalSection;        /* Layout for global vectors */
248:   PetscLayout             map;
249:   /* Constraints */
250:   PetscSection            defaultConstraintSection;
251:   Mat                     defaultConstraintMat;
252:   /* Basis transformation */
253:   DM                      transformDM;          /* Layout for basis transformation */
254:   Vec                     transform;            /* Basis transformation matrices */
255:   void                   *transformCtx;         /* Basis transformation context */
256:   PetscErrorCode        (*transformSetUp)(DM, void *);
257:   PetscErrorCode        (*transformDestroy)(DM, void *);
258:   PetscErrorCode        (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
259:   /* Coordinates */
260:   PetscInt                dimEmbed;             /* The dimension of the embedding space */
261:   DM                      coordinateDM;         /* Layout for coordinates */
262:   Vec                     coordinates;          /* Coordinate values in global vector */
263:   Vec                     coordinatesLocal;     /* Coordinate values in local  vector */
264:   PetscBool               periodic;             /* Is the DM periodic? */
265:   DMField                 coordinateField;      /* Coordinates as an abstract field */
266:   PetscReal              *L, *maxCell;          /* Size of periodic box and max cell size for determining periodicity */
267:   DMBoundaryType         *bdtype;               /* Indicates type of topological boundary */
268:   /* Null spaces -- of course I should make this have a variable number of fields */
269:   NullSpaceFunc           nullspaceConstructors[10];
270:   NullSpaceFunc           nearnullspaceConstructors[10];
271:   /* Fields are represented by objects */
272:   PetscInt                Nf;                   /* Number of fields defined on the total domain */
273:   RegionField            *fields;               /* Array of discretization fields with regions of validity */
274:   DMBoundary              boundary;             /* List of boundary conditions */
275:   PetscInt                Nds;                  /* Number of discrete systems defined on the total domain */
276:   DMSpace                *probs;                /* Array of discrete systems */
277:   /* Output structures */
278:   DM                      dmBC;                 /* The DM with boundary conditions in the global DM */
279:   PetscInt                outputSequenceNum;    /* The current sequence number for output */
280:   PetscReal               outputSequenceVal;    /* The current sequence value for output */
281:   PetscErrorCode        (*monitor[MAXDMMONITORS])(DM, void *);
282:   PetscErrorCode        (*monitordestroy[MAXDMMONITORS])(void **);
283:   void                   *monitorcontext[MAXDMMONITORS];
284:   PetscInt                numbermonitors;

286:   PetscObject             dmksp,dmsnes,dmts;
287: #ifdef PETSC_HAVE_LIBCEED
288:   Ceed                    ceed;                 /* LibCEED context */
289:   CeedElemRestriction     ceedERestrict;        /* Map from the local vector (Lvector) to the cells (Evector) */
290: #endif
291: };

293: PETSC_EXTERN PetscLogEvent DM_Convert;
294: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
295: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
296: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
297: PETSC_EXTERN PetscLogEvent DM_Coarsen;
298: PETSC_EXTERN PetscLogEvent DM_Refine;
299: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
300: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
301: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
302: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
303: PETSC_EXTERN PetscLogEvent DM_Load;
304: PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;

306: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
307: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);

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

311: /*

313:           Composite Vectors

315:       Single global representation
316:       Individual global representations
317:       Single local representation
318:       Individual local representations

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

322:        DM da_u, da_v, da_p

324:        DM dm_velocities

326:        DM dm

328:        DMDACreate(,&da_u);
329:        DMDACreate(,&da_v);
330:        DMCompositeCreate(,&dm_velocities);
331:        DMCompositeAddDM(dm_velocities,(DM)du);
332:        DMCompositeAddDM(dm_velocities,(DM)dv);

334:        DMDACreate(,&da_p);
335:        DMCompositeCreate(,&dm_velocities);
336:        DMCompositeAddDM(dm,(DM)dm_velocities);
337:        DMCompositeAddDM(dm,(DM)dm_p);

339:     Access parts of composite vectors (Composite only)
340:     ---------------------------------
341:       DMCompositeGetAccess  - access the global vector as subvectors and array (for redundant arrays)
342:       ADD for local vector -

344:     Element access
345:     --------------
346:       From global vectors
347:          -DAVecGetArray   - for DMDA
348:          -VecGetArray - for DMSliced
349:          ADD for DMComposite???  maybe

351:       From individual vector
352:           -DAVecGetArray - for DMDA
353:           -VecGetArray -for sliced
354:          ADD for DMComposite??? maybe

356:       From single local vector
357:           ADD         * single local vector as arrays?

359:    Communication
360:    -------------
361:       DMGlobalToLocal - global vector to single local vector

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

365:    Obtaining vectors
366:    -----------------
367:       DMCreateGlobal/Local
368:       DMGetGlobal/Local
369:       DMCompositeGetLocalVectors   - gives individual local work vectors and arrays

371:     ?????   individual global vectors   ????

373: */

375: #if defined(PETSC_HAVE_HDF5)
376: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
377: #endif

379: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
380: {
382: #if defined(PETSC_USE_DEBUG)
383:   {
384:     PetscInt       dof;
386:     *start = *end = 0; /* Silence overzealous compiler warning */
387:     if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
388:     PetscSectionGetOffset(dm->localSection, point, start);
389:     PetscSectionGetDof(dm->localSection, point, &dof);
390:     *end = *start + dof;
391:   }
392: #else
393:   {
394:     const PetscSection s = dm->localSection;
395:     *start = s->atlasOff[point - s->pStart];
396:     *end   = *start + s->atlasDof[point - s->pStart];
397:   }
398: #endif
399:   return(0);
400: }

402: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
403: {
405: #if defined(PETSC_USE_DEBUG)
406:   {
407:     PetscInt       dof;
409:     *start = *end = 0; /* Silence overzealous compiler warning */
410:     if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
411:     PetscSectionGetFieldOffset(dm->localSection, point, field, start);
412:     PetscSectionGetFieldDof(dm->localSection, point, field, &dof);
413:     *end = *start + dof;
414:   }
415: #else
416:   {
417:     const PetscSection s = dm->localSection->field[field];
418:     *start = s->atlasOff[point - s->pStart];
419:     *end   = *start + s->atlasDof[point - s->pStart];
420:   }
421: #endif
422:   return(0);
423: }

425: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
426: {
428: #if defined(PETSC_USE_DEBUG)
429:   {
431:     PetscInt       dof,cdof;
432:     *start = *end = 0; /* Silence overzealous compiler warning */
433:     if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
434:     if (!dm->globalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
435:     PetscSectionGetOffset(dm->globalSection, point, start);
436:     PetscSectionGetDof(dm->globalSection, point, &dof);
437:     PetscSectionGetConstraintDof(dm->globalSection, point, &cdof);
438:     *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
439:   }
440: #else
441:   {
442:     const PetscSection s    = dm->globalSection;
443:     const PetscInt     dof  = s->atlasDof[point - s->pStart];
444:     const PetscInt     cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
445:     *start = s->atlasOff[point - s->pStart];
446:     *end   = *start + dof - cdof + (dof < 0 ? 1 : 0);
447:   }
448: #endif
449:   return(0);
450: }

452: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
453: {
455: #if defined(PETSC_USE_DEBUG)
456:   {
457:     PetscInt       loff, lfoff, fdof, fcdof, ffcdof, f;
459:     *start = *end = 0; /* Silence overzealous compiler warning */
460:     if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
461:     if (!dm->globalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be crated automatically by DMGetGlobalSection()");
462:     PetscSectionGetOffset(dm->globalSection, point, start);
463:     PetscSectionGetOffset(dm->localSection, point, &loff);
464:     PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff);
465:     PetscSectionGetFieldDof(dm->localSection, point, field, &fdof);
466:     PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof);
467:     *start = *start < 0 ? *start - (lfoff-loff) : *start + lfoff-loff;
468:     for (f = 0; f < field; ++f) {
469:       PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof);
470:       *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
471:     }
472:     *end   = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
473:   }
474: #else
475:   {
476:     const PetscSection s     = dm->localSection;
477:     const PetscSection fs    = dm->localSection->field[field];
478:     const PetscSection gs    = dm->globalSection;
479:     const PetscInt     loff  = s->atlasOff[point - s->pStart];
480:     const PetscInt     goff  = gs->atlasOff[point - s->pStart];
481:     const PetscInt     lfoff = fs->atlasOff[point - s->pStart];
482:     const PetscInt     fdof  = fs->atlasDof[point - s->pStart];
483:     const PetscInt     fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
484:     PetscInt           ffcdof = 0, f;

486:     for (f = 0; f < field; ++f) {
487:       const PetscSection ffs = dm->localSection->field[f];
488:       ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
489:     }
490:     *start = goff + (goff < 0 ? loff-lfoff + ffcdof : lfoff-loff - ffcdof);
491:     *end   = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
492:   }
493: #endif
494:   return(0);
495: }

497: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
498: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
499: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);

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

504: PETSC_INTERN PetscErrorCode DMSetLabelValue_Fast(DM, DMLabel*, const char[], PetscInt, PetscInt);

506: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
507: PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
508: PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
509: PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
510: PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);

512: #endif