1
2
3 #if !defined(_DMIMPL_H)
4 #define _DMIMPL_H
5
6 #include <petscdm.h>
7 #include <petsc/private/petscimpl.h>
8 #include <petsc/private/petscdsimpl.h>
9 #include <petsc/private/sectionimpl.h> /* for inline access to atlasOff */
10
11 PETSC_EXTERN PetscBool DMRegisterAllCalled;
12 PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
13 typedef PetscErrorCode (*NullSpaceFunc)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace);
14
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*);
30
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*);
38
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*);
45
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);
52
53 PetscErrorCode (*destroy)(DM);
54
55 PetscErrorCode (*computevariablebounds)(DM,Vec,Vec);
56
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**);
62
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**);
69
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 (*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);
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 *);
78
79 PetscErrorCode (*getcompatibility)(DM,DM,PetscBool*,PetscBool*);
80 };
81
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[]);
85
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 };
93
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 };
101
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 };
109
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 };
117
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 };
125
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 };
134
135 typedef struct _DMWorkLink *DMWorkLink;
136 struct _DMWorkLink {
137 size_t bytes;
138 void *mem;
139 DMWorkLink next;
140 };
141
142 #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
143
144 struct _n_DMLabelLink {
145 DMLabel label;
146 PetscBool output;
147 struct _n_DMLabelLink *next;
148 };
149 typedef struct _n_DMLabelLink *DMLabelLink;
150
151 typedef struct _n_Boundary *DMBoundary;
152
153 struct _n_Boundary {
154 DSBoundary dsboundary;
155 DMLabel label;
156 DMBoundary next;
157 };
158
159 typedef struct _n_Field {
160 PetscObject disc; /* Field discretization, or a PetscContainer with the field name */
161 DMLabel label; /* Label defining the domain of definition of the field */
162 PetscBool adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
163 } RegionField;
164
165 typedef struct _n_Space {
166 PetscDS ds; /* Approximation space in this domain */
167 DMLabel label; /* Label defining the domain of definition of the discretization */
168 IS fields; /* Map from DS field numbers to original field numbers in the DM */
169 } DMSpace;
170
171 PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);
172
173 #define MAXDMMONITORS 5
174
175 struct _p_DM {
176 PETSCHEADER(struct _DMOps);
177 Vec localin[DM_MAX_WORK_VECTORS],localout[DM_MAX_WORK_VECTORS];
178 Vec globalin[DM_MAX_WORK_VECTORS],globalout[DM_MAX_WORK_VECTORS];
179 DMNamedVecLink namedglobal;
180 DMNamedVecLink namedlocal;
181 DMWorkLink workin,workout;
182 DMLabelLink labels; /* Linked list of labels */
183 DMLabel depthLabel; /* Optimized access to depth label */
184 DMLabel celltypeLabel; /* Optimized access to celltype label */
185 void *ctx; /* a user context */
186 PetscErrorCode (*ctxdestroy)(void**);
187 ISColoringType coloringtype;
188 MatFDColoring fd;
189 VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
190 MatType mattype; /* type of matrix created with DMCreateMatrix() */
191 PetscInt bs;
192 ISLocalToGlobalMapping ltogmap;
193 PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
194 PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix structure without values */
195 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 */
196 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 */
197 PetscBool setfromoptionscalled;
198 void *data;
199 /* Hierarchy / Submeshes */
200 DM coarseMesh;
201 DM fineMesh;
202 DMCoarsenHookLink coarsenhook; /* For transfering auxiliary problem data to coarser grids */
203 DMRefineHookLink refinehook;
204 DMSubDomainHookLink subdomainhook;
205 DMGlobalToLocalHookLink gtolhook;
206 DMLocalToGlobalHookLink ltoghook;
207 /* Topology */
208 PetscInt dim; /* The topological dimension */
209 /* Flexible communication */
210 PetscSF sfMigration; /* SF for point distribution created during distribution */
211 PetscSF sf; /* SF for parallel point overlap */
212 PetscSF sectionSF; /* SF for parallel dof overlap using section */
213 PetscSF sfNatural; /* SF mapping to the "natural" ordering */
214 PetscBool useNatural; /* Create the natural SF */
215 /* Allows a non-standard data layout */
216 PetscBool adjacency[2]; /* [use cone() or support() first, use the transitive closure] */
217 PetscSection localSection; /* Layout for local vectors */
218 PetscSection globalSection; /* Layout for global vectors */
219 PetscLayout map;
220 /* Constraints */
221 PetscSection defaultConstraintSection;
222 Mat defaultConstraintMat;
223 /* Basis transformation */
224 DM transformDM; /* Layout for basis transformation */
225 Vec transform; /* Basis transformation matrices */
226 void *transformCtx; /* Basis transformation context */
227 PetscErrorCode (*transformSetUp)(DM, void *);
228 PetscErrorCode (*transformDestroy)(DM, void *);
229 PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
230 /* Coordinates */
231 PetscInt dimEmbed; /* The dimension of the embedding space */
232 DM coordinateDM; /* Layout for coordinates */
233 Vec coordinates; /* Coordinate values in global vector */
234 Vec coordinatesLocal; /* Coordinate values in local vector */
235 PetscBool periodic; /* Is the DM periodic? */
236 DMField coordinateField; /* Coordinates as an abstract field */
237 PetscReal *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
238 DMBoundaryType *bdtype; /* Indicates type of topological boundary */
239 /* Null spaces -- of course I should make this have a variable number of fields */
240 /* I now believe this might not be the right way: see below */
241 NullSpaceFunc nullspaceConstructors[10];
242 NullSpaceFunc nearnullspaceConstructors[10];
243 /* Fields are represented by objects */
244 PetscInt Nf; /* Number of fields defined on the total domain */
245 RegionField *fields; /* Array of discretization fields with regions of validity */
246 DMBoundary boundary; /* List of boundary conditions */
247 PetscInt Nds; /* Number of discrete systems defined on the total domain */
248 DMSpace *probs; /* Array of discrete systems */
249 /* Output structures */
250 DM dmBC; /* The DM with boundary conditions in the global DM */
251 PetscInt outputSequenceNum; /* The current sequence number for output */
252 PetscReal outputSequenceVal; /* The current sequence value for output */
253 PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
254 PetscErrorCode (*monitordestroy[MAXDMMONITORS])(void **);
255 void *monitorcontext[MAXDMMONITORS];
256 PetscInt numbermonitors;
257
258 PetscObject dmksp,dmsnes,dmts;
259 };
260
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;
269 PETSC_EXTERN PetscLogEvent DM_CreateInjection;
270 PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
271 PETSC_EXTERN PetscLogEvent DM_Load;
272 PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
273
274 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM,Vec*);
275 PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM,Vec*);
276
277 PETSC_EXTERN PetscErrorCode DMView_GLVis(DM,PetscViewer,PetscErrorCode(*)(DM,PetscViewer));
278
279 /*
280
281 Composite Vectors
282
283 Single global representation
284 Individual global representations
285 Single local representation
286 Individual local representations
287
288 Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
289
290 DM da_u, da_v, da_p
291
292 DM dm_velocities
293
294 DM dm
295
296 DMDACreate(,&da_u);
297 DMDACreate(,&da_v);
298 DMCompositeCreate(,&dm_velocities);
299 DMCompositeAddDM(dm_velocities,(DM)du);
300 DMCompositeAddDM(dm_velocities,(DM)dv);
301
302 DMDACreate(,&da_p);
303 DMCompositeCreate(,&dm_velocities);
304 DMCompositeAddDM(dm,(DM)dm_velocities);
305 DMCompositeAddDM(dm,(DM)dm_p);
306
307
308 Access parts of composite vectors (Composite only)
309 ---------------------------------
310 DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
311 ADD for local vector -
312
313 Element access
314 --------------
315 From global vectors
316 -DAVecGetArray - for DMDA
317 -VecGetArray - for DMSliced
318 ADD for DMComposite??? maybe
319
320 From individual vector
321 -DAVecGetArray - for DMDA
322 -VecGetArray -for sliced
323 ADD for DMComposite??? maybe
324
325 From single local vector
326 ADD * single local vector as arrays?
327
328 Communication
329 -------------
330 DMGlobalToLocal - global vector to single local vector
331
332 DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
333
334 Obtaining vectors
335 -----------------
336 DMCreateGlobal/Local
337 DMGetGlobal/Local
338 DMCompositeGetLocalVectors - gives individual local work vectors and arrays
339
340
341 ????? individual global vectors ????
342
343 */
344
345 #if defined(PETSC_HAVE_HDF5)
346 PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char *, PetscInt, PetscScalar *, PetscViewer);
347 #endif
348
DMGetLocalOffset_Private(DM dm,PetscInt point,PetscInt * start,PetscInt * end)349 PETSC_STATIC_INLINE PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
350 {
351 PetscFunctionBeginHot;
352 #if defined(PETSC_USE_DEBUG)
353 {
354 PetscInt dof;
355 PetscErrorCode ierr;
356 *start = *end = 0; /* Silence overzealous compiler warning */
357 if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
358 ierr = PetscSectionGetOffset(dm->localSection, point, start);CHKERRQ(ierr);
359 ierr = PetscSectionGetDof(dm->localSection, point, &dof);CHKERRQ(ierr);
360 *end = *start + dof;
361 }
362 #else
363 {
364 const PetscSection s = dm->localSection;
365 *start = s->atlasOff[point - s->pStart];
366 *end = *start + s->atlasDof[point - s->pStart];
367 }
368 #endif
369 PetscFunctionReturn(0);
370 }
371
DMGetLocalFieldOffset_Private(DM dm,PetscInt point,PetscInt field,PetscInt * start,PetscInt * end)372 PETSC_STATIC_INLINE PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
373 {
374 PetscFunctionBegin;
375 #if defined(PETSC_USE_DEBUG)
376 {
377 PetscInt dof;
378 PetscErrorCode ierr;
379 *start = *end = 0; /* Silence overzealous compiler warning */
380 if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
381 ierr = PetscSectionGetFieldOffset(dm->localSection, point, field, start);CHKERRQ(ierr);
382 ierr = PetscSectionGetFieldDof(dm->localSection, point, field, &dof);CHKERRQ(ierr);
383 *end = *start + dof;
384 }
385 #else
386 {
387 const PetscSection s = dm->localSection->field[field];
388 *start = s->atlasOff[point - s->pStart];
389 *end = *start + s->atlasDof[point - s->pStart];
390 }
391 #endif
392 PetscFunctionReturn(0);
393 }
394
DMGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt * start,PetscInt * end)395 PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
396 {
397 PetscFunctionBegin;
398 #if defined(PETSC_USE_DEBUG)
399 {
400 PetscErrorCode ierr;
401 PetscInt dof,cdof;
402 *start = *end = 0; /* Silence overzealous compiler warning */
403 if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
404 if (!dm->globalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
405 ierr = PetscSectionGetOffset(dm->globalSection, point, start);CHKERRQ(ierr);
406 ierr = PetscSectionGetDof(dm->globalSection, point, &dof);CHKERRQ(ierr);
407 ierr = PetscSectionGetConstraintDof(dm->globalSection, point, &cdof);CHKERRQ(ierr);
408 *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
409 }
410 #else
411 {
412 const PetscSection s = dm->globalSection;
413 const PetscInt dof = s->atlasDof[point - s->pStart];
414 const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
415 *start = s->atlasOff[point - s->pStart];
416 *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
417 }
418 #endif
419 PetscFunctionReturn(0);
420 }
421
DMGetGlobalFieldOffset_Private(DM dm,PetscInt point,PetscInt field,PetscInt * start,PetscInt * end)422 PETSC_STATIC_INLINE PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
423 {
424 PetscFunctionBegin;
425 #if defined(PETSC_USE_DEBUG)
426 {
427 PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
428 PetscErrorCode ierr;
429 *start = *end = 0; /* Silence overzealous compiler warning */
430 if (!dm->localSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
431 if (!dm->globalSection) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be crated automatically by DMGetGlobalSection()");
432 ierr = PetscSectionGetOffset(dm->globalSection, point, start);CHKERRQ(ierr);
433 ierr = PetscSectionGetOffset(dm->localSection, point, &loff);CHKERRQ(ierr);
434 ierr = PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff);CHKERRQ(ierr);
435 ierr = PetscSectionGetFieldDof(dm->localSection, point, field, &fdof);CHKERRQ(ierr);
436 ierr = PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof);CHKERRQ(ierr);
437 *start = *start < 0 ? *start - (lfoff-loff) : *start + lfoff-loff;
438 for (f = 0; f < field; ++f) {
439 ierr = PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof);CHKERRQ(ierr);
440 *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
441 }
442 *end = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
443 }
444 #else
445 {
446 const PetscSection s = dm->localSection;
447 const PetscSection fs = dm->localSection->field[field];
448 const PetscSection gs = dm->globalSection;
449 const PetscInt loff = s->atlasOff[point - s->pStart];
450 const PetscInt goff = gs->atlasOff[point - s->pStart];
451 const PetscInt lfoff = fs->atlasOff[point - s->pStart];
452 const PetscInt fdof = fs->atlasDof[point - s->pStart];
453 const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
454 PetscInt ffcdof = 0, f;
455
456 for (f = 0; f < field; ++f) {
457 const PetscSection ffs = dm->localSection->field[f];
458 ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
459 }
460 *start = goff + (goff < 0 ? loff-lfoff + ffcdof : lfoff-loff - ffcdof);
461 *end = *start < 0 ? *start - (fdof-fcdof) : *start + fdof-fcdof;
462 }
463 #endif
464 PetscFunctionReturn(0);
465 }
466
467 PETSC_EXTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
468 PETSC_EXTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
469 PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);
470
471 PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
472 PETSC_INTERN PetscErrorCode DMSetField_Internal(DM, PetscInt, DMLabel, PetscObject);
473
474 #endif
475