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