xref: /petsc/src/dm/field/interface/dmfield.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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 
53da551e6SToby Isaac const char *const DMFieldContinuities[] = {
63da551e6SToby Isaac   "VERTEX",
73da551e6SToby Isaac   "EDGE",
83da551e6SToby Isaac   "FACET",
93da551e6SToby Isaac   "CELL",
10ea78f98cSLisandro Dalcin   NULL
113da551e6SToby Isaac };
123da551e6SToby Isaac 
133da551e6SToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm,PetscInt numComponents,DMFieldContinuity continuity,DMField *field)
143da551e6SToby Isaac {
153da551e6SToby Isaac   PetscErrorCode ierr;
163da551e6SToby Isaac   DMField        b;
173da551e6SToby Isaac 
183da551e6SToby Isaac   PetscFunctionBegin;
193da551e6SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
203da551e6SToby Isaac   PetscValidPointer(field,2);
213da551e6SToby Isaac   ierr = DMFieldInitializePackage();CHKERRQ(ierr);
223da551e6SToby Isaac 
233da551e6SToby Isaac   ierr = PetscHeaderCreate(b,DMFIELD_CLASSID,"DMField","Field over DM","DM",PetscObjectComm((PetscObject)dm),DMFieldDestroy,DMFieldView);CHKERRQ(ierr);
243da551e6SToby Isaac   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
253da551e6SToby Isaac   b->dm = dm;
263da551e6SToby Isaac   b->continuity = continuity;
273da551e6SToby Isaac   b->numComponents = numComponents;
283da551e6SToby Isaac   *field = b;
293da551e6SToby Isaac   PetscFunctionReturn(0);
303da551e6SToby Isaac }
313da551e6SToby Isaac 
323da551e6SToby Isaac /*@
333da551e6SToby Isaac    DMFieldDestroy - destroy a DMField
343da551e6SToby Isaac 
353da551e6SToby Isaac    Collective
363da551e6SToby Isaac 
373da551e6SToby Isaac    Input Arguments:
383da551e6SToby Isaac .  field - address of DMField
393da551e6SToby Isaac 
403da551e6SToby Isaac    Level: advanced
413da551e6SToby Isaac 
423da551e6SToby Isaac .seealso: DMFieldCreate()
433da551e6SToby Isaac @*/
443da551e6SToby Isaac PetscErrorCode DMFieldDestroy(DMField *field)
453da551e6SToby Isaac {
463da551e6SToby Isaac   PetscErrorCode ierr;
473da551e6SToby Isaac 
483da551e6SToby Isaac   PetscFunctionBegin;
493da551e6SToby Isaac   if (!*field) PetscFunctionReturn(0);
503da551e6SToby Isaac   PetscValidHeaderSpecific((*field),DMFIELD_CLASSID,1);
51ea78f98cSLisandro Dalcin   if (--((PetscObject)(*field))->refct > 0) {*field = NULL; PetscFunctionReturn(0);}
523da551e6SToby Isaac   if ((*field)->ops->destroy) {ierr = (*(*field)->ops->destroy)(*field);CHKERRQ(ierr);}
533da551e6SToby Isaac   ierr = DMDestroy(&((*field)->dm));CHKERRQ(ierr);
543da551e6SToby Isaac   ierr = PetscHeaderDestroy(field);CHKERRQ(ierr);
553da551e6SToby Isaac   PetscFunctionReturn(0);
563da551e6SToby Isaac }
573da551e6SToby Isaac 
583da551e6SToby Isaac /*@C
593da551e6SToby Isaac    DMFieldView - view a DMField
603da551e6SToby Isaac 
613da551e6SToby Isaac    Collective
623da551e6SToby Isaac 
633da551e6SToby Isaac    Input Arguments:
643da551e6SToby Isaac +  field - DMField
653da551e6SToby Isaac -  viewer - viewer to display field, for example PETSC_VIEWER_STDOUT_WORLD
663da551e6SToby Isaac 
673da551e6SToby Isaac    Level: advanced
683da551e6SToby Isaac 
693da551e6SToby Isaac .seealso: DMFieldCreate()
703da551e6SToby Isaac @*/
713da551e6SToby Isaac PetscErrorCode DMFieldView(DMField field,PetscViewer viewer)
723da551e6SToby Isaac {
733da551e6SToby Isaac   PetscErrorCode    ierr;
743da551e6SToby Isaac   PetscBool         iascii;
753da551e6SToby Isaac 
763da551e6SToby Isaac   PetscFunctionBegin;
773da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
783da551e6SToby Isaac   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field),&viewer);CHKERRQ(ierr);}
793da551e6SToby Isaac   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
803da551e6SToby Isaac   PetscCheckSameComm(field,1,viewer,2);
813da551e6SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
823da551e6SToby Isaac   if (iascii) {
833da551e6SToby Isaac     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)field,viewer);CHKERRQ(ierr);
843da551e6SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
853da551e6SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer,"%D components\n",field->numComponents);CHKERRQ(ierr);
863da551e6SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer,"%s continuity\n",DMFieldContinuities[field->continuity]);CHKERRQ(ierr);
873da551e6SToby Isaac     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
883da551e6SToby Isaac     ierr = DMView(field->dm,viewer);CHKERRQ(ierr);
893da551e6SToby Isaac     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
903da551e6SToby Isaac   }
913da551e6SToby Isaac   if (field->ops->view) {ierr = (*field->ops->view)(field,viewer);CHKERRQ(ierr);}
923da551e6SToby Isaac   if (iascii) {
933da551e6SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
943da551e6SToby Isaac   }
953da551e6SToby Isaac   PetscFunctionReturn(0);
963da551e6SToby Isaac }
973da551e6SToby Isaac 
983da551e6SToby Isaac /*@C
993da551e6SToby Isaac    DMFieldSetType - set the DMField implementation
1003da551e6SToby Isaac 
101d083f849SBarry Smith    Collective on field
1023da551e6SToby Isaac 
1033da551e6SToby Isaac    Input Parameters:
1043da551e6SToby Isaac +  field - the DMField context
1053da551e6SToby Isaac -  type - a known method
1063da551e6SToby Isaac 
1073da551e6SToby Isaac    Notes:
1083da551e6SToby Isaac    See "include/petscvec.h" for available methods (for instance)
1093da551e6SToby Isaac +    DMFIELDDA    - a field defined only by its values at the corners of a DMDA
1103da551e6SToby Isaac .    DMFIELDDS    - a field defined by a discretization over a mesh set with DMSetField()
1113da551e6SToby Isaac -    DMFIELDSHELL - a field defined by arbitrary callbacks
1123da551e6SToby Isaac 
1133da551e6SToby Isaac   Level: advanced
1143da551e6SToby Isaac 
1153da551e6SToby Isaac .seealso: DMFieldType,
1163da551e6SToby Isaac @*/
1173da551e6SToby Isaac PetscErrorCode DMFieldSetType(DMField field,DMFieldType type)
1183da551e6SToby Isaac {
1193da551e6SToby Isaac   PetscErrorCode ierr,(*r)(DMField);
1203da551e6SToby Isaac   PetscBool      match;
1213da551e6SToby Isaac 
1223da551e6SToby Isaac   PetscFunctionBegin;
1233da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1243da551e6SToby Isaac   PetscValidCharPointer(type,2);
1253da551e6SToby Isaac 
1263da551e6SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)field,type,&match);CHKERRQ(ierr);
1273da551e6SToby Isaac   if (match) PetscFunctionReturn(0);
1283da551e6SToby Isaac 
1293da551e6SToby Isaac   ierr = PetscFunctionListFind(DMFieldList,type,&r);CHKERRQ(ierr);
1303da551e6SToby Isaac   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DMField type %s",type);
1313da551e6SToby Isaac   /* Destroy the previous private DMField context */
1323da551e6SToby Isaac   if (field->ops->destroy) {
1333da551e6SToby Isaac     ierr = (*(field)->ops->destroy)(field);CHKERRQ(ierr);
1343da551e6SToby Isaac   }
1353da551e6SToby Isaac   ierr = PetscMemzero(field->ops,sizeof(*field->ops));CHKERRQ(ierr);
1363da551e6SToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject)field,type);CHKERRQ(ierr);
1373da551e6SToby Isaac   field->ops->create = r;
1383da551e6SToby Isaac   ierr = (*r)(field);CHKERRQ(ierr);
1393da551e6SToby Isaac   PetscFunctionReturn(0);
1403da551e6SToby Isaac }
1413da551e6SToby Isaac 
1423da551e6SToby Isaac /*@C
1433da551e6SToby Isaac   DMFieldGetType - Gets the DMField type name (as a string) from the DMField.
1443da551e6SToby Isaac 
1453da551e6SToby Isaac   Not Collective
1463da551e6SToby Isaac 
1473da551e6SToby Isaac   Input Parameter:
1483da551e6SToby Isaac . field  - The DMField context
1493da551e6SToby Isaac 
1503da551e6SToby Isaac   Output Parameter:
1513da551e6SToby Isaac . type - The DMField type name
1523da551e6SToby Isaac 
1533da551e6SToby Isaac   Level: advanced
1543da551e6SToby Isaac 
1553da551e6SToby Isaac .seealso: DMFieldSetType()
1563da551e6SToby Isaac @*/
1573da551e6SToby Isaac PetscErrorCode  DMFieldGetType(DMField field, DMFieldType *type)
1583da551e6SToby Isaac {
1593da551e6SToby Isaac   PetscErrorCode ierr;
1603da551e6SToby Isaac 
1613da551e6SToby Isaac   PetscFunctionBegin;
162064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(field, DMFIELD_CLASSID,1);
1633da551e6SToby Isaac   PetscValidPointer(type,2);
1643da551e6SToby Isaac   ierr = DMFieldRegisterAll();CHKERRQ(ierr);
1653da551e6SToby Isaac   *type = ((PetscObject)field)->type_name;
1663da551e6SToby Isaac   PetscFunctionReturn(0);
1673da551e6SToby Isaac }
1683da551e6SToby Isaac 
169da311338SToby Isaac /*@
170da311338SToby Isaac   DMFieldGetNumComponents - Returns the number of components in the field
171da311338SToby Isaac 
172da311338SToby Isaac   Not collective
173da311338SToby Isaac 
174da311338SToby Isaac   Input Parameter:
175da311338SToby Isaac . field - The DMField object
176da311338SToby Isaac 
177da311338SToby Isaac   Output Parameter:
178da311338SToby Isaac . nc - The number of field components
179da311338SToby Isaac 
180da311338SToby Isaac   Level: intermediate
181da311338SToby Isaac 
182da311338SToby Isaac .seealso: DMFieldEvaluate()
183da311338SToby Isaac @*/
1843da551e6SToby Isaac PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc)
1853da551e6SToby Isaac {
1863da551e6SToby Isaac   PetscFunctionBegin;
1873da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1883da551e6SToby Isaac   PetscValidIntPointer(nc,2);
1893da551e6SToby Isaac   *nc = field->numComponents;
1903da551e6SToby Isaac   PetscFunctionReturn(0);
1913da551e6SToby Isaac }
1923da551e6SToby Isaac 
193da311338SToby Isaac /*@
194da311338SToby Isaac   DMFieldGetDM - Returns the DM for the manifold over which the field is defined.
195da311338SToby Isaac 
196da311338SToby Isaac   Not collective
197da311338SToby Isaac 
198da311338SToby Isaac   Input Parameter:
199da311338SToby Isaac . field - The DMField object
200da311338SToby Isaac 
201da311338SToby Isaac   Output Parameter:
202da311338SToby Isaac . dm - The DM object
203da311338SToby Isaac 
204da311338SToby Isaac   Level: intermediate
205da311338SToby Isaac 
206da311338SToby Isaac .seealso: DMFieldEvaluate()
207da311338SToby Isaac @*/
2083da551e6SToby Isaac PetscErrorCode DMFieldGetDM(DMField field, DM *dm)
2093da551e6SToby Isaac {
2103da551e6SToby Isaac   PetscFunctionBegin;
2113da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2123da551e6SToby Isaac   PetscValidPointer(dm,2);
2133da551e6SToby Isaac   *dm = field->dm;
2143da551e6SToby Isaac   PetscFunctionReturn(0);
2153da551e6SToby Isaac }
2163da551e6SToby Isaac 
217da311338SToby Isaac /*@
218da311338SToby Isaac   DMFieldEvaluate - Evaluate the field and its derivatives on a set of points
219da311338SToby Isaac 
220d083f849SBarry Smith   Collective on points
221da311338SToby Isaac 
222*d8d19677SJose E. Roman   Input Parameters:
223da311338SToby Isaac + field - The DMField object
224da311338SToby Isaac . points - The points at which to evaluate the field.  Should have size d x n,
225da311338SToby Isaac            where d is the coordinate dimension of the manifold and n is the number
226da311338SToby Isaac            of points
227da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR.
228da311338SToby Isaac              If the field is complex and datatype is PETSC_REAL, the real part of the
229da311338SToby Isaac              field is returned.
230da311338SToby Isaac 
231*d8d19677SJose E. Roman   Output Parameters:
232da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
233da311338SToby Isaac       If B is not NULL, the values of the field are written in this array, varying first by component,
234da311338SToby Isaac       then by point.
235da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype).
236da311338SToby Isaac       If D is not NULL, the values of the field's spatial derivatives are written in this array,
237da311338SToby Isaac       varying first by the partial derivative component, then by field component, then by point.
238da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype).
239da311338SToby Isaac       If H is not NULL, the values of the field's second spatial derivatives are written in this array,
240da311338SToby Isaac       varying first by the second partial derivative component, then by field component, then by point.
241da311338SToby Isaac 
242da311338SToby Isaac   Level: intermediate
243da311338SToby Isaac 
244da311338SToby Isaac .seealso: DMFieldGetDM(), DMFieldGetNumComponents(), DMFieldEvaluateFE(), DMFieldEvaluateFV()
245da311338SToby Isaac @*/
2463da551e6SToby Isaac PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
2473da551e6SToby Isaac {
2483da551e6SToby Isaac   PetscErrorCode ierr;
2493da551e6SToby Isaac 
2503da551e6SToby Isaac   PetscFunctionBegin;
2513da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2523da551e6SToby Isaac   PetscValidHeaderSpecific(points,VEC_CLASSID,2);
253064a246eSJacob Faibussowitsch   if (B) PetscValidPointer(B,4);
254064a246eSJacob Faibussowitsch   if (D) PetscValidPointer(D,5);
255064a246eSJacob Faibussowitsch   if (H) PetscValidPointer(H,6);
2563da551e6SToby Isaac   if (field->ops->evaluate) {
2573da551e6SToby Isaac     ierr = (*field->ops->evaluate) (field, points, datatype, B, D, H);CHKERRQ(ierr);
2583da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
2593da551e6SToby Isaac   PetscFunctionReturn(0);
2603da551e6SToby Isaac }
2613da551e6SToby Isaac 
262da311338SToby Isaac /*@
263da311338SToby Isaac   DMFieldEvaluateFE - Evaluate the field and its derivatives on a set of points mapped from
264da311338SToby Isaac   quadrature points on a reference point.  The derivatives are taken with respect to the
265da311338SToby Isaac   reference coordinates.
266da311338SToby Isaac 
267da311338SToby Isaac   Not collective
268da311338SToby Isaac 
269*d8d19677SJose E. Roman   Input Parameters:
270da311338SToby Isaac + field - The DMField object
271da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field
272da311338SToby Isaac . points - The quadature containing the points in the reference cell at which to evaluate the field.
273da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR.
274da311338SToby Isaac              If the field is complex and datatype is PETSC_REAL, the real part of the
275da311338SToby Isaac              field is returned.
276da311338SToby Isaac 
277*d8d19677SJose E. Roman   Output Parameters:
278da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
279da311338SToby Isaac       If B is not NULL, the values of the field are written in this array, varying first by component,
280da311338SToby Isaac       then by point.
281da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype).
282da311338SToby Isaac       If D is not NULL, the values of the field's spatial derivatives are written in this array,
283da311338SToby Isaac       varying first by the partial derivative component, then by field component, then by point.
284da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype).
285da311338SToby Isaac       If H is not NULL, the values of the field's second spatial derivatives are written in this array,
286da311338SToby Isaac       varying first by the second partial derivative component, then by field component, then by point.
287da311338SToby Isaac 
288da311338SToby Isaac   Level: intermediate
289da311338SToby Isaac 
290da311338SToby Isaac .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFV()
291da311338SToby Isaac @*/
2923da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
2933da551e6SToby Isaac {
2943da551e6SToby Isaac   PetscErrorCode ierr;
2953da551e6SToby Isaac 
2963da551e6SToby Isaac   PetscFunctionBegin;
2973da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2983da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
2993da551e6SToby Isaac   PetscValidHeader(points,3);
300064a246eSJacob Faibussowitsch   if (B) PetscValidPointer(B,5);
301064a246eSJacob Faibussowitsch   if (D) PetscValidPointer(D,6);
302064a246eSJacob Faibussowitsch   if (H) PetscValidPointer(H,7);
3033da551e6SToby Isaac   if (field->ops->evaluateFE) {
3043da551e6SToby Isaac     ierr = (*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H);CHKERRQ(ierr);
3053da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
3063da551e6SToby Isaac   PetscFunctionReturn(0);
3073da551e6SToby Isaac }
3083da551e6SToby Isaac 
309da311338SToby Isaac /*@
310da311338SToby Isaac   DMFieldEvaluateFV - Evaluate the mean of a field and its finite volume derivatives on a set of points.
311da311338SToby Isaac 
312da311338SToby Isaac   Not collective
313da311338SToby Isaac 
314*d8d19677SJose E. Roman   Input Parameters:
315da311338SToby Isaac + field - The DMField object
316da311338SToby Isaac . cellIS - Index set for cells on which to evaluate the field
317da311338SToby Isaac - datatype - The PetscDataType of the output arrays: either PETSC_REAL or PETSC_SCALAR.
318da311338SToby Isaac              If the field is complex and datatype is PETSC_REAL, the real part of the
319da311338SToby Isaac              field is returned.
320da311338SToby Isaac 
321*d8d19677SJose E. Roman   Output Parameters:
322da311338SToby Isaac + B - pointer to data of size c * n * sizeof(datatype), where c is the number of components in the field.
323da311338SToby Isaac       If B is not NULL, the values of the field are written in this array, varying first by component,
324da311338SToby Isaac       then by point.
325da311338SToby Isaac . D - pointer to data of size d * c * n * sizeof(datatype).
326da311338SToby Isaac       If D is not NULL, the values of the field's spatial derivatives are written in this array,
327da311338SToby Isaac       varying first by the partial derivative component, then by field component, then by point.
328da311338SToby Isaac - H - pointer to data of size d * d * c * n * sizeof(datatype).
329da311338SToby Isaac       If H is not NULL, the values of the field's second spatial derivatives are written in this array,
330da311338SToby Isaac       varying first by the second partial derivative component, then by field component, then by point.
331da311338SToby Isaac 
332da311338SToby Isaac   Level: intermediate
333da311338SToby Isaac 
334da311338SToby Isaac .seealso: DMFieldGetNumComponents(), DMFieldEvaluate(), DMFieldEvaluateFE()
335da311338SToby Isaac @*/
3363da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
3373da551e6SToby Isaac {
3383da551e6SToby Isaac   PetscErrorCode ierr;
3393da551e6SToby Isaac 
3403da551e6SToby Isaac   PetscFunctionBegin;
3413da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
3423da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
343064a246eSJacob Faibussowitsch   if (B) PetscValidPointer(B,4);
344064a246eSJacob Faibussowitsch   if (D) PetscValidPointer(D,5);
345064a246eSJacob Faibussowitsch   if (H) PetscValidPointer(H,6);
3463da551e6SToby Isaac   if (field->ops->evaluateFV) {
3473da551e6SToby Isaac     ierr = (*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H);CHKERRQ(ierr);
3483da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
3493da551e6SToby Isaac   PetscFunctionReturn(0);
3503da551e6SToby Isaac }
3513da551e6SToby Isaac 
352da311338SToby Isaac /*@
353b7260050SToby Isaac   DMFieldGetDegree - Get the polynomial degree of a field when pulled back onto the
354da311338SToby Isaac   reference element
355da311338SToby Isaac 
356da311338SToby Isaac   Not collective
357da311338SToby Isaac 
358da311338SToby Isaac   Input Arguments:
359da311338SToby Isaac + field - the DMField object
360da311338SToby Isaac - cellIS - the index set of points over which we want know the invariance
361da311338SToby Isaac 
362da311338SToby Isaac   Output Arguments:
363b7260050SToby Isaac + minDegree - the degree of the largest polynomial space contained in the field on each element
364b7260050SToby Isaac - maxDegree - the largest degree of the smallest polynomial space containing the field on any element
365da311338SToby Isaac 
366da311338SToby Isaac   Level: intermediate
367da311338SToby Isaac 
368da311338SToby Isaac .seealso: DMFieldEvaluateFE()
369da311338SToby Isaac @*/
370b7260050SToby Isaac PetscErrorCode DMFieldGetDegree(DMField field, IS cellIS, PetscInt *minDegree, PetscInt *maxDegree)
3713da551e6SToby Isaac {
3723da551e6SToby Isaac   PetscErrorCode ierr;
3733da551e6SToby Isaac 
3743da551e6SToby Isaac   PetscFunctionBegin;
3753da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
3763da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
377b7260050SToby Isaac   if (minDegree) PetscValidPointer(minDegree,3);
378b7260050SToby Isaac   if (maxDegree) PetscValidPointer(maxDegree,4);
3793da551e6SToby Isaac 
380b7260050SToby Isaac   if (minDegree) *minDegree = -1;
381b7260050SToby Isaac   if (maxDegree) *maxDegree = PETSC_MAX_INT;
3823da551e6SToby Isaac 
383b7260050SToby Isaac   if (field->ops->getDegree) {
384b7260050SToby Isaac     ierr = (*field->ops->getDegree) (field,cellIS,minDegree,maxDegree);CHKERRQ(ierr);
3853da551e6SToby Isaac   }
3863da551e6SToby Isaac   PetscFunctionReturn(0);
3873da551e6SToby Isaac }
3883da551e6SToby Isaac 
389da311338SToby Isaac /*@
390da311338SToby Isaac   DMFieldCreateDefaultQuadrature - Creates a quadrature sufficient to integrate the field on the selected
391da311338SToby Isaac   points via pullback onto the reference element
392da311338SToby Isaac 
393da311338SToby Isaac   Not collective
394da311338SToby Isaac 
395da311338SToby Isaac   Input Arguments:
396da311338SToby Isaac + field - the DMField object
397da311338SToby Isaac - pointIS - the index set of points over which we wish to integrate the field
398da311338SToby Isaac 
399da311338SToby Isaac   Output Arguments:
400da311338SToby Isaac . quad - a PetscQuadrature object
401da311338SToby Isaac 
402da311338SToby Isaac   Level: developer
403da311338SToby Isaac 
404b7260050SToby Isaac .seealso: DMFieldEvaluteFE(), DMFieldGetDegree()
405da311338SToby Isaac @*/
4063da551e6SToby Isaac PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad)
4073da551e6SToby Isaac {
4083da551e6SToby Isaac   PetscErrorCode ierr;
4093da551e6SToby Isaac 
4103da551e6SToby Isaac   PetscFunctionBegin;
4113da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
4123da551e6SToby Isaac   PetscValidHeaderSpecific(pointIS,IS_CLASSID,2);
4133da551e6SToby Isaac   PetscValidPointer(quad,3);
4143da551e6SToby Isaac 
4153da551e6SToby Isaac   *quad = NULL;
4163da551e6SToby Isaac   if (field->ops->createDefaultQuadrature) {
4173da551e6SToby Isaac     ierr = (*field->ops->createDefaultQuadrature)(field, pointIS, quad);CHKERRQ(ierr);
4183da551e6SToby Isaac   }
4193da551e6SToby Isaac   PetscFunctionReturn(0);
4203da551e6SToby Isaac }
4213da551e6SToby Isaac 
422da311338SToby Isaac /*@C
423da311338SToby Isaac   DMFieldCreateFEGeom - Compute and create the geometric factors of a coordinate field
424da311338SToby Isaac 
425da311338SToby Isaac   Not collective
426da311338SToby Isaac 
427da311338SToby Isaac   Input Arguments:
428da311338SToby Isaac + field - the DMField object
429da311338SToby Isaac . pointIS - the index set of points over which we wish to integrate the field
430da311338SToby Isaac . quad - the quadrature points at which to evaluate the geometric factors
431da311338SToby Isaac - faceData - whether additional data for facets (the normal vectors and adjacent cells) should
432da311338SToby Isaac   be calculated
433da311338SToby Isaac 
434da311338SToby Isaac   Output Arguments:
435da311338SToby Isaac . geom - the geometric factors
436da311338SToby Isaac 
437da311338SToby Isaac   Level: developer
438da311338SToby Isaac 
439b7260050SToby Isaac .seealso: DMFieldEvaluateFE(), DMFieldCreateDefaulteQuadrature(), DMFieldGetDegree()
44006dd6b0eSSatish Balay @*/
4412bb1671aSToby Isaac PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
4422bb1671aSToby Isaac {
4432bb1671aSToby Isaac   PetscInt       dim, dE;
4442bb1671aSToby Isaac   PetscInt       nPoints;
445b7260050SToby Isaac   PetscInt       maxDegree;
4462bb1671aSToby Isaac   PetscFEGeom    *g;
4472bb1671aSToby Isaac   PetscErrorCode ierr;
4482bb1671aSToby Isaac 
4492bb1671aSToby Isaac   PetscFunctionBegin;
4502bb1671aSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
4512bb1671aSToby Isaac   PetscValidHeaderSpecific(pointIS,IS_CLASSID,2);
4522bb1671aSToby Isaac   PetscValidHeader(quad,3);
4532bb1671aSToby Isaac   ierr = ISGetLocalSize(pointIS,&nPoints);CHKERRQ(ierr);
4542bb1671aSToby Isaac   dE = field->numComponents;
4552bb1671aSToby Isaac   ierr = PetscFEGeomCreate(quad,nPoints,dE,faceData,&g);CHKERRQ(ierr);
4562bb1671aSToby Isaac   ierr = DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL);CHKERRQ(ierr);
4572bb1671aSToby Isaac   dim = g->dim;
4582bb1671aSToby Isaac   if (dE > dim) {
4592bb1671aSToby Isaac     /* space out J and make square Jacobians */
4602bb1671aSToby Isaac     PetscInt  i, j, k, N = g->numPoints * g->numCells;
4612bb1671aSToby Isaac 
4622bb1671aSToby Isaac     for (i = N-1; i >= 0; i--) {
4632bb1671aSToby Isaac       PetscReal   J[9];
4642bb1671aSToby Isaac 
4652bb1671aSToby Isaac       for (j = 0; j < dE; j++) {
4662bb1671aSToby Isaac         for (k = 0; k < dim; k++) {
4672bb1671aSToby Isaac           J[j*dE + k] = g->J[i*dE*dim + j*dim + k];
4682bb1671aSToby Isaac         }
4692bb1671aSToby Isaac       }
4702bb1671aSToby Isaac       switch (dim) {
4712bb1671aSToby Isaac       case 0:
4722bb1671aSToby Isaac         for (j = 0; j < dE; j++) {
4732bb1671aSToby Isaac           for (k = 0; k < dE; k++) {
4742bb1671aSToby Isaac             J[j * dE + k] = (j == k) ? 1. : 0.;
4752bb1671aSToby Isaac           }
4762bb1671aSToby Isaac         }
4772bb1671aSToby Isaac         break;
4782bb1671aSToby Isaac       case 1:
4792bb1671aSToby Isaac         if (dE == 2) {
4802bb1671aSToby Isaac           PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]);
4812bb1671aSToby Isaac 
4822bb1671aSToby Isaac           J[1] = -J[2] / norm;
4832bb1671aSToby Isaac           J[3] =  J[0] / norm;
4842bb1671aSToby Isaac         } else {
4852bb1671aSToby Isaac           PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]);
4862bb1671aSToby Isaac           PetscReal x = J[0] * inorm;
4872bb1671aSToby Isaac           PetscReal y = J[3] * inorm;
4882bb1671aSToby Isaac           PetscReal z = J[6] * inorm;
4892bb1671aSToby Isaac 
4902bb1671aSToby Isaac           if (x > 0.) {
4912bb1671aSToby Isaac             PetscReal inv1pX   = 1./ (1. + x);
4922bb1671aSToby Isaac 
4932bb1671aSToby Isaac             J[1] = -y;              J[2] = -z;
4942bb1671aSToby Isaac             J[4] = 1. - y*y*inv1pX; J[5] =     -y*z*inv1pX;
4952bb1671aSToby Isaac             J[7] =     -y*z*inv1pX; J[8] = 1. - z*z*inv1pX;
4962bb1671aSToby Isaac           } else {
4972bb1671aSToby Isaac             PetscReal inv1mX   = 1./ (1. - x);
4982bb1671aSToby Isaac 
4992bb1671aSToby Isaac             J[1] = z;               J[2] = y;
5002bb1671aSToby Isaac             J[4] =     -y*z*inv1mX; J[5] = 1. - y*y*inv1mX;
5012bb1671aSToby Isaac             J[7] = 1. - z*z*inv1mX; J[8] =     -y*z*inv1mX;
5022bb1671aSToby Isaac           }
5032bb1671aSToby Isaac         }
5042bb1671aSToby Isaac         break;
5052bb1671aSToby Isaac       case 2:
5062bb1671aSToby Isaac         {
5072bb1671aSToby Isaac           PetscReal inorm;
5082bb1671aSToby Isaac 
5092bb1671aSToby Isaac           J[2] = J[3] * J[7] - J[6] * J[4];
5102bb1671aSToby Isaac           J[5] = J[6] * J[1] - J[0] * J[7];
5112bb1671aSToby Isaac           J[8] = J[0] * J[4] - J[3] * J[1];
5122bb1671aSToby Isaac 
5132bb1671aSToby Isaac           inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]);
5142bb1671aSToby Isaac 
5152bb1671aSToby Isaac           J[2] *= inorm;
5162bb1671aSToby Isaac           J[5] *= inorm;
5172bb1671aSToby Isaac           J[8] *= inorm;
5182bb1671aSToby Isaac         }
5192bb1671aSToby Isaac         break;
5202bb1671aSToby Isaac       }
5212bb1671aSToby Isaac       for (j = 0; j < dE*dE; j++) {
5222bb1671aSToby Isaac         g->J[i*dE*dE + j] = J[j];
5232bb1671aSToby Isaac       }
5242bb1671aSToby Isaac     }
5252bb1671aSToby Isaac   }
5262bb1671aSToby Isaac   ierr = PetscFEGeomComplete(g);CHKERRQ(ierr);
527b7260050SToby Isaac   ierr = DMFieldGetDegree(field,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
528b7260050SToby Isaac   g->isAffine = (maxDegree <= 1) ? PETSC_TRUE : PETSC_FALSE;
5290145028aSToby Isaac   if (faceData) {
5300145028aSToby Isaac     if (!field->ops->computeFaceData) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMField implementation does not compute face data\n");
5310145028aSToby Isaac     ierr = (*field->ops->computeFaceData) (field, pointIS, quad, g);CHKERRQ(ierr);
5320145028aSToby Isaac   }
5332bb1671aSToby Isaac   *geom = g;
5342bb1671aSToby Isaac   PetscFunctionReturn(0);
5352bb1671aSToby Isaac }
536