Actual source code: dmimpl.h
petsc-3.13.4 2020-08-01
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/sectionimpl.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 (*createlocalsection)(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 (*hascreateinjection)(DM,PetscBool*);
37: PetscErrorCode (*createinjection)(DM,DM,Mat*);
39: PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
40: PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*);
41: PetscErrorCode (*refinehierarchy)(DM,PetscInt,DM*);
42: PetscErrorCode (*coarsenhierarchy)(DM,PetscInt,DM*);
43: PetscErrorCode (*adaptlabel)(DM,DMLabel,DM*);
44: PetscErrorCode (*adaptmetric)(DM,Vec,DMLabel,DM*);
46: PetscErrorCode (*globaltolocalbegin)(DM,Vec,InsertMode,Vec);
47: PetscErrorCode (*globaltolocalend)(DM,Vec,InsertMode,Vec);
48: PetscErrorCode (*localtoglobalbegin)(DM,Vec,InsertMode,Vec);
49: PetscErrorCode (*localtoglobalend)(DM,Vec,InsertMode,Vec);
50: PetscErrorCode (*localtolocalbegin)(DM,Vec,InsertMode,Vec);
51: PetscErrorCode (*localtolocalend)(DM,Vec,InsertMode,Vec);
53: PetscErrorCode (*destroy)(DM);
55: PetscErrorCode (*computevariablebounds)(DM,Vec,Vec);
57: PetscErrorCode (*createsubdm)(DM,PetscInt,const PetscInt*,IS*,DM*);
58: PetscErrorCode (*createsuperdm)(DM*,PetscInt,IS**,DM*);
59: PetscErrorCode (*createfielddecomposition)(DM,PetscInt*,char***,IS**,DM**);
60: PetscErrorCode (*createdomaindecomposition)(DM,PetscInt*,char***,IS**,IS**,DM**);
61: PetscErrorCode (*createddscatters)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
63: PetscErrorCode (*getdimpoints)(DM,PetscInt,PetscInt*,PetscInt*);
64: PetscErrorCode (*locatepoints)(DM,Vec,DMPointLocationType,PetscSF);
65: PetscErrorCode (*getneighbors)(DM,PetscInt*,const PetscMPIInt**);
66: PetscErrorCode (*getboundingbox)(DM,PetscReal*,PetscReal*);
67: PetscErrorCode (*getlocalboundingbox)(DM,PetscReal*,PetscReal*);
68: PetscErrorCode (*locatepointssubdomain)(DM,Vec,PetscMPIInt**);
70: PetscErrorCode (*projectfunctionlocal)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
71: PetscErrorCode (*projectfunctionlabellocal)(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
72: 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);
73: 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);
74: PetscErrorCode (*computel2diff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
75: PetscErrorCode (*computel2gradientdiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal [],const PetscReal[],PetscInt, PetscScalar *,void *),void **,Vec,const PetscReal[],PetscReal *);
76: PetscErrorCode (*computel2fielddiff)(DM,PetscReal,PetscErrorCode(**)(PetscInt, PetscReal,const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
78: PetscErrorCode (*getcompatibility)(DM,DM,PetscBool*,PetscBool*);
79: };
81: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
82: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
83: PETSC_EXTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
85: typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
86: struct _DMCoarsenHookLink {
87: PetscErrorCode (*coarsenhook)(DM,DM,void*); /* Run once, when coarse DM is created */
88: PetscErrorCode (*restricthook)(DM,Mat,Vec,Mat,DM,void*); /* Run each time a new problem is restricted to a coarse grid */
89: void *ctx;
90: DMCoarsenHookLink next;
91: };
93: typedef struct _DMRefineHookLink *DMRefineHookLink;
94: struct _DMRefineHookLink {
95: PetscErrorCode (*refinehook)(DM,DM,void*); /* Run once, when a fine DM is created */
96: PetscErrorCode (*interphook)(DM,Mat,DM,void*); /* Run each time a new problem is interpolated to a fine grid */
97: void *ctx;
98: DMRefineHookLink next;
99: };
101: typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
102: struct _DMSubDomainHookLink {
103: PetscErrorCode (*ddhook)(DM,DM,void*);
104: PetscErrorCode (*restricthook)(DM,VecScatter,VecScatter,DM,void*);
105: void *ctx;
106: DMSubDomainHookLink next;
107: };
109: typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
110: struct _DMGlobalToLocalHookLink {
111: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
112: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
113: void *ctx;
114: DMGlobalToLocalHookLink next;
115: };
117: typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
118: struct _DMLocalToGlobalHookLink {
119: PetscErrorCode (*beginhook)(DM,Vec,InsertMode,Vec,void*);
120: PetscErrorCode (*endhook)(DM,Vec,InsertMode,Vec,void*);
121: void *ctx;
122: DMLocalToGlobalHookLink next;
123: };
125: typedef enum {DMVEC_STATUS_IN,DMVEC_STATUS_OUT} DMVecStatus;
126: typedef struct _DMNamedVecLink *DMNamedVecLink;
127: struct _DMNamedVecLink {
128: Vec X;
129: char *name;
130: DMVecStatus status;
131: DMNamedVecLink next;
132: };
134: typedef struct _DMWorkLink *DMWorkLink;
135: struct _DMWorkLink {
136: size_t bytes;
137: void *mem;
138: DMWorkLink next;
139: };
141: #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
143: struct _n_DMLabelLink {
144: DMLabel label;
145: PetscBool output;
146: struct _n_DMLabelLink *next;
147: };
148: typedef struct _n_DMLabelLink *DMLabelLink;
150: typedef struct _n_Boundary *DMBoundary;
152: struct _n_Boundary {
153: DSBoundary dsboundary;
154: DMLabel label;
155: DMBoundary next;
156: };
158: typedef struct _n_Field {
159: PetscObject disc; /* Field discretization, or a PetscContainer with the field name */
160: DMLabel label; /* Label defining the domain of definition of the field */
161: PetscBool adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
162: } RegionField;
164: typedef struct _n_Space {
165: PetscDS ds; /* Approximation space in this domain */
166: DMLabel label; /* Label defining the domain of definition of the discretization */
167: IS fields; /* Map from DS field numbers to original field numbers in the DM */
168: } DMSpace;
170: PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);
172: #define MAXDMMONITORS 5
174: struct _p_DM {
175: PETSCHEADER(struct _DMOps);
176: Vec localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
177: Vec globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
178: DMNamedVecLink namedglobal;
179: DMNamedVecLink namedlocal;
180: DMWorkLink workin,workout;
181: DMLabelLink labels; /* Linked list of labels */
182: DMLabel depthLabel; /* Optimized access to depth label */
183: DMLabel celltypeLabel; /* Optimized access to celltype label */
184: void *ctx; /* a user context */
185: PetscErrorCode (*ctxdestroy)(void**);
186: ISColoringType coloringtype;
187: MatFDColoring fd;
188: VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
189: MatType mattype; /* type of matrix created with DMCreateMatrix() */
190: PetscInt bs;
191: ISLocalToGlobalMapping ltogmap;
192: PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
193: PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
194: 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 */
195: 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 */
196: PetscBool setfromoptionscalled;
197: void *data;
198: /* Hierarchy / Submeshes */
199: DM coarseMesh;
200: DM fineMesh;
201: DMCoarsenHookLink coarsenhook; /* For transfering auxiliary problem data to coarser grids */
202: DMRefineHookLink refinehook;
203: DMSubDomainHookLink subdomainhook;
204: DMGlobalToLocalHookLink gtolhook;
205: DMLocalToGlobalHookLink ltoghook;
206: /* Topology */
207: PetscInt dim; /* The topological dimension */
208: /* Flexible communication */
209: PetscSF sfMigration; /* SF for point distribution created during distribution */
210: PetscSF sf; /* SF for parallel point overlap */
211: PetscSF sectionSF; /* SF for parallel dof overlap using section */
212: PetscSF sfNatural; /* SF mapping to the "natural" ordering */
213: PetscBool useNatural; /* Create the natural SF */
214: /* Allows a non-standard data layout */
215: PetscBool adjacency[2]; /* [use cone() or support() first, use the transitive closure] */
216: PetscSection localSection; /* Layout for local vectors */
217: PetscSection globalSection; /* Layout for global vectors */
218: PetscLayout map;
219: /* Constraints */
220: PetscSection defaultConstraintSection;
221: Mat defaultConstraintMat;
222: /* Basis transformation */
223: DM transformDM; /* Layout for basis transformation */
224: Vec transform; /* Basis transformation matrices */
225: void *transformCtx; /* Basis transformation context */
226: PetscErrorCode (*transformSetUp)(DM, void *);
227: PetscErrorCode (*transformDestroy)(DM, void *);
228: PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
229: /* Coordinates */
230: PetscInt dimEmbed; /* The dimension of the embedding space */
231: DM coordinateDM; /* Layout for coordinates */
232: Vec coordinates; /* Coordinate values in global vector */
233: Vec coordinatesLocal; /* Coordinate values in local vector */
234: PetscBool periodic; /* Is the DM periodic? */
235: DMField coordinateField; /* Coordinates as an abstract field */
236: PetscReal *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
237: DMBoundaryType *bdtype; /* Indicates type of topological boundary */
238: /* Null spaces -- of course I should make this have a variable number of fields */
239: /* I now believe this might not be the right way: see below */
240: NullSpaceFunc nullspaceConstructors[10];
241: NullSpaceFunc nearnullspaceConstructors[10];
242: /* Fields are represented by objects */
243: PetscInt Nf; /* Number of fields defined on the total domain */
244: RegionField *fields; /* Array of discretization fields with regions of validity */
245: DMBoundary boundary; /* List of boundary conditions */
246: PetscInt Nds; /* Number of discrete systems defined on the total domain */
247: DMSpace *probs; /* Array of discrete systems */
248: /* Output structures */
249: DM dmBC; /* The DM with boundary conditions in the global DM */
250: PetscInt outputSequenceNum; /* The current sequence number for output */
251: PetscReal outputSequenceVal; /* The current sequence value for output */
252: PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
253: PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
254: void *monitorcontext[MAXDMMONITORS];
255: PetscInt numbermonitors;
257: PetscObject dmksp,dmsnes,dmts;
258: };
260: PETSC_EXTERN PetscLogEvent DM_Convert;
261: PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
262: PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
263: PETSC_EXTERN PetscLogEvent DM_LocatePoints;
264: PETSC_EXTERN PetscLogEvent DM_Coarsen;
265: PETSC_EXTERN PetscLogEvent DM_Refine;
266: PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
267: PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
268: PETSC_EXTERN PetscLogEvent DM_CreateInjection;
269: PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
270: PETSC_EXTERN PetscLogEvent DM_Load;
272: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
273: PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);
275: PETSC_EXTERN PetscErrorCode DMView_GLVis(DM,PetscViewer,PetscErrorCode(*)(DM,PetscViewer));
277: /*
279: Composite Vectors
281: Single global representation
282: Individual global representations
283: Single local representation
284: Individual local representations
286: Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
288: DM da_u, da_v, da_p
290: DM dm_velocities
292: DM dm
294: DMDACreate(,&da_u);
295: DMDACreate(,&da_v);
296: DMCompositeCreate(,&dm_velocities);
297: DMCompositeAddDM(dm_velocities,(DM)du);
298: DMCompositeAddDM(dm_velocities,(DM)dv);
300: DMDACreate(,&da_p);
301: DMCompositeCreate(,&dm_velocities);
302: DMCompositeAddDM(dm,(DM)dm_velocities);
303: DMCompositeAddDM(dm,(DM)dm_p);
306: Access parts of composite vectors (Composite only)
307: ---------------------------------
308: DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
309: ADD for local vector -
311: Element access
312: --------------
313: From global vectors
314: -DAVecGetArray - for DMDA
315: -VecGetArray - for DMSliced
316: ADD for DMComposite??? maybe
318: From individual vector
319: -DAVecGetArray - for DMDA
320: -VecGetArray -for sliced
321: ADD for DMComposite??? maybe
323: From single local vector
324: ADD * single local vector as arrays?
326: Communication
327: -------------
328: DMGlobalToLocal - global vector to single local vector
330: DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
332: Obtaining vectors
333: -----------------
334: DMCreateGlobal/Local
335: DMGetGlobal/Local
336: DMCompositeGetLocalVectors - gives individual local work vectors and arrays
339: ????? individual global vectors ????
341: */
343: #if defined(PETSC_HAVE_HDF5)
344: PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
345: #endif
347: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
348: {
350: #if defined(PETSC_USE_DEBUG)
351: {
352: PetscInt dof;
354: *start = *end = 0; /* Silence overzealous compiler warning */
355: if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
356: PetscSectionGetOffset(dm->localSection, point, start);
357: PetscSectionGetDof(dm->localSection, point, &dof);
358: *end = *start + dof;
359: }
360: #else
361: {
362: const PetscSection s = dm->localSection;
363: *start = s->atlasOff[point - s->pStart];
364: *end = *start + s->atlasDof[point - s->pStart];
365: }
366: #endif
367: return(0);
368: }
370: PETSC_STATIC_INLINE PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
371: {
373: #if defined(PETSC_USE_DEBUG)
374: {
375: PetscInt dof;
377: *start = *end = 0; /* Silence overzealous compiler warning */
378: if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
379: PetscSectionGetFieldOffset(dm->localSection, point, field, start);
380: PetscSectionGetFieldDof(dm->localSection, point, field, &dof);
381: *end = *start + dof;
382: }
383: #else
384: {
385: const PetscSection s = dm->localSection->field[field];
386: *start = s->atlasOff[point - s->pStart];
387: *end = *start + s->atlasDof[point - s->pStart];
388: }
389: #endif
390: return(0);
391: }
393: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
394: {
396: #if defined(PETSC_USE_DEBUG)
397: {
399: PetscInt dof,cdof;
400: *start = *end = 0; /* Silence overzealous compiler warning */
401: if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
402: if (!dm->globalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
403: PetscSectionGetOffset(dm->globalSection, point, start);
404: PetscSectionGetDof(dm->globalSection, point, &dof);
405: PetscSectionGetConstraintDof(dm->globalSection, point, &cdof);
406: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
407: }
408: #else
409: {
410: const PetscSection s = dm->globalSection;
411: const PetscInt dof = s->atlasDof[point - s->pStart];
412: const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
413: *start = s->atlasOff[point - s->pStart];
414: *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
415: }
416: #endif
417: return(0);
418: }
420: PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
421: {
423: #if defined(PETSC_USE_DEBUG)
424: {
425: PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
427: *start = *end = 0; /* Silence overzealous compiler warning */
428: if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
429: if (!dm->globalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be crated automatically by DMGetGlobalSection()");
430: PetscSectionGetOffset(dm->globalSection, point, start);
431: PetscSectionGetOffset(dm->localSection, point, &loff);
432: PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff);
433: PetscSectionGetFieldDof(dm->localSection, point, field, &fdof);
434: PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof);
435: *start = *start < 0 ? *start - (lfoff-loff) : *start + lfoff-loff;
436: for (f = 0; f < field; ++f) {
437: PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof);
438: *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
439: }
440: *end = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
441: }
442: #else
443: {
444: const PetscSection s = dm->localSection;
445: const PetscSection fs = dm->localSection->field[field];
446: const PetscSection gs = dm->globalSection;
447: const PetscInt loff = s->atlasOff[point - s->pStart];
448: const PetscInt goff = gs->atlasOff[point - s->pStart];
449: const PetscInt lfoff = fs->atlasOff[point - s->pStart];
450: const PetscInt fdof = fs->atlasDof[point - s->pStart];
451: const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
452: PetscInt ffcdof = 0, f;
454: for (f = 0; f < field; ++f) {
455: const PetscSection ffs = dm->localSection->field[f];
456: ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
457: }
458: *start = goff + (goff < 0 ? loff-lfoff + ffcdof : lfoff-loff - ffcdof);
459: *end = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
460: }
461: #endif
462: return(0);
463: }
465: PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
466: PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
467: PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);
469: PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
471: #endif