xref: /petsc/src/dm/field/interface/dmfield.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
13da551e6SToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
20145028aSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/
30145028aSToby Isaac #include <petscdmplex.h>
43da551e6SToby Isaac 
59371c9d4SSatish Balay const char *const DMFieldContinuities[] = {"VERTEX", "EDGE", "FACET", "CELL", NULL};
63da551e6SToby Isaac 
79371c9d4SSatish Balay PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm, PetscInt numComponents, DMFieldContinuity continuity, DMField *field) {
83da551e6SToby Isaac   DMField b;
93da551e6SToby Isaac 
103da551e6SToby Isaac   PetscFunctionBegin;
113da551e6SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
123da551e6SToby Isaac   PetscValidPointer(field, 2);
139566063dSJacob Faibussowitsch   PetscCall(DMFieldInitializePackage());
143da551e6SToby Isaac 
159566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(b, DMFIELD_CLASSID, "DMField", "Field over DM", "DM", PetscObjectComm((PetscObject)dm), DMFieldDestroy, DMFieldView));
169566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dm));
173da551e6SToby Isaac   b->dm            = dm;
183da551e6SToby Isaac   b->continuity    = continuity;
193da551e6SToby Isaac   b->numComponents = numComponents;
203da551e6SToby Isaac   *field           = b;
213da551e6SToby Isaac   PetscFunctionReturn(0);
223da551e6SToby Isaac }
233da551e6SToby Isaac 
243da551e6SToby Isaac /*@
253da551e6SToby Isaac    DMFieldDestroy - destroy a DMField
263da551e6SToby Isaac 
273da551e6SToby Isaac    Collective
283da551e6SToby Isaac 
294165533cSJose E. Roman    Input Parameter:
303da551e6SToby Isaac .  field - address of DMField
313da551e6SToby Isaac 
323da551e6SToby Isaac    Level: advanced
333da551e6SToby Isaac 
34db781477SPatrick Sanan .seealso: `DMFieldCreate()`
353da551e6SToby Isaac @*/
369371c9d4SSatish Balay PetscErrorCode DMFieldDestroy(DMField *field) {
373da551e6SToby Isaac   PetscFunctionBegin;
383da551e6SToby Isaac   if (!*field) PetscFunctionReturn(0);
393da551e6SToby Isaac   PetscValidHeaderSpecific((*field), DMFIELD_CLASSID, 1);
409371c9d4SSatish Balay   if (--((PetscObject)(*field))->refct > 0) {
419371c9d4SSatish Balay     *field = NULL;
429371c9d4SSatish Balay     PetscFunctionReturn(0);
439371c9d4SSatish Balay   }
44dbbe0bcdSBarry Smith   PetscTryTypeMethod((*field), destroy);
459566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&((*field)->dm)));
469566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(field));
473da551e6SToby Isaac   PetscFunctionReturn(0);
483da551e6SToby Isaac }
493da551e6SToby Isaac 
503da551e6SToby Isaac /*@C
513da551e6SToby Isaac    DMFieldView - view a DMField
523da551e6SToby Isaac 
533da551e6SToby Isaac    Collective
543da551e6SToby Isaac 
554165533cSJose E. Roman    Input Parameters:
563da551e6SToby Isaac +  field - DMField
573da551e6SToby Isaac -  viewer - viewer to display field, for example PETSC_VIEWER_STDOUT_WORLD
583da551e6SToby Isaac 
593da551e6SToby Isaac    Level: advanced
603da551e6SToby Isaac 
61db781477SPatrick Sanan .seealso: `DMFieldCreate()`
623da551e6SToby Isaac @*/
639371c9d4SSatish Balay PetscErrorCode DMFieldView(DMField field, PetscViewer viewer) {
643da551e6SToby Isaac   PetscBool iascii;
653da551e6SToby Isaac 
663da551e6SToby Isaac   PetscFunctionBegin;
673da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
689566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field), &viewer));
693da551e6SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
703da551e6SToby Isaac   PetscCheckSameComm(field, 1, viewer, 2);
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
723da551e6SToby Isaac   if (iascii) {
739566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)field, viewer));
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
7563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " components\n", field->numComponents));
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%s continuity\n", DMFieldContinuities[field->continuity]));
779566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_DEFAULT));
789566063dSJacob Faibussowitsch     PetscCall(DMView(field->dm, viewer));
799566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
803da551e6SToby Isaac   }
81dbbe0bcdSBarry Smith   PetscTryTypeMethod(field, view, viewer);
821baa6e33SBarry Smith   if (iascii) PetscCall(PetscViewerASCIIPopTab(viewer));
833da551e6SToby Isaac   PetscFunctionReturn(0);
843da551e6SToby Isaac }
853da551e6SToby Isaac 
863da551e6SToby Isaac /*@C
873da551e6SToby Isaac    DMFieldSetType - set the DMField implementation
883da551e6SToby Isaac 
89d083f849SBarry Smith    Collective on field
903da551e6SToby Isaac 
913da551e6SToby Isaac    Input Parameters:
923da551e6SToby Isaac +  field - the DMField context
933da551e6SToby Isaac -  type - a known method
943da551e6SToby Isaac 
953da551e6SToby Isaac    Notes:
963da551e6SToby Isaac    See "include/petscvec.h" for available methods (for instance)
973da551e6SToby Isaac +    DMFIELDDA    - a field defined only by its values at the corners of a DMDA
983da551e6SToby Isaac .    DMFIELDDS    - a field defined by a discretization over a mesh set with DMSetField()
993da551e6SToby Isaac -    DMFIELDSHELL - a field defined by arbitrary callbacks
1003da551e6SToby Isaac 
1013da551e6SToby Isaac   Level: advanced
1023da551e6SToby Isaac 
103db781477SPatrick Sanan .seealso: `DMFieldType`,
1043da551e6SToby Isaac @*/
1059371c9d4SSatish Balay PetscErrorCode DMFieldSetType(DMField field, DMFieldType type) {
1063da551e6SToby Isaac   PetscBool match;
1075f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(DMField);
1083da551e6SToby Isaac 
1093da551e6SToby Isaac   PetscFunctionBegin;
1103da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
1113da551e6SToby Isaac   PetscValidCharPointer(type, 2);
1123da551e6SToby Isaac 
1139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)field, type, &match));
1143da551e6SToby Isaac   if (match) PetscFunctionReturn(0);
1153da551e6SToby Isaac 
1169566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(DMFieldList, type, &r));
1175f80ce2aSJacob Faibussowitsch   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested DMField type %s", type);
1183da551e6SToby Isaac   /* Destroy the previous private DMField context */
119dbbe0bcdSBarry Smith   PetscTryTypeMethod(field, destroy);
1205f80ce2aSJacob Faibussowitsch 
1219566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(field->ops, sizeof(*field->ops)));
1229566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)field, type));
1233da551e6SToby Isaac   field->ops->create = r;
1249566063dSJacob Faibussowitsch   PetscCall((*r)(field));
1253da551e6SToby Isaac   PetscFunctionReturn(0);
1263da551e6SToby Isaac }
1273da551e6SToby Isaac 
1283da551e6SToby Isaac /*@C
1293da551e6SToby Isaac   DMFieldGetType - Gets the DMField type name (as a string) from the DMField.
1303da551e6SToby Isaac 
1313da551e6SToby Isaac   Not Collective
1323da551e6SToby Isaac 
1333da551e6SToby Isaac   Input Parameter:
1343da551e6SToby Isaac . field  - The DMField context
1353da551e6SToby Isaac 
1363da551e6SToby Isaac   Output Parameter:
1373da551e6SToby Isaac . type - The DMField type name
1383da551e6SToby Isaac 
1393da551e6SToby Isaac   Level: advanced
1403da551e6SToby Isaac 
141db781477SPatrick Sanan .seealso: `DMFieldSetType()`
1423da551e6SToby Isaac @*/
1439371c9d4SSatish Balay PetscErrorCode DMFieldGetType(DMField field, DMFieldType *type) {
1443da551e6SToby Isaac   PetscFunctionBegin;
145064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
1463da551e6SToby Isaac   PetscValidPointer(type, 2);
1479566063dSJacob Faibussowitsch   PetscCall(DMFieldRegisterAll());
1483da551e6SToby Isaac   *type = ((PetscObject)field)->type_name;
1493da551e6SToby Isaac   PetscFunctionReturn(0);
1503da551e6SToby Isaac }
1513da551e6SToby Isaac 
152da311338SToby Isaac /*@
153da311338SToby Isaac   DMFieldGetNumComponents - Returns the number of components in the field
154da311338SToby Isaac 
155da311338SToby Isaac   Not collective
156da311338SToby Isaac 
157da311338SToby Isaac   Input Parameter:
158da311338SToby Isaac . field - The DMField object
159da311338SToby Isaac 
160da311338SToby Isaac   Output Parameter:
161da311338SToby Isaac . nc - The number of field components
162da311338SToby Isaac 
163da311338SToby Isaac   Level: intermediate
164da311338SToby Isaac 
165db781477SPatrick Sanan .seealso: `DMFieldEvaluate()`
166da311338SToby Isaac @*/
1679371c9d4SSatish Balay PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc) {
1683da551e6SToby Isaac   PetscFunctionBegin;
1693da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
1703da551e6SToby Isaac   PetscValidIntPointer(nc, 2);
1713da551e6SToby Isaac   *nc = field->numComponents;
1723da551e6SToby Isaac   PetscFunctionReturn(0);
1733da551e6SToby Isaac }
1743da551e6SToby Isaac 
175da311338SToby Isaac /*@
176da311338SToby Isaac   DMFieldGetDM - Returns the DM for the manifold over which the field is defined.
177da311338SToby Isaac 
178da311338SToby Isaac   Not collective
179da311338SToby Isaac 
180da311338SToby Isaac   Input Parameter:
181da311338SToby Isaac . field - The DMField object
182da311338SToby Isaac 
183da311338SToby Isaac   Output Parameter:
184da311338SToby Isaac . dm - The DM object
185da311338SToby Isaac 
186da311338SToby Isaac   Level: intermediate
187da311338SToby Isaac 
188db781477SPatrick Sanan .seealso: `DMFieldEvaluate()`
189da311338SToby Isaac @*/
1909371c9d4SSatish Balay PetscErrorCode DMFieldGetDM(DMField field, DM *dm) {
1913da551e6SToby Isaac   PetscFunctionBegin;
1923da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
1933da551e6SToby Isaac   PetscValidPointer(dm, 2);
1943da551e6SToby Isaac   *dm = field->dm;
1953da551e6SToby Isaac   PetscFunctionReturn(0);
1963da551e6SToby Isaac }
1973da551e6SToby Isaac 
198da311338SToby Isaac /*@
199da311338SToby Isaac   DMFieldEvaluate - Evaluate the field and its derivatives on a set of points
200da311338SToby Isaac 
201d083f849SBarry Smith   Collective on points
202da311338SToby Isaac 
203d8d19677SJose E. Roman   Input Parameters:
204da311338SToby Isaac + field - The DMField object
205da311338SToby Isaac . points - The points at which to evaluate the field.  Should have size d x n,
206da311338SToby Isaac            where d is the coordinate dimension of the manifold and n is the number
207da311338SToby Isaac            of points
208da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR.
209da311338SToby Isaac              If the field is complex and datatype is PETSC_REAL, the real part of the
210da311338SToby Isaac              field is returned.
211da311338SToby Isaac 
212d8d19677SJose E. Roman   Output Parameters:
213da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
214da311338SToby Isaac       If B is not NULL, the values of the field are written in this array, varying first by component,
215da311338SToby Isaac       then by point.
216da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype).
217da311338SToby Isaac       If D is not NULL, the values of the field's spatial derivatives are written in this array,
218da311338SToby Isaac       varying first by the partial derivative component, then by field component, then by point.
219da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype).
220da311338SToby Isaac       If H is not NULL, the values of the field's second spatial derivatives are written in this array,
221da311338SToby Isaac       varying first by the second partial derivative component, then by field component, then by point.
222da311338SToby Isaac 
223da311338SToby Isaac   Level: intermediate
224da311338SToby Isaac 
225db781477SPatrick Sanan .seealso: `DMFieldGetDM()`, `DMFieldGetNumComponents()`, `DMFieldEvaluateFE()`, `DMFieldEvaluateFV()`
226da311338SToby Isaac @*/
2279371c9d4SSatish Balay PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H) {
2283da551e6SToby Isaac   PetscFunctionBegin;
2293da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
2303da551e6SToby Isaac   PetscValidHeaderSpecific(points, VEC_CLASSID, 2);
231064a246eSJacob Faibussowitsch   if (B) PetscValidPointer(B, 4);
232064a246eSJacob Faibussowitsch   if (D) PetscValidPointer(D, 5);
233064a246eSJacob Faibussowitsch   if (H) PetscValidPointer(H, 6);
2343da551e6SToby Isaac   if (field->ops->evaluate) {
2359566063dSJacob Faibussowitsch     PetscCall((*field->ops->evaluate)(field, points, datatype, B, D, H));
2363da551e6SToby Isaac   } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented for this type");
2373da551e6SToby Isaac   PetscFunctionReturn(0);
2383da551e6SToby Isaac }
2393da551e6SToby Isaac 
240da311338SToby Isaac /*@
241da311338SToby Isaac   DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from
242da311338SToby Isaac   quadrature points on a reference point.  The derivatives are taken with respect to the
243da311338SToby Isaac   reference coordinates.
244da311338SToby Isaac 
245da311338SToby Isaac   Not collective
246da311338SToby Isaac 
247d8d19677SJose E. Roman   Input Parameters:
248da311338SToby Isaac + field - The DMField object
249da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field
250da311338SToby Isaac . points - The quadature containing the points in the reference cell at which to evaluate the field.
251da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR.
252da311338SToby Isaac              If the field is complex and datatype is PETSC_REAL, the real part of the
253da311338SToby Isaac              field is returned.
254da311338SToby Isaac 
255d8d19677SJose E. Roman   Output Parameters:
256da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
257da311338SToby Isaac       If B is not NULL, the values of the field are written in this array, varying first by component,
258da311338SToby Isaac       then by point.
259da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype).
260da311338SToby Isaac       If D is not NULL, the values of the field's spatial derivatives are written in this array,
261da311338SToby Isaac       varying first by the partial derivative component, then by field component, then by point.
262da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype).
263da311338SToby Isaac       If H is not NULL, the values of the field's second spatial derivatives are written in this array,
264da311338SToby Isaac       varying first by the second partial derivative component, then by field component, then by point.
265da311338SToby Isaac 
266da311338SToby Isaac   Level: intermediate
267da311338SToby Isaac 
268db781477SPatrick Sanan .seealso: `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFV()`
269da311338SToby Isaac @*/
2709371c9d4SSatish Balay PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H) {
2713da551e6SToby Isaac   PetscFunctionBegin;
2723da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
2733da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2);
2743da551e6SToby Isaac   PetscValidHeader(points, 3);
275064a246eSJacob Faibussowitsch   if (B) PetscValidPointer(B, 5);
276064a246eSJacob Faibussowitsch   if (D) PetscValidPointer(D, 6);
277064a246eSJacob Faibussowitsch   if (H) PetscValidPointer(H, 7);
2783da551e6SToby Isaac   if (field->ops->evaluateFE) {
2799566063dSJacob Faibussowitsch     PetscCall((*field->ops->evaluateFE)(field, cellIS, points, datatype, B, D, H));
2803da551e6SToby Isaac   } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented for this type");
2813da551e6SToby Isaac   PetscFunctionReturn(0);
2823da551e6SToby Isaac }
2833da551e6SToby Isaac 
284da311338SToby Isaac /*@
285da311338SToby Isaac   DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points.
286da311338SToby Isaac 
287da311338SToby Isaac   Not collective
288da311338SToby Isaac 
289d8d19677SJose E. Roman   Input Parameters:
290da311338SToby Isaac + field - The DMField object
291da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field
292da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR.
293da311338SToby Isaac              If the field is complex and datatype is PETSC_REAL, the real part of the
294da311338SToby Isaac              field is returned.
295da311338SToby Isaac 
296d8d19677SJose E. Roman   Output Parameters:
297da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
298da311338SToby Isaac       If B is not NULL, the values of the field are written in this array, varying first by component,
299da311338SToby Isaac       then by point.
300da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype).
301da311338SToby Isaac       If D is not NULL, the values of the field's spatial derivatives are written in this array,
302da311338SToby Isaac       varying first by the partial derivative component, then by field component, then by point.
303da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype).
304da311338SToby Isaac       If H is not NULL, the values of the field's second spatial derivatives are written in this array,
305da311338SToby Isaac       varying first by the second partial derivative component, then by field component, then by point.
306da311338SToby Isaac 
307da311338SToby Isaac   Level: intermediate
308da311338SToby Isaac 
309db781477SPatrick Sanan .seealso: `DMFieldGetNumComponents()`, `DMFieldEvaluate()`, `DMFieldEvaluateFE()`
310da311338SToby Isaac @*/
3119371c9d4SSatish Balay PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H) {
3123da551e6SToby Isaac   PetscFunctionBegin;
3133da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
3143da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2);
315064a246eSJacob Faibussowitsch   if (B) PetscValidPointer(B, 4);
316064a246eSJacob Faibussowitsch   if (D) PetscValidPointer(D, 5);
317064a246eSJacob Faibussowitsch   if (H) PetscValidPointer(H, 6);
3183da551e6SToby Isaac   if (field->ops->evaluateFV) {
3199566063dSJacob Faibussowitsch     PetscCall((*field->ops->evaluateFV)(field, cellIS, datatype, B, D, H));
3203da551e6SToby Isaac   } else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Not implemented for this type");
3213da551e6SToby Isaac   PetscFunctionReturn(0);
3223da551e6SToby Isaac }
3233da551e6SToby Isaac 
324da311338SToby Isaac /*@
325b7260050SToby Isaac   DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the
326da311338SToby Isaac   reference element
327da311338SToby Isaac 
328da311338SToby Isaac   Not collective
329da311338SToby Isaac 
3304165533cSJose E. Roman   Input Parameters:
331da311338SToby Isaac + field - the DMField object
332da311338SToby Isaac - cellIS - the index set of points over which we want know the invariance
333da311338SToby Isaac 
3344165533cSJose E. Roman   Output Parameters:
335b7260050SToby Isaac + minDegree - the degree of the largest polynomial space contained in the field on each element
336b7260050SToby Isaac - maxDegree - the largest degree of the smallest polynomial space containing the field on any element
337da311338SToby Isaac 
338da311338SToby Isaac   Level: intermediate
339da311338SToby Isaac 
340db781477SPatrick Sanan .seealso: `DMFieldEvaluateFE()`
341da311338SToby Isaac @*/
3429371c9d4SSatish Balay PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PetscInt *minDegree, PetscInt *maxDegree) {
3433da551e6SToby Isaac   PetscFunctionBegin;
3443da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
3453da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS, IS_CLASSID, 2);
346dadcf809SJacob Faibussowitsch   if (minDegree) PetscValidIntPointer(minDegree, 3);
347dadcf809SJacob Faibussowitsch   if (maxDegree) PetscValidIntPointer(maxDegree, 4);
3483da551e6SToby Isaac 
349b7260050SToby Isaac   if (minDegree) *minDegree = -1;
350b7260050SToby Isaac   if (maxDegree) *maxDegree = PETSC_MAX_INT;
3513da551e6SToby Isaac 
3521baa6e33SBarry Smith   if (field->ops->getDegree) PetscCall((*field->ops->getDegree)(field, cellIS, minDegree, maxDegree));
3533da551e6SToby Isaac   PetscFunctionReturn(0);
3543da551e6SToby Isaac }
3553da551e6SToby Isaac 
356da311338SToby Isaac /*@
357da311338SToby Isaac   DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected
358da311338SToby Isaac   points via pullback onto the reference element
359da311338SToby Isaac 
360da311338SToby Isaac   Not collective
361da311338SToby Isaac 
3624165533cSJose E. Roman   Input Parameters:
363da311338SToby Isaac + field - the DMField object
364da311338SToby Isaac - pointIS - the index set of points over which we wish to integrate the field
365da311338SToby Isaac 
3664165533cSJose E. Roman   Output Parameter:
367da311338SToby Isaac . quad - a PetscQuadrature object
368da311338SToby Isaac 
369da311338SToby Isaac   Level: developer
370da311338SToby Isaac 
371db781477SPatrick Sanan .seealso: `DMFieldEvaluteFE()`, `DMFieldGetDegree()`
372da311338SToby Isaac @*/
3739371c9d4SSatish Balay PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad) {
3743da551e6SToby Isaac   PetscFunctionBegin;
3753da551e6SToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
3763da551e6SToby Isaac   PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2);
3773da551e6SToby Isaac   PetscValidPointer(quad, 3);
3783da551e6SToby Isaac 
3793da551e6SToby Isaac   *quad = NULL;
380dbbe0bcdSBarry Smith   PetscTryTypeMethod(field, createDefaultQuadrature, pointIS, quad);
3813da551e6SToby Isaac   PetscFunctionReturn(0);
3823da551e6SToby Isaac }
3833da551e6SToby Isaac 
384da311338SToby Isaac /*@C
385da311338SToby Isaac   DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field
386da311338SToby Isaac 
387da311338SToby Isaac   Not collective
388da311338SToby Isaac 
3894165533cSJose E. Roman   Input Parameters:
390da311338SToby Isaac + field - the DMField object
391da311338SToby Isaac . pointIS - the index set of points over which we wish to integrate the field
392da311338SToby Isaac . quad - the quadrature points at which to evaluate the geometric factors
393da311338SToby Isaac - faceData - whether additional data for facets (the normal vectors and adjacent cells) should
394da311338SToby Isaac   be calculated
395da311338SToby Isaac 
3964165533cSJose E. Roman   Output Parameter:
397da311338SToby Isaac . geom - the geometric factors
398da311338SToby Isaac 
399da311338SToby Isaac   Level: developer
400da311338SToby Isaac 
401db781477SPatrick Sanan .seealso: `DMFieldEvaluateFE()`, `DMFieldCreateDefaulteQuadrature()`, `DMFieldGetDegree()`
40206dd6b0eSSatish Balay @*/
4039371c9d4SSatish Balay PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) {
4042bb1671aSToby Isaac   PetscInt     dim, dE;
4052bb1671aSToby Isaac   PetscInt     nPoints;
406b7260050SToby Isaac   PetscInt     maxDegree;
4072bb1671aSToby Isaac   PetscFEGeom *g;
4082bb1671aSToby Isaac 
4092bb1671aSToby Isaac   PetscFunctionBegin;
4102bb1671aSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
4112bb1671aSToby Isaac   PetscValidHeaderSpecific(pointIS, IS_CLASSID, 2);
4122bb1671aSToby Isaac   PetscValidHeader(quad, 3);
4139566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &nPoints));
4142bb1671aSToby Isaac   dE = field->numComponents;
4159566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomCreate(quad, nPoints, dE, faceData, &g));
4169566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluateFE(field, pointIS, quad, PETSC_REAL, g->v, g->J, NULL));
4172bb1671aSToby Isaac   dim = g->dim;
4182bb1671aSToby Isaac   if (dE > dim) {
4192bb1671aSToby Isaac     /* space out J and make square Jacobians */
4202bb1671aSToby Isaac     PetscInt i, j, k, N = g->numPoints * g->numCells;
4212bb1671aSToby Isaac 
4222bb1671aSToby Isaac     for (i = N - 1; i >= 0; i--) {
4232bb1671aSToby Isaac       PetscReal J[9];
4242bb1671aSToby Isaac 
4252bb1671aSToby Isaac       for (j = 0; j < dE; j++) {
4269371c9d4SSatish Balay         for (k = 0; k < dim; k++) { J[j * dE + k] = g->J[i * dE * dim + j * dim + k]; }
4272bb1671aSToby Isaac       }
4282bb1671aSToby Isaac       switch (dim) {
4292bb1671aSToby Isaac       case 0:
4302bb1671aSToby Isaac         for (j = 0; j < dE; j++) {
4319371c9d4SSatish Balay           for (k = 0; k < dE; k++) { J[j * dE + k] = (j == k) ? 1. : 0.; }
4322bb1671aSToby Isaac         }
4332bb1671aSToby Isaac         break;
4342bb1671aSToby Isaac       case 1:
4352bb1671aSToby Isaac         if (dE == 2) {
4362bb1671aSToby Isaac           PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]);
4372bb1671aSToby Isaac 
4382bb1671aSToby Isaac           J[1] = -J[2] / norm;
4392bb1671aSToby Isaac           J[3] = J[0] / norm;
4402bb1671aSToby Isaac         } else {
4412bb1671aSToby Isaac           PetscReal inorm = 1. / PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]);
4422bb1671aSToby Isaac           PetscReal x     = J[0] * inorm;
4432bb1671aSToby Isaac           PetscReal y     = J[3] * inorm;
4442bb1671aSToby Isaac           PetscReal z     = J[6] * inorm;
4452bb1671aSToby Isaac 
4462bb1671aSToby Isaac           if (x > 0.) {
4472bb1671aSToby Isaac             PetscReal inv1pX = 1. / (1. + x);
4482bb1671aSToby Isaac 
4499371c9d4SSatish Balay             J[1] = -y;
4509371c9d4SSatish Balay             J[2] = -z;
4519371c9d4SSatish Balay             J[4] = 1. - y * y * inv1pX;
4529371c9d4SSatish Balay             J[5] = -y * z * inv1pX;
4539371c9d4SSatish Balay             J[7] = -y * z * inv1pX;
4549371c9d4SSatish Balay             J[8] = 1. - z * z * inv1pX;
4552bb1671aSToby Isaac           } else {
4562bb1671aSToby Isaac             PetscReal inv1mX = 1. / (1. - x);
4572bb1671aSToby Isaac 
4589371c9d4SSatish Balay             J[1] = z;
4599371c9d4SSatish Balay             J[2] = y;
4609371c9d4SSatish Balay             J[4] = -y * z * inv1mX;
4619371c9d4SSatish Balay             J[5] = 1. - y * y * inv1mX;
4629371c9d4SSatish Balay             J[7] = 1. - z * z * inv1mX;
4639371c9d4SSatish Balay             J[8] = -y * z * inv1mX;
4642bb1671aSToby Isaac           }
4652bb1671aSToby Isaac         }
4662bb1671aSToby Isaac         break;
4679371c9d4SSatish Balay       case 2: {
4682bb1671aSToby Isaac         PetscReal inorm;
4692bb1671aSToby Isaac 
4702bb1671aSToby Isaac         J[2] = J[3] * J[7] - J[6] * J[4];
4712bb1671aSToby Isaac         J[5] = J[6] * J[1] - J[0] * J[7];
4722bb1671aSToby Isaac         J[8] = J[0] * J[4] - J[3] * J[1];
4732bb1671aSToby Isaac 
4742bb1671aSToby Isaac         inorm = 1. / PetscSqrtReal(J[2] * J[2] + J[5] * J[5] + J[8] * J[8]);
4752bb1671aSToby Isaac 
4762bb1671aSToby Isaac         J[2] *= inorm;
4772bb1671aSToby Isaac         J[5] *= inorm;
4782bb1671aSToby Isaac         J[8] *= inorm;
4799371c9d4SSatish Balay       } break;
4802bb1671aSToby Isaac       }
4819371c9d4SSatish Balay       for (j = 0; j < dE * dE; j++) { g->J[i * dE * dE + j] = J[j]; }
4822bb1671aSToby Isaac     }
4832bb1671aSToby Isaac   }
4849566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomComplete(g));
4859566063dSJacob Faibussowitsch   PetscCall(DMFieldGetDegree(field, pointIS, NULL, &maxDegree));
486b7260050SToby Isaac   g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE;
487*48a46eb9SPierre Jolivet   if (faceData) PetscCall((*field->ops->computeFaceData)(field, pointIS, quad, g));
4882bb1671aSToby Isaac   *geom = g;
4892bb1671aSToby Isaac   PetscFunctionReturn(0);
4902bb1671aSToby Isaac }
491