1 #pragma once
2
3 #include <petscdm.h>
4 #ifdef PETSC_HAVE_LIBCEED
5 #include <petscdmceed.h>
6 #endif
7 #include <petsc/private/petscimpl.h>
8 #include <petsc/private/petscdsimpl.h>
9 #include <petsc/private/sectionimpl.h> /* for inline access to atlasOff */
10
11 PETSC_EXTERN PetscBool DMRegisterAllCalled;
12 PETSC_EXTERN PetscErrorCode DMRegisterAll(void);
13
14 typedef struct _PetscHashAuxKey {
15 DMLabel label;
16 PetscInt value;
17 PetscInt part;
18 } PetscHashAuxKey;
19
20 #define PetscHashAuxKeyHash(key) PetscHashCombine(PetscHashCombine(PetscHashPointer((key).label), PetscHashInt((key).value)), PetscHashInt((key).part))
21
22 #define PetscHashAuxKeyEqual(k1, k2) (((k1).label == (k2).label) ? (((k1).value == (k2).value) ? ((k1).part == (k2).part) : 0) : 0)
23
24 PETSC_HASH_MAP(HMapAux, PetscHashAuxKey, Vec, PetscHashAuxKeyHash, PetscHashAuxKeyEqual, PETSC_NULLPTR)
25
26 struct _n_DMGeneratorFunctionList {
27 PetscErrorCode (*generate)(DM, PetscBool, DM *);
28 PetscErrorCode (*refine)(DM, PetscReal *, DM *);
29 PetscErrorCode (*adapt)(DM, Vec, DMLabel, DMLabel, DM *);
30 char *name;
31 PetscInt dim;
32 DMGeneratorFunctionList next;
33 };
34
35 typedef struct _DMOps *DMOps;
36 struct _DMOps {
37 PetscErrorCode (*view)(DM, PetscViewer);
38 PetscErrorCode (*load)(DM, PetscViewer);
39 PetscErrorCode (*clone)(DM, DM *);
40 PetscErrorCode (*setfromoptions)(DM, PetscOptionItems);
41 PetscErrorCode (*setup)(DM);
42 PetscErrorCode (*createlocalsection)(DM);
43 PetscErrorCode (*createsectionpermutation)(DM, IS *, PetscBT *);
44 PetscErrorCode (*createdefaultconstraints)(DM);
45 PetscErrorCode (*createglobalvector)(DM, Vec *);
46 PetscErrorCode (*createlocalvector)(DM, Vec *);
47 PetscErrorCode (*getlocaltoglobalmapping)(DM);
48 PetscErrorCode (*createfieldis)(DM, PetscInt *, char ***, IS **);
49 PetscErrorCode (*createcoordinatedm)(DM, DM *);
50 PetscErrorCode (*createcellcoordinatedm)(DM, DM *);
51 PetscErrorCode (*createcoordinatefield)(DM, DMField *);
52
53 PetscErrorCode (*getcoloring)(DM, ISColoringType, ISColoring *);
54 PetscErrorCode (*creatematrix)(DM, Mat *);
55 PetscErrorCode (*createinterpolation)(DM, DM, Mat *, Vec *);
56 PetscErrorCode (*createrestriction)(DM, DM, Mat *);
57 PetscErrorCode (*createmassmatrix)(DM, DM, Mat *);
58 PetscErrorCode (*createmassmatrixlumped)(DM, Vec *, Vec *);
59 PetscErrorCode (*creategradientmatrix)(DM, DM, Mat *);
60 PetscErrorCode (*hascreateinjection)(DM, PetscBool *);
61 PetscErrorCode (*createinjection)(DM, DM, Mat *);
62
63 PetscErrorCode (*refine)(DM, MPI_Comm, DM *);
64 PetscErrorCode (*coarsen)(DM, MPI_Comm, DM *);
65 PetscErrorCode (*refinehierarchy)(DM, PetscInt, DM *);
66 PetscErrorCode (*coarsenhierarchy)(DM, PetscInt, DM *);
67 PetscErrorCode (*extrude)(DM, PetscInt, DM *);
68
69 PetscErrorCode (*globaltolocalbegin)(DM, Vec, InsertMode, Vec);
70 PetscErrorCode (*globaltolocalend)(DM, Vec, InsertMode, Vec);
71 PetscErrorCode (*localtoglobalbegin)(DM, Vec, InsertMode, Vec);
72 PetscErrorCode (*localtoglobalend)(DM, Vec, InsertMode, Vec);
73 PetscErrorCode (*localtolocalbegin)(DM, Vec, InsertMode, Vec);
74 PetscErrorCode (*localtolocalend)(DM, Vec, InsertMode, Vec);
75
76 PetscErrorCode (*destroy)(DM);
77
78 PetscErrorCode (*computevariablebounds)(DM, Vec, Vec);
79
80 PetscErrorCode (*createsubdm)(DM, PetscInt, const PetscInt *, IS *, DM *);
81 PetscErrorCode (*createsuperdm)(DM *, PetscInt, IS **, DM *);
82 PetscErrorCode (*createfielddecomposition)(DM, PetscInt *, char ***, IS **, DM **);
83 PetscErrorCode (*createdomaindecomposition)(DM, PetscInt *, char ***, IS **, IS **, DM **);
84 PetscErrorCode (*createddscatters)(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
85
86 PetscErrorCode (*getdimpoints)(DM, PetscInt, PetscInt *, PetscInt *);
87 PetscErrorCode (*locatepoints)(DM, Vec, DMPointLocationType, PetscSF);
88 PetscErrorCode (*getneighbors)(DM, PetscInt *, const PetscMPIInt **);
89 PetscErrorCode (*getboundingbox)(DM, PetscReal *, PetscReal *);
90 PetscErrorCode (*getlocalboundingbox)(DM, PetscReal[], PetscReal[], PetscInt[], PetscInt[]);
91 PetscErrorCode (*locatepointssubdomain)(DM, Vec, PetscMPIInt **);
92 PetscErrorCode (*snaptogeommodel)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
93
94 PetscErrorCode (*projectfunctionlocal)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
95 PetscErrorCode (*projectfunctionlabellocal)(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
96 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);
97 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);
98 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);
99 PetscErrorCode (*computel2diff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
100 PetscErrorCode (*computel2gradientdiff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
101 PetscErrorCode (*computel2fielddiff)(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
102
103 PetscErrorCode (*getcompatibility)(DM, DM, PetscBool *, PetscBool *);
104 };
105
106 PETSC_INTERN PetscErrorCode DMLocalizeCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
107 PETSC_INTERN PetscErrorCode DMLocalizeCoordinateReal_Internal(DM, PetscInt, const PetscReal[], const PetscReal[], PetscReal[]);
108 PETSC_INTERN PetscErrorCode DMLocalizeAddCoordinate_Internal(DM, PetscInt, const PetscScalar[], const PetscScalar[], PetscScalar[]);
109
110 PETSC_INTERN PetscErrorCode DMCountNonCyclicReferences(PetscObject, PetscInt *);
111
112 typedef struct _DMCoarsenHookLink *DMCoarsenHookLink;
113 struct _DMCoarsenHookLink {
114 PetscErrorCode (*coarsenhook)(DM, DM, void *); /* Run once, when coarse DM is created */
115 PetscErrorCode (*restricthook)(DM, Mat, Vec, Mat, DM, void *); /* Run each time a new problem is restricted to a coarse grid */
116 void *ctx;
117 DMCoarsenHookLink next;
118 };
119
120 typedef struct _DMRefineHookLink *DMRefineHookLink;
121 struct _DMRefineHookLink {
122 PetscErrorCode (*refinehook)(DM, DM, void *); /* Run once, when a fine DM is created */
123 PetscErrorCode (*interphook)(DM, Mat, DM, void *); /* Run each time a new problem is interpolated to a fine grid */
124 void *ctx;
125 DMRefineHookLink next;
126 };
127
128 typedef struct _DMSubDomainHookLink *DMSubDomainHookLink;
129 struct _DMSubDomainHookLink {
130 PetscErrorCode (*ddhook)(DM, DM, void *);
131 PetscErrorCode (*restricthook)(DM, VecScatter, VecScatter, DM, void *);
132 void *ctx;
133 DMSubDomainHookLink next;
134 };
135
136 typedef struct _DMGlobalToLocalHookLink *DMGlobalToLocalHookLink;
137 struct _DMGlobalToLocalHookLink {
138 PetscErrorCode (*beginhook)(DM, Vec, InsertMode, Vec, void *);
139 PetscErrorCode (*endhook)(DM, Vec, InsertMode, Vec, void *);
140 void *ctx;
141 DMGlobalToLocalHookLink next;
142 };
143
144 typedef struct _DMLocalToGlobalHookLink *DMLocalToGlobalHookLink;
145 struct _DMLocalToGlobalHookLink {
146 PetscErrorCode (*beginhook)(DM, Vec, InsertMode, Vec, void *);
147 PetscErrorCode (*endhook)(DM, Vec, InsertMode, Vec, void *);
148 void *ctx;
149 DMLocalToGlobalHookLink next;
150 };
151
152 typedef enum {
153 DMVEC_STATUS_IN,
154 DMVEC_STATUS_OUT
155 } DMVecStatus;
156 typedef struct _DMNamedVecLink *DMNamedVecLink;
157 struct _DMNamedVecLink {
158 Vec X;
159 char *name;
160 DMVecStatus status;
161 DMNamedVecLink next;
162 };
163
164 typedef struct _DMWorkLink *DMWorkLink;
165 struct _DMWorkLink {
166 size_t bytes;
167 void *mem;
168 DMWorkLink next;
169 };
170
171 #define DM_MAX_WORK_VECTORS 100 /* work vectors available to users via DMGetGlobalVector(), DMGetLocalVector() */
172
173 struct _n_DMLabelLink {
174 DMLabel label;
175 PetscBool output;
176 struct _n_DMLabelLink *next;
177 };
178 typedef struct _n_DMLabelLink *DMLabelLink;
179
180 typedef struct _n_Boundary *DMBoundary;
181
182 struct _n_Boundary {
183 DSBoundary dsboundary;
184 DMLabel label;
185 DMBoundary next;
186 };
187
188 typedef struct _n_Field {
189 PetscObject disc; /* Field discretization, or a PetscContainer with the field name */
190 DMLabel label; /* Label defining the domain of definition of the field */
191 PetscBool adjacency[2]; /* Flags for defining variable influence (adjacency) for each field [use cone() or support() first, use the transitive closure] */
192 PetscBool avoidTensor; /* Flag to avoid defining field over tensor cells */
193 } RegionField;
194
195 typedef struct _n_Space {
196 PetscDS ds; /* Approximation space in this domain */
197 PetscDS dsIn; /* Approximation space for input to this domain (now only used for cohesive cells) */
198 DMLabel label; /* Label defining the domain of definition of the discretization */
199 IS fields; /* Map from DS field numbers to original field numbers in the DM */
200 } DMSpace;
201
202 struct _p_UniversalLabel {
203 DMLabel label; /* The universal label */
204 PetscInt Nl; /* Number of labels encoded */
205 char **names; /* The label names */
206 PetscInt *indices; /* The original indices in the input DM */
207 PetscInt Nv; /* Total number of values in all the labels */
208 PetscInt *bits; /* Starting bit for values of each label */
209 PetscInt *masks; /* Masks to pull out label value bits for each label */
210 PetscInt *offsets; /* Starting offset for label values for each label */
211 PetscInt *values; /* Original label values before renumbering */
212 };
213
214 PETSC_INTERN PetscErrorCode DMDestroyLabelLinkList_Internal(DM);
215
216 #define MAXDMMONITORS 5
217
218 typedef struct {
219 PetscInt dim; /* The dimension of the embedding space */
220 DM dm; /* Layout for coordinates */
221 Vec x; /* Coordinate values in global vector */
222 Vec xl; /* Coordinate values in local vector */
223 DMField field; /* Coordinates as an abstract field */
224 } DMCoordinates;
225
226 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*NullSpaceFn)(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullSpace);
227
228 struct _p_DM {
229 PETSCHEADER(struct _DMOps);
230 Vec localin[DM_MAX_WORK_VECTORS], localout[DM_MAX_WORK_VECTORS];
231 Vec globalin[DM_MAX_WORK_VECTORS], globalout[DM_MAX_WORK_VECTORS];
232 DMNamedVecLink namedglobal;
233 DMNamedVecLink namedlocal;
234 DMWorkLink workin, workout;
235 DMLabelLink labels; /* Linked list of labels */
236 DMLabel depthLabel; /* Optimized access to depth label */
237 DMLabel celltypeLabel; /* Optimized access to celltype label */
238 void *ctx; /* a user context */
239 PetscCtxDestroyFn *ctxdestroy;
240 ISColoringType coloringtype;
241 MatFDColoring fd;
242 VecType vectype; /* type of vector created with DMCreateLocalVector() and DMCreateGlobalVector() */
243 MatType mattype; /* type of matrix created with DMCreateMatrix() */
244 PetscInt bind_below; /* Local size threshold (in entries/rows) below which Vec/Mat objects are bound to CPU */
245 PetscInt bs;
246 DMBlockingType blocking_type;
247 ISLocalToGlobalMapping ltogmap;
248 PetscBool prealloc_skip; // Flag indicating the DMCreateMatrix() should not preallocate (only set sizes and local-to-global)
249 PetscBool prealloc_only; /* Flag indicating the DMCreateMatrix() should only preallocate, not fill the matrix */
250 PetscBool structure_only; /* Flag indicating the DMCreateMatrix() create matrix nonzero structure without values */
251 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 */
252 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 */
253 PetscBool setfromoptionscalled;
254 void *data;
255 /* Hierarchy / Submeshes */
256 DM coarseMesh;
257 DM fineMesh;
258 DMCoarsenHookLink coarsenhook; /* For transferring auxiliary problem data to coarser grids */
259 DMRefineHookLink refinehook;
260 DMSubDomainHookLink subdomainhook;
261 DMGlobalToLocalHookLink gtolhook;
262 DMLocalToGlobalHookLink ltoghook;
263 /* Topology */
264 PetscInt dim; /* The topological dimension */
265 /* Auxiliary data */
266 PetscHMapAux auxData; /* Auxiliary DM and Vec for region denoted by the key */
267 /* Flexible communication */
268 PetscSF sfMigration; /* SF for point distribution created during distribution */
269 PetscSF sf; /* SF for parallel point overlap */
270 PetscSF sectionSF; /* SF for parallel dof overlap using section */
271 PetscSF sfNatural; /* SF mapping to the "natural" ordering */
272 PetscBool useNatural; /* Create the natural SF */
273 /* Allows a non-standard data layout */
274 PetscBool adjacency[2]; /* [use cone() or support() first, use the transitive closure] */
275 PetscSection localSection; /* Layout for local vectors */
276 PetscSection globalSection; /* Layout for global vectors */
277 PetscLayout map; /* Parallel division of unknowns across processes */
278 DMReorderDefaultFlag reorderSection; /* Reorder the local section by default */
279 MatOrderingType reorderSectionType; /* The type of reordering */
280
281 // Affine transform applied in DMGlobalToLocal
282 struct {
283 PetscInt num_affines;
284 VecScatter *affine_to_local;
285 Vec *affine;
286 PetscErrorCode (*setup)(DM);
287 } periodic;
288 /* Constraints */
289 struct {
290 PetscSection section;
291 Mat mat;
292 Vec bias;
293 } defaultConstraint;
294 /* Basis transformation */
295 DM transformDM; /* Layout for basis transformation */
296 Vec transform; /* Basis transformation matrices */
297 void *transformCtx; /* Basis transformation context */
298 PetscErrorCode (*transformSetUp)(DM, void *);
299 PetscErrorCode (*transformDestroy)(DM, void *);
300 PetscErrorCode (*transformGetMatrix)(DM, const PetscReal[], PetscBool, const PetscScalar **, void *);
301 /* Coordinates */
302 DMCoordinates coordinates[2]; /* Coordinates, 0 is default, 1 is possible DG coordinate field for periodicity */
303 /* Periodicity */
304 PetscReal *Lstart, *L, *maxCell; /* Size of periodic box and max cell size for determining periodicity */
305 PetscBool sparseLocalize; /* Localize coordinates only for cells near periodic boundary */
306 /* Null spaces */
307 NullSpaceFn *nullspaceConstructors;
308 NullSpaceFn *nearnullspaceConstructors;
309 /* Fields are represented by objects */
310 PetscInt Nf; /* Number of fields defined on the total domain */
311 RegionField *fields; /* Array of discretization fields with regions of validity */
312 DMBoundary boundary; /* List of boundary conditions */
313 PetscInt Nds; /* Number of discrete systems defined on the total domain */
314 DMSpace *probs; /* Array of discrete systems */
315 /* Output structures */
316 DM dmBC; /* The DM with boundary conditions in the global DM */
317 PetscBool ignorePermOutput; /* Ignore the local section permutation on output */
318 PetscInt outputSequenceNum; /* The current sequence number for output */
319 PetscReal outputSequenceVal; /* The current sequence value for output */
320 PetscErrorCode (*monitor[MAXDMMONITORS])(DM, void *);
321 PetscCtxDestroyFn *monitordestroy[MAXDMMONITORS];
322 void *monitorcontext[MAXDMMONITORS];
323 PetscInt numbermonitors;
324 /* Configuration */
325 PetscBool cloneOpts; /* Flag indicating that this is a linked clone and should not respond to some options. This is currently used to prevent transformations from also affecting the coordinate DM */
326
327 PetscObject dmksp, dmsnes, dmts;
328 #ifdef PETSC_HAVE_LIBCEED
329 Ceed ceed; // LibCEED context
330 CeedElemRestriction ceedERestrict; // Map from the local vector (Lvector) to the cells (Evector)
331 #endif
332 DMCeed dmceed; // CEED operator and data for this problem
333 };
334
335 PETSC_EXTERN PetscLogEvent DM_Convert;
336 PETSC_EXTERN PetscLogEvent DM_GlobalToLocal;
337 PETSC_EXTERN PetscLogEvent DM_LocalToGlobal;
338 PETSC_EXTERN PetscLogEvent DM_LocatePoints;
339 PETSC_EXTERN PetscLogEvent DM_Coarsen;
340 PETSC_EXTERN PetscLogEvent DM_Refine;
341 PETSC_EXTERN PetscLogEvent DM_CreateInterpolation;
342 PETSC_EXTERN PetscLogEvent DM_CreateRestriction;
343 PETSC_EXTERN PetscLogEvent DM_CreateInjection;
344 PETSC_EXTERN PetscLogEvent DM_CreateMatrix;
345 PETSC_EXTERN PetscLogEvent DM_CreateMassMatrix;
346 PETSC_EXTERN PetscLogEvent DM_Load;
347 PETSC_EXTERN PetscLogEvent DM_View;
348 PETSC_EXTERN PetscLogEvent DM_AdaptInterpolator;
349 PETSC_EXTERN PetscLogEvent DM_ProjectFunction;
350
351 PETSC_INTERN PetscErrorCode DMCreateGlobalVector_Section_Private(DM, Vec *);
352 PETSC_INTERN PetscErrorCode DMCreateLocalVector_Section_Private(DM, Vec *);
353
354 PETSC_INTERN PetscErrorCode DMView_GLVis(DM, PetscViewer, PetscErrorCode (*)(DM, PetscViewer));
355
356 /*
357
358 Composite Vectors
359
360 Single global representation
361 Individual global representations
362 Single local representation
363 Individual local representations
364
365 Subsets of individual as a single????? Do we handle this by having DMComposite inside composite??????
366
367 DM da_u, da_v, da_p
368
369 DM dm_velocities
370
371 DM dm
372
373 DMDACreate(,&da_u);
374 DMDACreate(,&da_v);
375 DMCompositeCreate(,&dm_velocities);
376 DMCompositeAddDM(dm_velocities,(DM)du);
377 DMCompositeAddDM(dm_velocities,(DM)dv);
378
379 DMDACreate(,&da_p);
380 DMCompositeCreate(,&dm_velocities);
381 DMCompositeAddDM(dm,(DM)dm_velocities);
382 DMCompositeAddDM(dm,(DM)dm_p);
383
384 Access parts of composite vectors (Composite only)
385 ---------------------------------
386 DMCompositeGetAccess - access the global vector as subvectors and array (for redundant arrays)
387 ADD for local vector -
388
389 Element access
390 --------------
391 From global vectors
392 -DAVecGetArray - for DMDA
393 -VecGetArray - for DMSliced
394 ADD for DMComposite??? maybe
395
396 From individual vector
397 -DAVecGetArray - for DMDA
398 -VecGetArray -for sliced
399 ADD for DMComposite??? maybe
400
401 From single local vector
402 ADD * single local vector as arrays?
403
404 Communication
405 -------------
406 DMGlobalToLocal - global vector to single local vector
407
408 DMCompositeScatter/Gather - direct to individual local vectors and arrays CHANGE name closer to GlobalToLocal?
409
410 Obtaining vectors
411 -----------------
412 DMCreateGlobal/Local
413 DMGetGlobal/Local
414 DMCompositeGetLocalVectors - gives individual local work vectors and arrays
415
416 ????? individual global vectors ????
417
418 */
419
420 #if defined(PETSC_HAVE_HDF5)
421 PETSC_EXTERN PetscErrorCode DMSequenceLoad_HDF5_Internal(DM, const char[], PetscInt, PetscScalar *, PetscViewer);
422 PETSC_EXTERN PetscErrorCode DMSequenceGetLength_HDF5_Internal(DM, const char[], PetscInt *, PetscViewer);
423 #endif
424
DMGetLocalOffset_Private(DM dm,PetscInt point,PetscInt * start,PetscInt * end)425 static inline PetscErrorCode DMGetLocalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
426 {
427 PetscFunctionBeginHot;
428 #if defined(PETSC_USE_DEBUG)
429 {
430 PetscInt dof;
431
432 *start = *end = 0; /* Silence overzealous compiler warning */
433 PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
434 PetscCall(PetscSectionGetOffset(dm->localSection, point, start));
435 PetscCall(PetscSectionGetDof(dm->localSection, point, &dof));
436 *end = *start + dof;
437 }
438 #else
439 {
440 const PetscSection s = dm->localSection;
441 *start = s->atlasOff[point - s->pStart];
442 *end = *start + s->atlasDof[point - s->pStart];
443 }
444 #endif
445 PetscFunctionReturn(PETSC_SUCCESS);
446 }
447
DMGetLocalFieldOffset_Private(DM dm,PetscInt point,PetscInt field,PetscInt * start,PetscInt * end)448 static inline PetscErrorCode DMGetLocalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
449 {
450 PetscFunctionBegin;
451 #if defined(PETSC_USE_DEBUG)
452 {
453 PetscInt dof;
454 *start = *end = 0; /* Silence overzealous compiler warning */
455 PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
456 PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, start));
457 PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &dof));
458 *end = *start + dof;
459 }
460 #else
461 {
462 const PetscSection s = dm->localSection->field[field];
463 *start = s->atlasOff[point - s->pStart];
464 *end = *start + s->atlasDof[point - s->pStart];
465 }
466 #endif
467 PetscFunctionReturn(PETSC_SUCCESS);
468 }
469
DMGetGlobalOffset_Private(DM dm,PetscInt point,PetscInt * start,PetscInt * end)470 static inline PetscErrorCode DMGetGlobalOffset_Private(DM dm, PetscInt point, PetscInt *start, PetscInt *end)
471 {
472 PetscFunctionBegin;
473 #if defined(PETSC_USE_DEBUG)
474 {
475 PetscInt dof, cdof;
476 *start = *end = 0; /* Silence overzealous compiler warning */
477 PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
478 PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
479 PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
480 PetscCall(PetscSectionGetDof(dm->globalSection, point, &dof));
481 PetscCall(PetscSectionGetConstraintDof(dm->globalSection, point, &cdof));
482 *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
483 }
484 #else
485 {
486 const PetscSection s = dm->globalSection;
487 const PetscInt dof = s->atlasDof[point - s->pStart];
488 const PetscInt cdof = s->bc ? s->bc->atlasDof[point - s->bc->pStart] : 0;
489 *start = s->atlasOff[point - s->pStart];
490 *end = *start + dof - cdof + (dof < 0 ? 1 : 0);
491 }
492 #endif
493 PetscFunctionReturn(PETSC_SUCCESS);
494 }
495
DMGetGlobalFieldOffset_Private(DM dm,PetscInt point,PetscInt field,PetscInt * start,PetscInt * end)496 static inline PetscErrorCode DMGetGlobalFieldOffset_Private(DM dm, PetscInt point, PetscInt field, PetscInt *start, PetscInt *end)
497 {
498 PetscFunctionBegin;
499 #if defined(PETSC_USE_DEBUG)
500 {
501 PetscInt loff, lfoff, fdof, fcdof, ffcdof, f;
502 *start = *end = 0; /* Silence overzealous compiler warning */
503 PetscCheck(dm->localSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a local section, see DMSetLocalSection()");
504 PetscCheck(dm->globalSection, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "DM must have a global section. It will be created automatically by DMGetGlobalSection()");
505 PetscCall(PetscSectionGetOffset(dm->globalSection, point, start));
506 PetscCall(PetscSectionGetOffset(dm->localSection, point, &loff));
507 PetscCall(PetscSectionGetFieldOffset(dm->localSection, point, field, &lfoff));
508 PetscCall(PetscSectionGetFieldDof(dm->localSection, point, field, &fdof));
509 PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, field, &fcdof));
510 *start = *start < 0 ? *start - (lfoff - loff) : *start + lfoff - loff;
511 for (f = 0; f < field; ++f) {
512 PetscCall(PetscSectionGetFieldConstraintDof(dm->localSection, point, f, &ffcdof));
513 *start = *start < 0 ? *start + ffcdof : *start - ffcdof;
514 }
515 *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
516 }
517 #else
518 {
519 const PetscSection s = dm->localSection;
520 const PetscSection fs = dm->localSection->field[field];
521 const PetscSection gs = dm->globalSection;
522 const PetscInt loff = s->atlasOff[point - s->pStart];
523 const PetscInt goff = gs->atlasOff[point - s->pStart];
524 const PetscInt lfoff = fs->atlasOff[point - s->pStart];
525 const PetscInt fdof = fs->atlasDof[point - s->pStart];
526 const PetscInt fcdof = fs->bc ? fs->bc->atlasDof[point - fs->bc->pStart] : 0;
527 PetscInt ffcdof = 0, f;
528
529 for (f = 0; f < field; ++f) {
530 const PetscSection ffs = dm->localSection->field[f];
531 ffcdof += ffs->bc ? ffs->bc->atlasDof[point - ffs->bc->pStart] : 0;
532 }
533 *start = goff + (goff < 0 ? loff - lfoff + ffcdof : lfoff - loff - ffcdof);
534 *end = *start < 0 ? *start - (fdof - fcdof) : *start + fdof - fcdof;
535 }
536 #endif
537 PetscFunctionReturn(PETSC_SUCCESS);
538 }
539
540 PETSC_INTERN PetscErrorCode DMGetCoordinateDegree_Internal(DM, PetscInt *);
541 PETSC_INTERN PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM, PetscReal[], PetscReal[], PetscInt[], PetscInt[]);
542 PETSC_INTERN PetscErrorCode DMGetIsoperiodicPointSF_Internal(DM dm, PetscSF *sf);
543
544 PETSC_INTERN PetscErrorCode DMGetBasisTransformDM_Internal(DM, DM *);
545 PETSC_INTERN PetscErrorCode DMGetBasisTransformVec_Internal(DM, Vec *);
546 PETSC_INTERN PetscErrorCode DMConstructBasisTransform_Internal(DM);
547
548 PETSC_INTERN PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM, PetscReal[], PetscReal[]);
549 PETSC_INTERN PetscErrorCode DMSetField_Internal(DM, PetscInt, DMLabel, PetscObject);
550
551 PETSC_INTERN PetscErrorCode DMSetLabelValue_Fast(DM, DMLabel *, const char[], PetscInt, PetscInt);
552 PETSC_INTERN PetscErrorCode DMGetPoints_Internal(DM, DMLabel, PetscInt, PetscInt, IS *);
553
554 PETSC_INTERN PetscErrorCode DMCompleteBCLabels_Internal(DM dm);
555 PETSC_EXTERN PetscErrorCode DMUniversalLabelCreate(DM, DMUniversalLabel *);
556 PETSC_EXTERN PetscErrorCode DMUniversalLabelDestroy(DMUniversalLabel *);
557 PETSC_EXTERN PetscErrorCode DMUniversalLabelGetLabel(DMUniversalLabel, DMLabel *);
558 PETSC_EXTERN PetscErrorCode DMUniversalLabelCreateLabels(DMUniversalLabel, PetscBool, DM);
559 PETSC_EXTERN PetscErrorCode DMUniversalLabelSetLabelValue(DMUniversalLabel, DM, PetscBool, PetscInt, PetscInt);
560