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