xref: /petsc/src/dm/field/interface/dmfield.c (revision 0145028a27feb5d014cb6ee82515d46f0128efc1)
13da551e6SToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
2*0145028aSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/
3*0145028aSToby Isaac #include <petscdmplex.h>
43da551e6SToby Isaac 
53da551e6SToby Isaac const char *const DMFieldContinuities[] = {
63da551e6SToby Isaac   "VERTEX",
73da551e6SToby Isaac   "EDGE",
83da551e6SToby Isaac   "FACET",
93da551e6SToby Isaac   "CELL",
103da551e6SToby Isaac   0
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);
513da551e6SToby Isaac   if (--((PetscObject)(*field))->refct > 0) {*field = 0; 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 
1013da551e6SToby Isaac    Collective on DMField
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 .keywords: DMField, set, type
1163da551e6SToby Isaac 
1173da551e6SToby Isaac .seealso: DMFieldType,
1183da551e6SToby Isaac @*/
1193da551e6SToby Isaac PetscErrorCode DMFieldSetType(DMField field,DMFieldType type)
1203da551e6SToby Isaac {
1213da551e6SToby Isaac   PetscErrorCode ierr,(*r)(DMField);
1223da551e6SToby Isaac   PetscBool      match;
1233da551e6SToby Isaac 
1243da551e6SToby Isaac   PetscFunctionBegin;
1253da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1263da551e6SToby Isaac   PetscValidCharPointer(type,2);
1273da551e6SToby Isaac 
1283da551e6SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)field,type,&match);CHKERRQ(ierr);
1293da551e6SToby Isaac   if (match) PetscFunctionReturn(0);
1303da551e6SToby Isaac 
1313da551e6SToby Isaac   ierr = PetscFunctionListFind(DMFieldList,type,&r);CHKERRQ(ierr);
1323da551e6SToby Isaac   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DMField type %s",type);
1333da551e6SToby Isaac   /* Destroy the previous private DMField context */
1343da551e6SToby Isaac   if (field->ops->destroy) {
1353da551e6SToby Isaac     ierr = (*(field)->ops->destroy)(field);CHKERRQ(ierr);
1363da551e6SToby Isaac   }
1373da551e6SToby Isaac   ierr = PetscMemzero(field->ops,sizeof(*field->ops));CHKERRQ(ierr);
1383da551e6SToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject)field,type);CHKERRQ(ierr);
1393da551e6SToby Isaac   field->ops->create = r;
1403da551e6SToby Isaac   ierr = (*r)(field);CHKERRQ(ierr);
1413da551e6SToby Isaac   PetscFunctionReturn(0);
1423da551e6SToby Isaac }
1433da551e6SToby Isaac 
1443da551e6SToby Isaac /*@C
1453da551e6SToby Isaac   DMFieldGetType - Gets the DMField type name (as a string) from the DMField.
1463da551e6SToby Isaac 
1473da551e6SToby Isaac   Not Collective
1483da551e6SToby Isaac 
1493da551e6SToby Isaac   Input Parameter:
1503da551e6SToby Isaac . field  - The DMField context
1513da551e6SToby Isaac 
1523da551e6SToby Isaac   Output Parameter:
1533da551e6SToby Isaac . type - The DMField type name
1543da551e6SToby Isaac 
1553da551e6SToby Isaac   Level: advanced
1563da551e6SToby Isaac 
1573da551e6SToby Isaac .keywords: DMField, get, type, name
1583da551e6SToby Isaac .seealso: DMFieldSetType()
1593da551e6SToby Isaac @*/
1603da551e6SToby Isaac PetscErrorCode  DMFieldGetType(DMField field, DMFieldType *type)
1613da551e6SToby Isaac {
1623da551e6SToby Isaac   PetscErrorCode ierr;
1633da551e6SToby Isaac 
1643da551e6SToby Isaac   PetscFunctionBegin;
1653da551e6SToby Isaac   PetscValidHeaderSpecific(field, VEC_TAGGER_CLASSID,1);
1663da551e6SToby Isaac   PetscValidPointer(type,2);
1673da551e6SToby Isaac   ierr = DMFieldRegisterAll();CHKERRQ(ierr);
1683da551e6SToby Isaac   *type = ((PetscObject)field)->type_name;
1693da551e6SToby Isaac   PetscFunctionReturn(0);
1703da551e6SToby Isaac }
1713da551e6SToby Isaac 
1723da551e6SToby Isaac PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc)
1733da551e6SToby Isaac {
1743da551e6SToby Isaac   PetscFunctionBegin;
1753da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1763da551e6SToby Isaac   PetscValidIntPointer(nc,2);
1773da551e6SToby Isaac   *nc = field->numComponents;
1783da551e6SToby Isaac   PetscFunctionReturn(0);
1793da551e6SToby Isaac }
1803da551e6SToby Isaac 
1813da551e6SToby Isaac PetscErrorCode DMFieldGetDM(DMField field, DM *dm)
1823da551e6SToby Isaac {
1833da551e6SToby Isaac   PetscFunctionBegin;
1843da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1853da551e6SToby Isaac   PetscValidPointer(dm,2);
1863da551e6SToby Isaac   *dm = field->dm;
1873da551e6SToby Isaac   PetscFunctionReturn(0);
1883da551e6SToby Isaac }
1893da551e6SToby Isaac 
1903da551e6SToby Isaac PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
1913da551e6SToby Isaac {
1923da551e6SToby Isaac   PetscErrorCode ierr;
1933da551e6SToby Isaac 
1943da551e6SToby Isaac   PetscFunctionBegin;
1953da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1963da551e6SToby Isaac   PetscValidHeaderSpecific(points,VEC_CLASSID,2);
1973da551e6SToby Isaac   if (B) PetscValidPointer(B,3);
1983da551e6SToby Isaac   if (D) PetscValidPointer(D,4);
1993da551e6SToby Isaac   if (H) PetscValidPointer(H,5);
2003da551e6SToby Isaac   if (field->ops->evaluate) {
2013da551e6SToby Isaac     ierr = (*field->ops->evaluate) (field, points, datatype, B, D, H);CHKERRQ(ierr);
2023da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
2033da551e6SToby Isaac   PetscFunctionReturn(0);
2043da551e6SToby Isaac }
2053da551e6SToby Isaac 
2063da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
2073da551e6SToby Isaac {
2083da551e6SToby Isaac   PetscErrorCode ierr;
2093da551e6SToby Isaac 
2103da551e6SToby Isaac   PetscFunctionBegin;
2113da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2123da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
2133da551e6SToby Isaac   PetscValidHeader(points,3);
2143da551e6SToby Isaac   if (B) PetscValidPointer(B,4);
2153da551e6SToby Isaac   if (D) PetscValidPointer(D,5);
2163da551e6SToby Isaac   if (H) PetscValidPointer(H,6);
2173da551e6SToby Isaac   if (field->ops->evaluateFE) {
2183da551e6SToby Isaac     ierr = (*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H);CHKERRQ(ierr);
2193da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
2203da551e6SToby Isaac   PetscFunctionReturn(0);
2213da551e6SToby Isaac }
2223da551e6SToby Isaac 
2233da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
2243da551e6SToby Isaac {
2253da551e6SToby Isaac   PetscErrorCode ierr;
2263da551e6SToby Isaac 
2273da551e6SToby Isaac   PetscFunctionBegin;
2283da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2293da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
2303da551e6SToby Isaac   if (B) PetscValidPointer(B,3);
2313da551e6SToby Isaac   if (D) PetscValidPointer(D,4);
2323da551e6SToby Isaac   if (H) PetscValidPointer(H,5);
2333da551e6SToby Isaac   if (field->ops->evaluateFV) {
2343da551e6SToby Isaac     ierr = (*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H);CHKERRQ(ierr);
2353da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
2363da551e6SToby Isaac   PetscFunctionReturn(0);
2373da551e6SToby Isaac }
2383da551e6SToby Isaac 
2393da551e6SToby Isaac PetscErrorCode DMFieldGetFEInvariance(DMField field, IS cellIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic)
2403da551e6SToby Isaac {
2413da551e6SToby Isaac   PetscErrorCode ierr;
2423da551e6SToby Isaac 
2433da551e6SToby Isaac   PetscFunctionBegin;
2443da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2453da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
2463da551e6SToby Isaac   if (isConstant) PetscValidPointer(isConstant,3);
2473da551e6SToby Isaac   if (isAffine) PetscValidPointer(isAffine,4);
2483da551e6SToby Isaac   if (isQuadratic) PetscValidPointer(isQuadratic,5);
2493da551e6SToby Isaac 
2503da551e6SToby Isaac   if (isConstant)  *isConstant  = PETSC_FALSE;
2513da551e6SToby Isaac   if (isAffine)    *isAffine    = PETSC_FALSE;
2523da551e6SToby Isaac   if (isQuadratic) *isQuadratic = PETSC_FALSE;
2533da551e6SToby Isaac 
2543da551e6SToby Isaac   if (field->ops->getFEInvariance) {
2553da551e6SToby Isaac     ierr = (*field->ops->getFEInvariance) (field,cellIS,isConstant,isAffine,isQuadratic);CHKERRQ(ierr);
2563da551e6SToby Isaac   }
2573da551e6SToby Isaac   PetscFunctionReturn(0);
2583da551e6SToby Isaac }
2593da551e6SToby Isaac 
2603da551e6SToby Isaac PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad)
2613da551e6SToby Isaac {
2623da551e6SToby Isaac   PetscErrorCode ierr;
2633da551e6SToby Isaac 
2643da551e6SToby Isaac   PetscFunctionBegin;
2653da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2663da551e6SToby Isaac   PetscValidHeaderSpecific(pointIS,IS_CLASSID,2);
2673da551e6SToby Isaac   PetscValidPointer(quad,3);
2683da551e6SToby Isaac 
2693da551e6SToby Isaac   *quad = NULL;
2703da551e6SToby Isaac   if (field->ops->createDefaultQuadrature) {
2713da551e6SToby Isaac     ierr = (*field->ops->createDefaultQuadrature)(field, pointIS, quad);CHKERRQ(ierr);
2723da551e6SToby Isaac   }
2733da551e6SToby Isaac   PetscFunctionReturn(0);
2743da551e6SToby Isaac }
2753da551e6SToby Isaac 
2762bb1671aSToby Isaac PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
2772bb1671aSToby Isaac {
2782bb1671aSToby Isaac   PetscInt       dim, dE;
2792bb1671aSToby Isaac   PetscInt       nPoints;
2802bb1671aSToby Isaac   PetscFEGeom    *g;
2812bb1671aSToby Isaac   PetscErrorCode ierr;
2822bb1671aSToby Isaac 
2832bb1671aSToby Isaac   PetscFunctionBegin;
2842bb1671aSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2852bb1671aSToby Isaac   PetscValidHeaderSpecific(pointIS,IS_CLASSID,2);
2862bb1671aSToby Isaac   PetscValidHeader(quad,3);
2872bb1671aSToby Isaac   ierr = ISGetLocalSize(pointIS,&nPoints);CHKERRQ(ierr);
2882bb1671aSToby Isaac   dE = field->numComponents;
2892bb1671aSToby Isaac   ierr = PetscFEGeomCreate(quad,nPoints,dE,faceData,&g);CHKERRQ(ierr);
2902bb1671aSToby Isaac   ierr = DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL);CHKERRQ(ierr);
2912bb1671aSToby Isaac   dim = g->dim;
2922bb1671aSToby Isaac   if (dE > dim) {
2932bb1671aSToby Isaac     /* space out J and make square Jacobians */
2942bb1671aSToby Isaac     PetscInt  i, j, k, N = g->numPoints * g->numCells;
2952bb1671aSToby Isaac 
2962bb1671aSToby Isaac     for (i = N-1; i >= 0; i--) {
2972bb1671aSToby Isaac       PetscReal   J[9];
2982bb1671aSToby Isaac 
2992bb1671aSToby Isaac       for (j = 0; j < dE; j++) {
3002bb1671aSToby Isaac         for (k = 0; k < dim; k++) {
3012bb1671aSToby Isaac           J[j*dE + k] = g->J[i*dE*dim + j*dim + k];
3022bb1671aSToby Isaac         }
3032bb1671aSToby Isaac       }
3042bb1671aSToby Isaac       switch (dim) {
3052bb1671aSToby Isaac       case 0:
3062bb1671aSToby Isaac         for (j = 0; j < dE; j++) {
3072bb1671aSToby Isaac           for (k = 0; k < dE; k++) {
3082bb1671aSToby Isaac             J[j * dE + k] = (j == k) ? 1. : 0.;
3092bb1671aSToby Isaac           }
3102bb1671aSToby Isaac         }
3112bb1671aSToby Isaac         break;
3122bb1671aSToby Isaac       case 1:
3132bb1671aSToby Isaac         if (dE == 2) {
3142bb1671aSToby Isaac           PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]);
3152bb1671aSToby Isaac 
3162bb1671aSToby Isaac           J[1] = -J[2] / norm;
3172bb1671aSToby Isaac           J[3] =  J[0] / norm;
3182bb1671aSToby Isaac         } else {
3192bb1671aSToby Isaac           PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]);
3202bb1671aSToby Isaac           PetscReal x = J[0] * inorm;
3212bb1671aSToby Isaac           PetscReal y = J[3] * inorm;
3222bb1671aSToby Isaac           PetscReal z = J[6] * inorm;
3232bb1671aSToby Isaac 
3242bb1671aSToby Isaac           if (x > 0.) {
3252bb1671aSToby Isaac             PetscReal inv1pX   = 1./ (1. + x);
3262bb1671aSToby Isaac 
3272bb1671aSToby Isaac             J[1] = -y;              J[2] = -z;
3282bb1671aSToby Isaac             J[4] = 1. - y*y*inv1pX; J[5] =     -y*z*inv1pX;
3292bb1671aSToby Isaac             J[7] =     -y*z*inv1pX; J[8] = 1. - z*z*inv1pX;
3302bb1671aSToby Isaac           } else {
3312bb1671aSToby Isaac             PetscReal inv1mX   = 1./ (1. - x);
3322bb1671aSToby Isaac 
3332bb1671aSToby Isaac             J[1] = z;               J[2] = y;
3342bb1671aSToby Isaac             J[4] =     -y*z*inv1mX; J[5] = 1. - y*y*inv1mX;
3352bb1671aSToby Isaac             J[7] = 1. - z*z*inv1mX; J[8] =     -y*z*inv1mX;
3362bb1671aSToby Isaac           }
3372bb1671aSToby Isaac         }
3382bb1671aSToby Isaac         break;
3392bb1671aSToby Isaac       case 2:
3402bb1671aSToby Isaac         {
3412bb1671aSToby Isaac           PetscReal inorm;
3422bb1671aSToby Isaac 
3432bb1671aSToby Isaac           J[2] = J[3] * J[7] - J[6] * J[4];
3442bb1671aSToby Isaac           J[5] = J[6] * J[1] - J[0] * J[7];
3452bb1671aSToby Isaac           J[8] = J[0] * J[4] - J[3] * J[1];
3462bb1671aSToby Isaac 
3472bb1671aSToby Isaac           inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]);
3482bb1671aSToby Isaac 
3492bb1671aSToby Isaac           J[2] *= inorm;
3502bb1671aSToby Isaac           J[5] *= inorm;
3512bb1671aSToby Isaac           J[8] *= inorm;
3522bb1671aSToby Isaac         }
3532bb1671aSToby Isaac         break;
3542bb1671aSToby Isaac       }
3552bb1671aSToby Isaac       for (j = 0; j < dE*dE; j++) {
3562bb1671aSToby Isaac         g->J[i*dE*dE + j] = J[j];
3572bb1671aSToby Isaac       }
3582bb1671aSToby Isaac     }
3592bb1671aSToby Isaac   }
3602bb1671aSToby Isaac   ierr = PetscFEGeomComplete(g);CHKERRQ(ierr);
3612bb1671aSToby Isaac   ierr = DMFieldGetFEInvariance(field,pointIS,NULL,&g->isAffine,NULL);CHKERRQ(ierr);
362*0145028aSToby Isaac   if (faceData) {
363*0145028aSToby Isaac     if (!field->ops->computeFaceData) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMField implementation does not compute face data\n");
364*0145028aSToby Isaac     ierr = (*field->ops->computeFaceData) (field, pointIS, quad, g);CHKERRQ(ierr);
365*0145028aSToby Isaac   }
3662bb1671aSToby Isaac   *geom = g;
3672bb1671aSToby Isaac   PetscFunctionReturn(0);
3682bb1671aSToby Isaac }
369