1 #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
2 #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/
3 #include <petscdmplex.h>
4
5 const char *const DMFieldContinuities[] = {"VERTEX", "EDGE", "FACET", "CELL", NULL};
6
DMFieldCreate(DM dm,PetscInt numComponents,DMFieldContinuity continuity,DMField * field)7 PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm, PetscInt numComponents, DMFieldContinuity continuity, DMField *field)
8 {
9 DMField b;
10
11 PetscFunctionBegin;
12 PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13 PetscAssertPointer(field, 4);
14 PetscCall(DMFieldInitializePackage());
15
16 PetscCall(PetscHeaderCreate(b, DMFIELD_CLASSID, "DMField", "Field over DM", "DM", PetscObjectComm((PetscObject)dm), DMFieldDestroy, DMFieldView));
17 PetscCall(PetscObjectReference((PetscObject)dm));
18 b->dm = dm;
19 b->continuity = continuity;
20 b->numComponents = numComponents;
21 *field = b;
22 PetscFunctionReturn(PETSC_SUCCESS);
23 }
24
25 /*@
26 DMFieldDestroy - destroy a `DMField`
27
28 Collective
29
30 Input Parameter:
31 . field - address of `DMField`
32
33 Level: advanced
34
35 .seealso: `DMField`, `DMFieldCreate()`
36 @*/
DMFieldDestroy(DMField * field)37 PetscErrorCode DMFieldDestroy(DMField *field)
38 {
39 PetscFunctionBegin;
40 if (!*field) PetscFunctionReturn(PETSC_SUCCESS);
41 PetscValidHeaderSpecific(*field, DMFIELD_CLASSID, 1);
42 if (--((PetscObject)*field)->refct > 0) {
43 *field = NULL;
44 PetscFunctionReturn(PETSC_SUCCESS);
45 }
46 PetscTryTypeMethod(*field, destroy);
47 PetscCall(DMDestroy(&(*field)->dm));
48 PetscCall(PetscHeaderDestroy(field));
49 PetscFunctionReturn(PETSC_SUCCESS);
50 }
51
52 /*@
53 DMFieldView - view a `DMField`
54
55 Collective
56
57 Input Parameters:
58 + field - `DMField`
59 - viewer - viewer to display field, for example `PETSC_VIEWER_STDOUT_WORLD`
60
61 Level: advanced
62
63 .seealso: `DMField`, `DMFieldCreate()`
64 @*/
DMFieldView(DMField field,PetscViewer viewer)65 PetscErrorCode DMFieldView(DMField field, PetscViewer viewer)
66 {
67 PetscBool isascii;
68
69 PetscFunctionBegin;
70 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
71 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field), &viewer));
72 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
73 PetscCheckSameComm(field, 1, viewer, 2);
74 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
75 if (isascii) {
76 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)field, viewer));
77 PetscCall(PetscViewerASCIIPushTab(viewer));
78 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " components\n", field->numComponents));
79 PetscCall(PetscViewerASCIIPrintf(viewer, "%s continuity\n", DMFieldContinuities[field->continuity]));
80 PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
81 PetscCall(DMView(field->dm, viewer));
82 PetscCall(PetscViewerPopFormat(viewer));
83 }
84 PetscTryTypeMethod(field, view, viewer);
85 if (isascii) PetscCall(PetscViewerASCIIPopTab(viewer));
86 PetscFunctionReturn(PETSC_SUCCESS);
87 }
88
89 /*@
90 DMFieldSetType - set the `DMField` implementation
91
92 Collective
93
94 Input Parameters:
95 + field - the `DMField` context
96 - type - a known method
97
98 Level: advanced
99
100 Notes:
101 See "include/petscvec.h" for available methods (for instance)
102 + `DMFIELDDA` - a field defined only by its values at the corners of a `DMDA`
103 . `DMFIELDDS` - a field defined by a discretization over a mesh set with `DMSetField()`
104 - `DMFIELDSHELL` - a field defined by arbitrary callbacks
105
106 .seealso: `DMField`, `DMFieldGetType()`, `DMFieldType`,
107 @*/
DMFieldSetType(DMField field,DMFieldType type)108 PetscErrorCode DMFieldSetType(DMField field, DMFieldType type)
109 {
110 PetscBool match;
111 PetscErrorCode (*r)(DMField);
112
113 PetscFunctionBegin;
114 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
115 PetscAssertPointer(type, 2);
116
117 PetscCall(PetscObjectTypeCompare((PetscObject)field, type, &match));
118 if (match) PetscFunctionReturn(PETSC_SUCCESS);
119
120 PetscCall(PetscFunctionListFind(DMFieldList, type, &r));
121 PetscCheck(r, PetscObjectComm((PetscObject)field), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested DMField type %s", type);
122 /* Destroy the previous private DMField context */
123 PetscTryTypeMethod(field, destroy);
124
125 PetscCall(PetscMemzero(field->ops, sizeof(*field->ops)));
126 PetscCall(PetscObjectChangeTypeName((PetscObject)field, type));
127 field->ops->create = r;
128 PetscCall((*r)(field));
129 PetscFunctionReturn(PETSC_SUCCESS);
130 }
131
132 /*@
133 DMFieldGetType - Gets the `DMFieldType` name (as a string) from the `DMField`.
134
135 Not Collective
136
137 Input Parameter:
138 . field - The `DMField` context
139
140 Output Parameter:
141 . type - The `DMFieldType` name
142
143 Level: advanced
144
145 .seealso: `DMField`, `DMFieldSetType()`, `DMFieldType`
146 @*/
DMFieldGetType(DMField field,DMFieldType * type)147 PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type)
148 {
149 PetscFunctionBegin;
150 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
151 PetscAssertPointer(type, 2);
152 PetscCall(DMFieldRegisterAll());
153 *type = ((PetscObject)field)->type_name;
154 PetscFunctionReturn(PETSC_SUCCESS);
155 }
156
157 /*@
158 DMFieldGetNumComponents - Returns the number of components in the field
159
160 Not Collective
161
162 Input Parameter:
163 . field - The `DMField` object
164
165 Output Parameter:
166 . nc - The number of field components
167
168 Level: intermediate
169
170 .seealso: `DMField`, `DMFieldEvaluate()`
171 @*/
DMFieldGetNumComponents(DMField field,PetscInt * nc)172 PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc)
173 {
174 PetscFunctionBegin;
175 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
176 PetscAssertPointer(nc, 2);
177 *nc = field->numComponents;
178 PetscFunctionReturn(PETSC_SUCCESS);
179 }
180
181 /*@
182 DMFieldGetDM - Returns the `DM` for the manifold over which the field is defined.
183
184 Not Collective
185
186 Input Parameter:
187 . field - The `DMField` object
188
189 Output Parameter:
190 . dm - The `DM` object
191
192 Level: intermediate
193
194 .seealso: `DMField`, `DM`, `DMFieldEvaluate()`
195 @*/
DMFieldGetDM(DMField field,DM * dm)196 PetscErrorCode DMFieldGetDM(DMField field, DM *dm)
197 {
198 PetscFunctionBegin;
199 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
200 PetscAssertPointer(dm, 2);
201 *dm = field->dm;
202 PetscFunctionReturn(PETSC_SUCCESS);
203 }
204
205 /*@
206 DMFieldEvaluate - Evaluate the field and its derivatives on a set of points
207
208 Collective
209
210 Input Parameters:
211 + field - The `DMField` object
212 . points - The points at which to evaluate the field. Should have size d x n,
213 where d is the coordinate dimension of the manifold and n is the number
214 of points
215 - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`.
216 If the field is complex and datatype is `PETSC_REAL`, the real part of the
217 field is returned.
218
219 Output Parameters:
220 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
221 If B is not `NULL`, the values of the field are written in this array, varying first by component,
222 then by point.
223 . D - pointer to data of size d * c * n * sizeof(datatype).
224 If `D` is not `NULL`, the values of the field's spatial derivatives are written in this array,
225 varying first by the partial derivative component, then by field component, then by point.
226 - H - pointer to data of size d * d * c * n * sizeof(datatype).
227 If `H` is not `NULL`, the values of the field's second spatial derivatives are written in this array,
228 varying first by the second partial derivative component, then by field component, then by point.
229
230 Level: intermediate
231
232 .seealso: `DMField`, `DMFieldGetDM()`, `DMFieldGetNumComponents()`, `DMFieldEvaluateFE()`, `DMFieldEvaluateFV()`, `PetscDataType`
233 @*/
DMFieldEvaluate(DMField field,Vec points,PetscDataType datatype,void * B,void * D,void * H)234 PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
235 {
236 PetscFunctionBegin;
237 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
238 PetscValidHeaderSpecific(points, VEC_CLASSID, 2);
239 if (B) PetscAssertPointer(B, 4);
240 if (D) PetscAssertPointer(D, 5);
241 if (H) PetscAssertPointer(H, 6);
242 PetscUseTypeMethod(field, evaluate, points, datatype, B, D, H);
243 PetscFunctionReturn(PETSC_SUCCESS);
244 }
245
246 /*@
247 DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from
248 quadrature points on a reference point. The derivatives are taken with respect to the
249 reference coordinates.
250
251 Not Collective
252
253 Input Parameters:
254 + field - The `DMField` object
255 . cellIS - Index set for cells on which to evaluate the field
256 . points - The quadature containing the points in the reference cell at which to evaluate the field.
257 - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`.
258 If the field is complex and datatype is `PETSC_REAL`, the real part of the
259 field is returned.
260
261 Output Parameters:
262 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
263 If B is not `NULL`, the values of the field are written in this array, varying first by component,
264 then by point.
265 . D - pointer to data of size d * c * n * sizeof(datatype).
266 If D is not `NULL`, the values of the field's spatial derivatives are written in this array,
267 varying first by the partial derivative component, then by field component, then by point.
268 - H - pointer to data of size d * d * c * n * sizeof(datatype).
269 If H is not `NULL`, the values of the field's second spatial derivatives are written in this array,
270 varying first by the second partial derivative component, then by field component, then by point.
271
272 Level: intermediate
273
274 .seealso: `DMField`, `DM`, `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFV()`
275 @*/
DMFieldEvaluateFE(DMField field,IS cellIS,PetscQuadrature points,PetscDataType datatype,void * B,void * D,void * H)276 PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
277 {
278 PetscFunctionBegin;
279 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
280 PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2);
281 PetscValidHeader(points, 3);
282 if (B) PetscAssertPointer(B, 5);
283 if (D) PetscAssertPointer(D, 6);
284 if (H) PetscAssertPointer(H, 7);
285 PetscUseTypeMethod(field, evaluateFE, cellIS, points, datatype, B, D, H);
286 PetscFunctionReturn(PETSC_SUCCESS);
287 }
288
289 /*@
290 DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points.
291
292 Not Collective
293
294 Input Parameters:
295 + field - The `DMField` object
296 . cellIS - Index set for cells on which to evaluate the field
297 - datatype - The PetscDataType of the output arrays: either `PETSC_REAL` or `PETSC_SCALAR`.
298 If the field is complex and datatype is `PETSC_REAL`, the real part of the
299 field is returned.
300
301 Output Parameters:
302 + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
303 If B is not `NULL`, the values of the field are written in this array, varying first by component,
304 then by point.
305 . D - pointer to data of size d * c * n * sizeof(datatype).
306 If D is not `NULL`, the values of the field's spatial derivatives are written in this array,
307 varying first by the partial derivative component, then by field component, then by point.
308 - H - pointer to data of size d * d * c * n * sizeof(datatype).
309 If H is not `NULL`, the values of the field's second spatial derivatives are written in this array,
310 varying first by the second partial derivative component, then by field component, then by point.
311
312 Level: intermediate
313
314 .seealso: `DMField`, `IS`, `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFE()`, `PetscDataType`
315 @*/
DMFieldEvaluateFV(DMField field,IS cellIS,PetscDataType datatype,void * B,void * D,void * H)316 PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
317 {
318 PetscFunctionBegin;
319 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
320 PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2);
321 if (B) PetscAssertPointer(B, 4);
322 if (D) PetscAssertPointer(D, 5);
323 if (H) PetscAssertPointer(H, 6);
324 PetscUseTypeMethod(field, evaluateFV, cellIS, datatype, B, D, H);
325 PetscFunctionReturn(PETSC_SUCCESS);
326 }
327
328 /*@
329 DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the
330 reference element
331
332 Not Collective
333
334 Input Parameters:
335 + field - the `DMField` object
336 - cellIS - the index set of points over which we want know the invariance
337
338 Output Parameters:
339 + minDegree - the degree of the largest polynomial space contained in the field on each element
340 - maxDegree - the largest degree of the smallest polynomial space containing the field on any element
341
342 Level: intermediate
343
344 .seealso: `DMField`, `IS`, `DMFieldEvaluateFE()`
345 @*/
DMFieldGetDegree(DMField field,IS cellIS,PeOp PetscInt * minDegree,PeOp PetscInt * maxDegree)346 PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PeOp PetscInt *minDegree, PeOp PetscInt *maxDegree)
347 {
348 PetscFunctionBegin;
349 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
350 PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2);
351 if (minDegree) PetscAssertPointer(minDegree, 3);
352 if (maxDegree) PetscAssertPointer(maxDegree, 4);
353
354 if (minDegree) *minDegree = -1;
355 if (maxDegree) *maxDegree = PETSC_INT_MAX;
356
357 PetscTryTypeMethod(field, getDegree, cellIS, minDegree, maxDegree);
358 PetscFunctionReturn(PETSC_SUCCESS);
359 }
360
361 /*@
362 DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected
363 points via pullback onto the reference element
364
365 Not Collective
366
367 Input Parameters:
368 + field - the `DMField` object
369 - pointIS - the index set of points over which we wish to integrate the field
370
371 Output Parameter:
372 . quad - a `PetscQuadrature` object
373
374 Level: developer
375
376 .seealso: `DMFieldCreateDefaultFaceQuadrature()`, `DMField`, `PetscQuadrature`, `IS`, `DMFieldEvaluteFE()`, `DMFieldGetDegree()`
377 @*/
DMFieldCreateDefaultQuadrature(DMField field,IS pointIS,PetscQuadrature * quad)378 PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad)
379 {
380 PetscFunctionBegin;
381 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
382 PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2);
383 PetscAssertPointer(quad, 3);
384
385 *quad = NULL;
386 PetscTryTypeMethod(field, createDefaultQuadrature, pointIS, quad);
387 PetscFunctionReturn(PETSC_SUCCESS);
388 }
389
390 /*@
391 DMFieldCreateDefaultFaceQuadrature - Creates a quadrature sufficient to integrate the field on all faces of the selected cells via pullback onto the reference element
392
393 Not Collective
394
395 Input Parameters:
396 + field - the `DMField` object
397 - pointIS - the index set of points over which we wish to integrate the field over faces
398
399 Output Parameter:
400 . quad - a `PetscQuadrature` object
401
402 Level: developer
403
404 .seealso: `DMFieldCreateDefaultQuadrature()`, `DMField`, `PetscQuadrature`, `IS`, `DMFieldEvaluteFE()`, `DMFieldGetDegree()`
405 @*/
DMFieldCreateDefaultFaceQuadrature(DMField field,IS pointIS,PetscQuadrature * quad)406 PetscErrorCode DMFieldCreateDefaultFaceQuadrature(DMField field, IS pointIS, PetscQuadrature *quad)
407 {
408 PetscFunctionBegin;
409 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
410 PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2);
411 PetscAssertPointer(quad, 3);
412
413 *quad = NULL;
414 PetscTryTypeMethod(field, createDefaultFaceQuadrature, pointIS, quad);
415 PetscFunctionReturn(PETSC_SUCCESS);
416 }
417
418 /*@C
419 DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field
420
421 Not Collective
422
423 Input Parameters:
424 + field - the `DMField` object
425 . pointIS - the index set of points over which we wish to integrate the field
426 . quad - the quadrature points at which to evaluate the geometric factors
427 - mode - Type of geometry data to store
428
429 Output Parameter:
430 . geom - the geometric factors
431
432 Level: developer
433
434 Note:
435 For some modes, the normal vectors and adjacent cells are calculated
436
437 .seealso: `DMField`, `PetscQuadrature`, `IS`, `PetscFEGeom`, `DMFieldEvaluateFE()`, `DMFieldCreateDefaulteQuadrature()`, `DMFieldGetDegree()`
438 @*/
DMFieldCreateFEGeom(DMField field,IS pointIS,PetscQuadrature quad,PetscFEGeomMode mode,PetscFEGeom ** geom)439 PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeomMode mode, PetscFEGeom **geom)
440 {
441 PetscBool faceData = mode == PETSC_FEGEOM_BOUNDARY || mode == PETSC_FEGEOM_COHESIVE ? PETSC_TRUE : PETSC_FALSE;
442 PetscInt dim, dE;
443 PetscInt nPoints;
444 PetscInt maxDegree;
445 PetscFEGeom *g;
446
447 PetscFunctionBegin;
448 PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
449 PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2);
450 PetscValidHeader(quad, 3);
451 PetscCall(ISGetLocalSize(pointIS, &nPoints));
452 dE = field->numComponents;
453 PetscCall(PetscFEGeomCreate(quad, nPoints, dE, mode, &g));
454 PetscCall(DMFieldEvaluateFE(field, pointIS, quad, PETSC_REAL, g->v, g->J, NULL));
455 dim = g->dim;
456 if (dE > dim) {
457 /* space out J and make square Jacobians */
458 PetscInt i, j, k, N = g->numPoints * g->numCells;
459
460 for (i = N - 1; i >= 0; i--) {
461 PetscReal J[16] = {0};
462
463 for (j = 0; j < dE; j++) {
464 for (k = 0; k < dim; k++) J[j * dE + k] = g->J[i * dE * dim + j * dim + k];
465 }
466 switch (dim) {
467 case 0:
468 for (j = 0; j < dE; j++) {
469 for (k = 0; k < dE; k++) J[j * dE + k] = (j == k) ? 1. : 0.;
470 }
471 break;
472 case 1:
473 if (dE == 2) {
474 PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]);
475
476 J[1] = -J[2] / norm;
477 J[3] = J[0] / norm;
478 } else {
479 PetscReal inorm = 1. / PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]);
480 PetscReal x = J[0] * inorm;
481 PetscReal y = J[3] * inorm;
482 PetscReal z = J[6] * inorm;
483
484 if (x > 0.) {
485 PetscReal inv1pX = 1. / (1. + x);
486
487 J[1] = -y;
488 J[2] = -z;
489 J[4] = 1. - y * y * inv1pX;
490 J[5] = -y * z * inv1pX;
491 J[7] = -y * z * inv1pX;
492 J[8] = 1. - z * z * inv1pX;
493 } else {
494 PetscReal inv1mX = 1. / (1. - x);
495
496 J[1] = z;
497 J[2] = y;
498 J[4] = -y * z * inv1mX;
499 J[5] = 1. - y * y * inv1mX;
500 J[7] = 1. - z * z * inv1mX;
501 J[8] = -y * z * inv1mX;
502 }
503 }
504 break;
505 case 2: {
506 PetscReal inorm;
507
508 J[2] = J[3] * J[7] - J[6] * J[4];
509 J[5] = J[6] * J[1] - J[0] * J[7];
510 J[8] = J[0] * J[4] - J[3] * J[1];
511
512 inorm = 1. / PetscSqrtReal(J[2] * J[2] + J[5] * J[5] + J[8] * J[8]);
513
514 J[2] *= inorm;
515 J[5] *= inorm;
516 J[8] *= inorm;
517 } break;
518 }
519 for (j = 0; j < dE * dE; j++) g->J[i * dE * dE + j] = J[j];
520 }
521 }
522 PetscCall(PetscFEGeomComplete(g));
523 PetscCall(DMFieldGetDegree(field, pointIS, NULL, &maxDegree));
524 g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE;
525 if (faceData) PetscUseTypeMethod(field, computeFaceData, pointIS, quad, g);
526 *geom = g;
527 PetscFunctionReturn(PETSC_SUCCESS);
528 }
529