Actual source code: dmimpl.h

petsc-master 2019-05-18
Report Typos and Errors


  3: #if !defined(_DMIMPL_H)
  4: #define _DMIMPL_H

  6:  #include <petscdm.h>
  7:  #include <petsc/private/petscimpl.h>
  8:  #include <petsc/private/petscdsimpl.h>
  9:  #include <petsc/private/isimpl.h>

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

 15: typedef struct _DMOps *DMOps;
 16: struct _DMOps {
 17:   PetscErrorCode (*view)(DM,PetscViewer);
 18:   PetscErrorCode (*load)(DM,PetscViewer);
 19:   PetscErrorCode (*clone)(DM,DM*);
 20:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,DM);
 21:   PetscErrorCode (*setup)(DM);
 22:   PetscErrorCode (*createdefaultsection)(DM);
 23:   PetscErrorCode (*createdefaultconstraints)(DM);
 24:   PetscErrorCode (*createglobalvector)(DM,Vec*);
 25:   PetscErrorCode (*createlocalvector)(DM,Vec*);
 26:   PetscErrorCode (*getlocaltoglobalmapping)(DM);
 27:   PetscErrorCode (*createfieldis)(DM,PetscInt*,char***,IS**);
 28:   PetscErrorCode (*createcoordinatedm)(DM,DM*);
 29:   PetscErrorCode (*createcoordinatefield)(DM,DMField*);

 31:   PetscErrorCode (*getcoloring)(DM,ISColoringType,ISColoring*);
 32:   PetscErrorCode (*creatematrix)(DM, Mat*);
 33:   PetscErrorCode (*createinterpolation)(DM,DM,Mat*,Vec*);
 34:   PetscErrorCode (*createrestriction)(DM,DM,Mat*);
 35:   PetscErrorCode (*createmassmatrix)(DM,DM,Mat*);
 36:   PetscErrorCode (*getaggregates)(DM,DM,Mat*);
 37:   PetscErrorCode (*hascreateinjection)(DM,PetscBool*);
 38:   PetscErrorCode (*getinjection)(DM,DM,Mat*);

 40:   PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
 41:   PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*);
 42:   PetscErrorCode (*refinehierarchy)(DM,PetscInt,DM*);
 43:   PetscErrorCode (*coarsenhierarchy)(DM,PetscInt,DM*);
 44:   PetscErrorCode (*adaptlabel)(DM,DMLabel,DM*);
 45:   PetscErrorCode (*adaptmetric)(DM,Vec,DMLabel,DM*);

 47:   PetscErrorCode (*globaltolocalbegin)(DM,Vec,InsertMode,Vec);
 48:   PetscErrorCode (*globaltolocalend)(DM,Vec,InsertMode,Vec);
 49:   PetscErrorCode (*localtoglobalbegin)(DM,Vec,InsertMode,Vec);
 50:   PetscErrorCode (*localtoglobalend)(DM,Vec,InsertMode,Vec);
 51:   PetscErrorCode (*localtolocalbegin)(DM,Vec,InsertMode,Vec);
 52:   PetscErrorCode (*localtolocalend)(DM,Vec,InsertMode,Vec);

 54:   PetscErrorCode (*destroy)(DM);

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

 58:   PetscErrorCode (*createsubdm)(DM,PetscInt,const PetscInt*,IS*,DM*);
 59:   PetscErrorCode (*createsuperdm)(DM*,PetscInt,IS**,DM*);
 60:   PetscErrorCode (*createfielddecomposition)(DM,PetscInt*,char***,IS**,DM**);
 61:   PetscErrorCode (*createdomaindecomposition)(DM,PetscInt*,char***,IS**,IS**,DM**);
 62:   PetscErrorCode (*createddscatters)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

 64:   PetscErrorCode (*getdimpoints)(DM,PetscInt,PetscInt*,PetscInt*);
 65:   PetscErrorCode (*locatepoints)(DM,Vec,DMPointLocationType,PetscSF);
 66:   PetscErrorCode (*getneighbors)(DM,PetscInt*,const PetscMPIInt**);
 67:   PetscErrorCode (*getboundingbox)(DM,PetscReal*,PetscReal*);
 68:   PetscErrorCode (*getlocalboundingbox)(DM,PetscReal*,PetscReal*);
 69:   PetscErrorCode (*locatepointssubdomain)(DM,Vec,PetscMPIInt**);

 71:   PetscErrorCode (*projectfunctionlocal)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
 72:   PetscErrorCode (*projectfunctionlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
 73:   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);
 74:   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);
 75:   PetscErrorCode (*computel2diff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
 76:   PetscErrorCode (*computel2gradientdiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal [],const PetscReal[],PetscInt, PetscScalar *,void *),void **,Vec,const PetscReal[],PetscReal *);
 77:   PetscErrorCode (*computel2fielddiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);

 79:   PetscErrorCode (*getcompatibility)(DM,DM,PetscBool*,PetscBool*);
 80: };

 82: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
 83: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
 84: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);

 86: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
 87: struct _DMCoarsenHookLink {
 88:   PetscErrorCode (*coarsenhook)(DM,DM,void*);              /* Run once, when coarse DM is created */
 89:   PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*); /* Run each time a new problem is restricted to a coarse grid */
 90:   void *ctx;
 91:   DMCoarsenHookLink next;
 92: };

 94: typedef struct _DMRefineHookLink *DMRefineHookLink;
 95: struct _DMRefineHookLink {
 96:   PetscErrorCode (*refinehook)(DM,DM,void*);     /* Run once, when a fine DM is created */
 97:   PetscErrorCode (*interphook)(DM,Mat,DM,void*); /* Run each time a new problem is interpolated to a fine grid */
 98:   void *ctx;
 99:   DMRefineHookLink next;
100: };

102: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
103: struct _DMSubDomainHookLink {
104:   PetscErrorCode (*ddhook)(DM,DM,void*);
105:   PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*);
106:   void *ctx;
107:   DMSubDomainHookLink next;
108: };

110: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
111: struct _DMGlobalToLocalHookLink {
112:   PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
113:   PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
114:   void *ctx;
115:   DMGlobalToLocalHookLink next;
116: };

118: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
119: struct _DMLocalToGlobalHookLink {
120:   PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
121:   PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
122:   void *ctx;
123:   DMLocalToGlobalHookLink next;
124: };

126: typedef enum {DMVEC_STATUS_IN,DMVEC_STATUS_OUT} DMVecStatus;
127: typedef struct _DMNamedVecLink *DMNamedVecLink;
128: struct _DMNamedVecLink {
129:   Vec X;
130:   char *name;
131:   DMVecStatus status;
132:   DMNamedVecLink next;
133: };

135: typedef struct _DMWorkLink *DMWorkLink;
136: struct _DMWorkLink {
137:   size_t     bytes;
138:   void       *mem;
139:   DMWorkLink next;
140: };

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

144: struct _n_DMLabelLink {
145:   DMLabel              label;
146:   PetscBool            output;
147:   struct _n_DMLabelLink *next;
148: };
149: typedef struct _n_DMLabelLink *DMLabelLink;

151: struct _n_DMLabelLinkList {
152:   PetscInt refct;
153:   DMLabelLink next;
154: };
155: typedef struct _n_DMLabelLinkList *DMLabelLinkList;

157: typedef struct _n_Boundary *DMBoundary;

159: struct _n_Boundary {
160:   DSBoundary  dsboundary;
161:   DMLabel     label;
162:   DMBoundary  next;
163: };

165: typedef struct _n_Field {
166:   PetscObject disc;         /* Field discretization, or a PetscContainer with the field name */
167:   DMLabel     label;        /* Label defining the domain of definition of the field */
168:   PetscBool   adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
169: } RegionField;

171: typedef struct _n_Space {
172:   PetscDS ds;     /* Approximation space in this domain */
173:   DMLabel label;  /* Label defining the domain of definition of the discretization */
174:   IS      fields; /* Map from DS field numbers to original field numbers in the DM */
175: } DMSpace;

177: PETSC_EXTERN PetscErrorCode DMDestroyLabelLinkList(DM);

179: struct _p_DM {
180:   PETSCHEADER(struct _DMOps);
181:   Vec                     localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
182:   Vec                     globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
183:   DMNamedVecLink          namedglobal;
184:   DMNamedVecLink          namedlocal;
185:   DMWorkLink              workin,workout;
186:   DMLabelLinkList         labels;            /* Linked list of labels */
187:   DMLabel                 depthLabel;        /* Optimized access to depth label */
188:   void                    *ctx;    /* a user context */
189:   PetscErrorCode          (*ctxdestroy)(void**);
190:   Vec                     x;       /* location at which the functions/Jacobian are computed */
191:   ISColoringType          coloringtype;
192:   MatFDColoring           fd;
193:   VecType                 vectype;  /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
194:   MatType                 mattype;  /* type of matrix created with DMCreateMatrix() */
195:   PetscInt                bs;
196:   ISLocalToGlobalMapping  ltogmap;
197:   PetscBool               prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
198:   PetscBool               structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
199:   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 */
200:   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 */
201:   PetscBool               setfromoptionscalled;
202:   void                    *data;
203:   /* Hierarchy / Submeshes */
204:   DM                      coarseMesh;
205:   DM                      fineMesh;
206:   DMCoarsenHookLink       coarsenhook; /* For transfering auxiliary problem data to coarser grids */
207:   DMRefineHookLink        refinehook;
208:   DMSubDomainHookLink     subdomainhook;
209:   DMGlobalToLocalHookLink gtolhook;
210:   DMLocalToGlobalHookLink ltoghook;
211:   /* Topology */
212:   PetscInt                dim;                  /* The topological dimension */
213:   /* Flexible communication */
214:   PetscSF                 sfMigration;          /* SF for point distribution created during distribution */
215:   PetscSF                 sf;                   /* SF for parallel point overlap */
216:   PetscSF                 defaultSF;            /* SF for parallel dof overlap using default section */
217:   PetscSF                 sfNatural;            /* SF mapping to the "natural" ordering */
218:   PetscBool               useNatural;           /* Create the natural SF */
219:   /* Allows a non-standard data layout */
220:   PetscBool               adjacency[2];         /* [use cone() or support() first, use the transitive closure] */
221:   PetscSection            defaultSection;       /* Layout for local vectors */
222:   PetscSection            defaultGlobalSection; /* Layout for global vectors */
223:   PetscLayout             map;
224:   /* Constraints */
225:   PetscSection            defaultConstraintSection;
226:   Mat                     defaultConstraintMat;
227:   /* Basis transformation */
228:   DM                      transformDM;          /* Layout for basis transformation */
229:   Vec                     transform;            /* Basis transformation matrices */
230:   void                   *transformCtx;         /* Basis transformation context */
231:   PetscErrorCode        (*transformSetUp)(DM, void *);
232:   PetscErrorCode        (*transformDestroy)(DM, void *);
233:   PetscErrorCode        (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
234:   /* Coordinates */
235:   PetscInt                dimEmbed;             /* The dimension of the embedding space */
236:   DM                      coordinateDM;         /* Layout for coordinates (default section) */
237:   Vec                     coordinates;          /* Coordinate values in global vector */
238:   Vec                     coordinatesLocal;     /* Coordinate values in local  vector */
239:   PetscBool               periodic;             /* Is the DM periodic? */
240:   DMField                 coordinateField;      /* Coordinates as an abstract field */
241:   PetscReal              *L, *maxCell;          /* Size of periodic box and max cell size for determining periodicity */
242:   DMBoundaryType         *bdtype;               /* Indicates type of topological boundary */
243:   /* Null spaces -- of course I should make this have a variable number of fields */
244:   /*   I now believe this might not be the right way: see below */
245:   NullSpaceFunc           nullspaceConstructors[10];
246:   NullSpaceFunc           nearnullspaceConstructors[10];
247:   /* Fields are represented by objects */
248:   PetscInt                Nf;                   /* Number of fields defined on the total domain */
249:   RegionField            *fields;               /* Array of discretization fields with regions of validity */
250:   DMBoundary              boundary;             /* List of boundary conditions */
251:   PetscInt                Nds;                  /* Number of discrete systems defined on the total domain */
252:   DMSpace                *probs;                /* Array of discrete systems */
253:   /* Output structures */
254:   DM                      dmBC;                 /* The DM with boundary conditions in the global DM */
255:   PetscInt                outputSequenceNum;    /* The current sequence number for output */
256:   PetscReal               outputSequenceVal;    /* The current sequence value for output */

258:   PetscObject             dmksp,dmsnes,dmts;
259: };

261: PETSC_EXTERN PetscLogEvent DM_Convert;
262: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
263: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
264: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
265: PETSC_EXTERN PetscLogEvent DM_Coarsen;
266: PETSC_EXTERN PetscLogEvent DM_Refine;
267: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
268: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;

270: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
271: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);

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

275: /*

277:           Composite Vectors

279:       Single global representation
280:       Individual global representations
281:       Single local representation
282:       Individual local representations

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

286:        DM da_u, da_v, da_p

288:        DM dm_velocities

290:        DM dm

292:        DMDACreate(,&da_u);
293:        DMDACreate(,&da_v);
294:        DMCompositeCreate(,&dm_velocities);
295:        DMCompositeAddDM(dm_velocities,(DM)du);
296:        DMCompositeAddDM(dm_velocities,(DM)dv);

298:        DMDACreate(,&da_p);
299:        DMCompositeCreate(,&dm_velocities);
300:        DMCompositeAddDM(dm,(DM)dm_velocities);
301:        DMCompositeAddDM(dm,(DM)dm_p);


304:     Access parts of composite vectors (Composite only)
305:     ---------------------------------
306:       DMCompositeGetAccess  - access the global vector as subvectors and array (for redundant arrays)
307:       ADD for local vector -

309:     Element access
310:     --------------
311:       From global vectors
312:          -DAVecGetArray   - for DMDA
313:          -VecGetArray - for DMSliced
314:          ADD for DMComposite???  maybe

316:       From individual vector
317:           -DAVecGetArray - for DMDA
318:           -VecGetArray -for sliced
319:          ADD for DMComposite??? maybe

321:       From single local vector
322:           ADD         * single local vector as arrays?

324:    Communication
325:    -------------
326:       DMGlobalToLocal - global vector to single local vector

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

330:    Obtaining vectors
331:    -----------------
332:       DMCreateGlobal/Local
333:       DMGetGlobal/Local
334:       DMCompositeGetLocalVectors   - gives individual local work vectors and arrays


337: ?????   individual global vectors   ????

339: */

341: #if defined(PETSC_HAVE_HDF5)
342: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
343: #endif

345: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
346: {
348: #if defined(PETSC_USE_DEBUG)
349:   {
350:     PetscInt       dof;
352:     *start = *end = 0; /* Silence overzealous compiler warning */
353:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
354:     PetscSectionGetOffset(dm->defaultSection, point, start);
355:     PetscSectionGetDof(dm->defaultSection, point, &dof);
356:     *end = *start + dof;
357:   }
358: #else
359:   {
360:     const PetscSection s = dm->defaultSection;
361:     *start = s->atlasOff[point - s->pStart];
362:     *end   = *start + s->atlasDof[point - s->pStart];
363:   }
364: #endif
365:   return(0);
366: }

368: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
369: {
371: #if defined(PETSC_USE_DEBUG)
372:   {
373:     PetscInt       dof;
375:     *start = *end = 0; /* Silence overzealous compiler warning */
376:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
377:     PetscSectionGetFieldOffset(dm->defaultSection, point, field, start);
378:     PetscSectionGetFieldDof(dm->defaultSection, point, field, &dof);
379:     *end = *start + dof;
380:   }
381: #else
382:   {
383:     const PetscSection s = dm->defaultSection->field[field];
384:     *start = s->atlasOff[point - s->pStart];
385:     *end   = *start + s->atlasDof[point - s->pStart];
386:   }
387: #endif
388:   return(0);
389: }

391: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
392: {
394: #if defined(PETSC_USE_DEBUG)
395:   {
397:     PetscInt       dof,cdof;
398:     *start = *end = 0; /* Silence overzealous compiler warning */
399:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
400:     if (!dm->defaultGlobalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default global Section. It will be created automatically by DMGetDefaultGlobalSection()");
401:     PetscSectionGetOffset(dm->defaultGlobalSection, point, start);
402:     PetscSectionGetDof(dm->defaultGlobalSection, point, &dof);
403:     PetscSectionGetConstraintDof(dm->defaultGlobalSection, point, &cdof);
404:     *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
405:   }
406: #else
407:   {
408:     const PetscSection s    = dm->defaultGlobalSection;
409:     const PetscInt     dof  = s->atlasDof[point - s->pStart];
410:     const PetscInt     cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
411:     *start = s->atlasOff[point - s->pStart];
412:     *end   = *start + dof - cdof + (dof < 0 ? 1 : 0);
413:   }
414: #endif
415:   return(0);
416: }

418: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
419: {
421: #if defined(PETSC_USE_DEBUG)
422:   {
423:     PetscInt       loff, lfoff, fdof, fcdof, ffcdof, f;
425:     *start = *end = 0; /* Silence overzealous compiler warning */
426:     if (!dm->defaultSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default Section, see DMSetDefaultSection()");
427:     if (!dm->defaultGlobalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a default global Section. It will be crated automatically by DMGetDefaultGlobalSection()");
428:     PetscSectionGetOffset(dm->defaultGlobalSection, point, start);
429:     PetscSectionGetOffset(dm->defaultSection, point, &loff);
430:     PetscSectionGetFieldOffset(dm->defaultSection, point, field, &lfoff);
431:     PetscSectionGetFieldDof(dm->defaultSection, point, field, &fdof);
432:     PetscSectionGetFieldConstraintDof(dm->defaultSection, point, field, &fcdof);
433:     *start = *start < 0 ? *start - (lfoff-loff) : *start + lfoff-loff;
434:     for (f = 0; f < field; ++f) {
435:       PetscSectionGetFieldConstraintDof(dm->defaultSection, point, f, &ffcdof);
436:       *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
437:     }
438:     *end   = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
439:   }
440: #else
441:   {
442:     const PetscSection s     = dm->defaultSection;
443:     const PetscSection fs    = dm->defaultSection->field[field];
444:     const PetscSection gs    = dm->defaultGlobalSection;
445:     const PetscInt     loff  = s->atlasOff[point - s->pStart];
446:     const PetscInt     goff  = gs->atlasOff[point - s->pStart];
447:     const PetscInt     lfoff = fs->atlasOff[point - s->pStart];
448:     const PetscInt     fdof  = fs->atlasDof[point - s->pStart];
449:     const PetscInt     fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
450:     PetscInt           ffcdof = 0, f;

452:     for (f = 0; f < field; ++f) {
453:       const PetscSection ffs = dm->defaultSection->field[f];
454:       ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
455:     }
456:     *start = goff + (goff < 0 ? loff-lfoff + ffcdof : lfoff-loff - ffcdof);
457:     *end   = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
458:   }
459: #endif
460:   return(0);
461: }

463: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
464: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
465: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);

467: #endif