xref: /petsc/src/dm/field/interface/dmfield.c (revision 2bb1671ad356c2a323b47ef3bea96e0958aa2d7a)
13da551e6SToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
23da551e6SToby Isaac 
33da551e6SToby Isaac const char *const DMFieldContinuities[] = {
43da551e6SToby Isaac   "VERTEX",
53da551e6SToby Isaac   "EDGE",
63da551e6SToby Isaac   "FACET",
73da551e6SToby Isaac   "CELL",
83da551e6SToby Isaac   0
93da551e6SToby Isaac };
103da551e6SToby Isaac 
113da551e6SToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate(DM dm,PetscInt numComponents,DMFieldContinuity continuity,DMField *field)
123da551e6SToby Isaac {
133da551e6SToby Isaac   PetscErrorCode ierr;
143da551e6SToby Isaac   DMField        b;
153da551e6SToby Isaac 
163da551e6SToby Isaac   PetscFunctionBegin;
173da551e6SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
183da551e6SToby Isaac   PetscValidPointer(field,2);
193da551e6SToby Isaac   ierr = DMFieldInitializePackage();CHKERRQ(ierr);
203da551e6SToby Isaac 
213da551e6SToby Isaac   ierr = PetscHeaderCreate(b,DMFIELD_CLASSID,"DMField","Field over DM","DM",PetscObjectComm((PetscObject)dm),DMFieldDestroy,DMFieldView);CHKERRQ(ierr);
223da551e6SToby Isaac   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
233da551e6SToby Isaac   b->dm = dm;
243da551e6SToby Isaac   b->continuity = continuity;
253da551e6SToby Isaac   b->numComponents = numComponents;
263da551e6SToby Isaac   *field = b;
273da551e6SToby Isaac   PetscFunctionReturn(0);
283da551e6SToby Isaac }
293da551e6SToby Isaac 
303da551e6SToby Isaac /*@
313da551e6SToby Isaac    DMFieldDestroy - destroy a DMField
323da551e6SToby Isaac 
333da551e6SToby Isaac    Collective
343da551e6SToby Isaac 
353da551e6SToby Isaac    Input Arguments:
363da551e6SToby Isaac .  field - address of DMField
373da551e6SToby Isaac 
383da551e6SToby Isaac    Level: advanced
393da551e6SToby Isaac 
403da551e6SToby Isaac .seealso: DMFieldCreate()
413da551e6SToby Isaac @*/
423da551e6SToby Isaac PetscErrorCode DMFieldDestroy(DMField *field)
433da551e6SToby Isaac {
443da551e6SToby Isaac   PetscErrorCode ierr;
453da551e6SToby Isaac 
463da551e6SToby Isaac   PetscFunctionBegin;
473da551e6SToby Isaac   if (!*field) PetscFunctionReturn(0);
483da551e6SToby Isaac   PetscValidHeaderSpecific((*field),DMFIELD_CLASSID,1);
493da551e6SToby Isaac   if (--((PetscObject)(*field))->refct > 0) {*field = 0; PetscFunctionReturn(0);}
503da551e6SToby Isaac   if ((*field)->ops->destroy) {ierr = (*(*field)->ops->destroy)(*field);CHKERRQ(ierr);}
513da551e6SToby Isaac   ierr = DMDestroy(&((*field)->dm));CHKERRQ(ierr);
523da551e6SToby Isaac   ierr = PetscHeaderDestroy(field);CHKERRQ(ierr);
533da551e6SToby Isaac   PetscFunctionReturn(0);
543da551e6SToby Isaac }
553da551e6SToby Isaac 
563da551e6SToby Isaac /*@C
573da551e6SToby Isaac    DMFieldView - view a DMField
583da551e6SToby Isaac 
593da551e6SToby Isaac    Collective
603da551e6SToby Isaac 
613da551e6SToby Isaac    Input Arguments:
623da551e6SToby Isaac +  field - DMField
633da551e6SToby Isaac -  viewer - viewer to display field, for example PETSC_VIEWER_STDOUT_WORLD
643da551e6SToby Isaac 
653da551e6SToby Isaac    Level: advanced
663da551e6SToby Isaac 
673da551e6SToby Isaac .seealso: DMFieldCreate()
683da551e6SToby Isaac @*/
693da551e6SToby Isaac PetscErrorCode DMFieldView(DMField field,PetscViewer viewer)
703da551e6SToby Isaac {
713da551e6SToby Isaac   PetscErrorCode    ierr;
723da551e6SToby Isaac   PetscBool         iascii;
733da551e6SToby Isaac 
743da551e6SToby Isaac   PetscFunctionBegin;
753da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
763da551e6SToby Isaac   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)field),&viewer);CHKERRQ(ierr);}
773da551e6SToby Isaac   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
783da551e6SToby Isaac   PetscCheckSameComm(field,1,viewer,2);
793da551e6SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
803da551e6SToby Isaac   if (iascii) {
813da551e6SToby Isaac     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)field,viewer);CHKERRQ(ierr);
823da551e6SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
833da551e6SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer,"%D components\n",field->numComponents);CHKERRQ(ierr);
843da551e6SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer,"%s continuity\n",DMFieldContinuities[field->continuity]);CHKERRQ(ierr);
853da551e6SToby Isaac     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);
863da551e6SToby Isaac     ierr = DMView(field->dm,viewer);CHKERRQ(ierr);
873da551e6SToby Isaac     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
883da551e6SToby Isaac   }
893da551e6SToby Isaac   if (field->ops->view) {ierr = (*field->ops->view)(field,viewer);CHKERRQ(ierr);}
903da551e6SToby Isaac   if (iascii) {
913da551e6SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
923da551e6SToby Isaac   }
933da551e6SToby Isaac   PetscFunctionReturn(0);
943da551e6SToby Isaac }
953da551e6SToby Isaac 
963da551e6SToby Isaac /*@C
973da551e6SToby Isaac    DMFieldSetType - set the DMField implementation
983da551e6SToby Isaac 
993da551e6SToby Isaac    Collective on DMField
1003da551e6SToby Isaac 
1013da551e6SToby Isaac    Input Parameters:
1023da551e6SToby Isaac +  field - the DMField context
1033da551e6SToby Isaac -  type - a known method
1043da551e6SToby Isaac 
1053da551e6SToby Isaac    Notes:
1063da551e6SToby Isaac    See "include/petscvec.h" for available methods (for instance)
1073da551e6SToby Isaac +    DMFIELDDA    - a field defined only by its values at the corners of a DMDA
1083da551e6SToby Isaac .    DMFIELDDS    - a field defined by a discretization over a mesh set with DMSetField()
1093da551e6SToby Isaac -    DMFIELDSHELL - a field defined by arbitrary callbacks
1103da551e6SToby Isaac 
1113da551e6SToby Isaac   Level: advanced
1123da551e6SToby Isaac 
1133da551e6SToby Isaac .keywords: DMField, set, type
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 .keywords: DMField, get, type, name
1563da551e6SToby Isaac .seealso: DMFieldSetType()
1573da551e6SToby Isaac @*/
1583da551e6SToby Isaac PetscErrorCode  DMFieldGetType(DMField field, DMFieldType *type)
1593da551e6SToby Isaac {
1603da551e6SToby Isaac   PetscErrorCode ierr;
1613da551e6SToby Isaac 
1623da551e6SToby Isaac   PetscFunctionBegin;
1633da551e6SToby Isaac   PetscValidHeaderSpecific(field, VEC_TAGGER_CLASSID,1);
1643da551e6SToby Isaac   PetscValidPointer(type,2);
1653da551e6SToby Isaac   ierr = DMFieldRegisterAll();CHKERRQ(ierr);
1663da551e6SToby Isaac   *type = ((PetscObject)field)->type_name;
1673da551e6SToby Isaac   PetscFunctionReturn(0);
1683da551e6SToby Isaac }
1693da551e6SToby Isaac 
1703da551e6SToby Isaac PetscErrorCode DMFieldGetNumComponents(DMField field, PetscInt *nc)
1713da551e6SToby Isaac {
1723da551e6SToby Isaac   PetscFunctionBegin;
1733da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1743da551e6SToby Isaac   PetscValidIntPointer(nc,2);
1753da551e6SToby Isaac   *nc = field->numComponents;
1763da551e6SToby Isaac   PetscFunctionReturn(0);
1773da551e6SToby Isaac }
1783da551e6SToby Isaac 
1793da551e6SToby Isaac PetscErrorCode DMFieldGetDM(DMField field, DM *dm)
1803da551e6SToby Isaac {
1813da551e6SToby Isaac   PetscFunctionBegin;
1823da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1833da551e6SToby Isaac   PetscValidPointer(dm,2);
1843da551e6SToby Isaac   *dm = field->dm;
1853da551e6SToby Isaac   PetscFunctionReturn(0);
1863da551e6SToby Isaac }
1873da551e6SToby Isaac 
1883da551e6SToby Isaac PetscErrorCode DMFieldEvaluate(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
1893da551e6SToby Isaac {
1903da551e6SToby Isaac   PetscErrorCode ierr;
1913da551e6SToby Isaac 
1923da551e6SToby Isaac   PetscFunctionBegin;
1933da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
1943da551e6SToby Isaac   PetscValidHeaderSpecific(points,VEC_CLASSID,2);
1953da551e6SToby Isaac   if (B) PetscValidPointer(B,3);
1963da551e6SToby Isaac   if (D) PetscValidPointer(D,4);
1973da551e6SToby Isaac   if (H) PetscValidPointer(H,5);
1983da551e6SToby Isaac   if (field->ops->evaluate) {
1993da551e6SToby Isaac     ierr = (*field->ops->evaluate) (field, points, datatype, B, D, H);CHKERRQ(ierr);
2003da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
2013da551e6SToby Isaac   PetscFunctionReturn(0);
2023da551e6SToby Isaac }
2033da551e6SToby Isaac 
2043da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFE(DMField field, IS cellIS, PetscQuadrature points, PetscDataType datatype, void *B, void *D, void *H)
2053da551e6SToby Isaac {
2063da551e6SToby Isaac   PetscErrorCode ierr;
2073da551e6SToby Isaac 
2083da551e6SToby Isaac   PetscFunctionBegin;
2093da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2103da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
2113da551e6SToby Isaac   PetscValidHeader(points,3);
2123da551e6SToby Isaac   if (B) PetscValidPointer(B,4);
2133da551e6SToby Isaac   if (D) PetscValidPointer(D,5);
2143da551e6SToby Isaac   if (H) PetscValidPointer(H,6);
2153da551e6SToby Isaac   if (field->ops->evaluateFE) {
2163da551e6SToby Isaac     ierr = (*field->ops->evaluateFE) (field, cellIS, points, datatype, B, D, H);CHKERRQ(ierr);
2173da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
2183da551e6SToby Isaac   PetscFunctionReturn(0);
2193da551e6SToby Isaac }
2203da551e6SToby Isaac 
2213da551e6SToby Isaac PetscErrorCode DMFieldEvaluateFV(DMField field, IS cellIS, PetscDataType datatype, void *B, void *D, void *H)
2223da551e6SToby Isaac {
2233da551e6SToby Isaac   PetscErrorCode ierr;
2243da551e6SToby Isaac 
2253da551e6SToby Isaac   PetscFunctionBegin;
2263da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2273da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
2283da551e6SToby Isaac   if (B) PetscValidPointer(B,3);
2293da551e6SToby Isaac   if (D) PetscValidPointer(D,4);
2303da551e6SToby Isaac   if (H) PetscValidPointer(H,5);
2313da551e6SToby Isaac   if (field->ops->evaluateFV) {
2323da551e6SToby Isaac     ierr = (*field->ops->evaluateFV) (field, cellIS, datatype, B, D, H);CHKERRQ(ierr);
2333da551e6SToby Isaac   } else SETERRQ (PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented for this type");
2343da551e6SToby Isaac   PetscFunctionReturn(0);
2353da551e6SToby Isaac }
2363da551e6SToby Isaac 
2373da551e6SToby Isaac PetscErrorCode DMFieldGetFEInvariance(DMField field, IS cellIS, PetscBool *isConstant, PetscBool *isAffine, PetscBool *isQuadratic)
2383da551e6SToby Isaac {
2393da551e6SToby Isaac   PetscErrorCode ierr;
2403da551e6SToby Isaac 
2413da551e6SToby Isaac   PetscFunctionBegin;
2423da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2433da551e6SToby Isaac   PetscValidHeaderSpecific(cellIS,IS_CLASSID,2);
2443da551e6SToby Isaac   if (isConstant) PetscValidPointer(isConstant,3);
2453da551e6SToby Isaac   if (isAffine) PetscValidPointer(isAffine,4);
2463da551e6SToby Isaac   if (isQuadratic) PetscValidPointer(isQuadratic,5);
2473da551e6SToby Isaac 
2483da551e6SToby Isaac   if (isConstant)  *isConstant  = PETSC_FALSE;
2493da551e6SToby Isaac   if (isAffine)    *isAffine    = PETSC_FALSE;
2503da551e6SToby Isaac   if (isQuadratic) *isQuadratic = PETSC_FALSE;
2513da551e6SToby Isaac 
2523da551e6SToby Isaac   if (field->ops->getFEInvariance) {
2533da551e6SToby Isaac     ierr = (*field->ops->getFEInvariance) (field,cellIS,isConstant,isAffine,isQuadratic);CHKERRQ(ierr);
2543da551e6SToby Isaac   }
2553da551e6SToby Isaac   PetscFunctionReturn(0);
2563da551e6SToby Isaac }
2573da551e6SToby Isaac 
2583da551e6SToby Isaac PetscErrorCode DMFieldCreateDefaultQuadrature(DMField field, IS pointIS, PetscQuadrature *quad)
2593da551e6SToby Isaac {
2603da551e6SToby Isaac   PetscErrorCode ierr;
2613da551e6SToby Isaac 
2623da551e6SToby Isaac   PetscFunctionBegin;
2633da551e6SToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
2643da551e6SToby Isaac   PetscValidHeaderSpecific(pointIS,IS_CLASSID,2);
2653da551e6SToby Isaac   PetscValidPointer(quad,3);
2663da551e6SToby Isaac 
2673da551e6SToby Isaac   *quad = NULL;
2683da551e6SToby Isaac   if (field->ops->createDefaultQuadrature) {
2693da551e6SToby Isaac     ierr = (*field->ops->createDefaultQuadrature)(field, pointIS, quad);CHKERRQ(ierr);
2703da551e6SToby Isaac   }
2713da551e6SToby Isaac   PetscFunctionReturn(0);
2723da551e6SToby Isaac }
2733da551e6SToby Isaac 
274*2bb1671aSToby Isaac PetscErrorCode DMFieldCreateFEGeom(DMField field, IS pointIS, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom)
275*2bb1671aSToby Isaac {
276*2bb1671aSToby Isaac   PetscInt       dim, dE;
277*2bb1671aSToby Isaac   PetscInt       nPoints;
278*2bb1671aSToby Isaac   PetscFEGeom    *g;
279*2bb1671aSToby Isaac   PetscErrorCode ierr;
280*2bb1671aSToby Isaac 
281*2bb1671aSToby Isaac   PetscFunctionBegin;
282*2bb1671aSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
283*2bb1671aSToby Isaac   PetscValidHeaderSpecific(pointIS,IS_CLASSID,2);
284*2bb1671aSToby Isaac   PetscValidHeader(quad,3);
285*2bb1671aSToby Isaac   ierr = ISGetLocalSize(pointIS,&nPoints);CHKERRQ(ierr);
286*2bb1671aSToby Isaac   dE = field->numComponents;
287*2bb1671aSToby Isaac   ierr = PetscFEGeomCreate(quad,nPoints,dE,faceData,&g);CHKERRQ(ierr);
288*2bb1671aSToby Isaac   ierr = DMFieldEvaluateFE(field,pointIS,quad,PETSC_REAL,g->v,g->J,NULL);CHKERRQ(ierr);
289*2bb1671aSToby Isaac   dim = g->dim;
290*2bb1671aSToby Isaac   if (dE > dim) {
291*2bb1671aSToby Isaac     /* space out J and make square Jacobians */
292*2bb1671aSToby Isaac     PetscInt  i, j, k, N = g->numPoints * g->numCells;
293*2bb1671aSToby Isaac 
294*2bb1671aSToby Isaac     for (i = N-1; i >= 0; i--) {
295*2bb1671aSToby Isaac       PetscReal   J[9];
296*2bb1671aSToby Isaac 
297*2bb1671aSToby Isaac       for (j = 0; j < dE; j++) {
298*2bb1671aSToby Isaac         for (k = 0; k < dim; k++) {
299*2bb1671aSToby Isaac           J[j*dE + k] = g->J[i*dE*dim + j*dim + k];
300*2bb1671aSToby Isaac         }
301*2bb1671aSToby Isaac       }
302*2bb1671aSToby Isaac       switch (dim) {
303*2bb1671aSToby Isaac       case 0:
304*2bb1671aSToby Isaac         for (j = 0; j < dE; j++) {
305*2bb1671aSToby Isaac           for (k = 0; k < dE; k++) {
306*2bb1671aSToby Isaac             J[j * dE + k] = (j == k) ? 1. : 0.;
307*2bb1671aSToby Isaac           }
308*2bb1671aSToby Isaac         }
309*2bb1671aSToby Isaac         break;
310*2bb1671aSToby Isaac       case 1:
311*2bb1671aSToby Isaac         if (dE == 2) {
312*2bb1671aSToby Isaac           PetscReal norm = PetscSqrtReal(J[0] * J[0] + J[2] * J[2]);
313*2bb1671aSToby Isaac 
314*2bb1671aSToby Isaac           J[1] = -J[2] / norm;
315*2bb1671aSToby Isaac           J[3] =  J[0] / norm;
316*2bb1671aSToby Isaac         } else {
317*2bb1671aSToby Isaac           PetscReal inorm = 1./PetscSqrtReal(J[0] * J[0] + J[3] * J[3] + J[6] * J[6]);
318*2bb1671aSToby Isaac           PetscReal x = J[0] * inorm;
319*2bb1671aSToby Isaac           PetscReal y = J[3] * inorm;
320*2bb1671aSToby Isaac           PetscReal z = J[6] * inorm;
321*2bb1671aSToby Isaac 
322*2bb1671aSToby Isaac           if (x > 0.) {
323*2bb1671aSToby Isaac             PetscReal inv1pX   = 1./ (1. + x);
324*2bb1671aSToby Isaac 
325*2bb1671aSToby Isaac             J[1] = -y;              J[2] = -z;
326*2bb1671aSToby Isaac             J[4] = 1. - y*y*inv1pX; J[5] =     -y*z*inv1pX;
327*2bb1671aSToby Isaac             J[7] =     -y*z*inv1pX; J[8] = 1. - z*z*inv1pX;
328*2bb1671aSToby Isaac           } else {
329*2bb1671aSToby Isaac             PetscReal inv1mX   = 1./ (1. - x);
330*2bb1671aSToby Isaac 
331*2bb1671aSToby Isaac             J[1] = z;               J[2] = y;
332*2bb1671aSToby Isaac             J[4] =     -y*z*inv1mX; J[5] = 1. - y*y*inv1mX;
333*2bb1671aSToby Isaac             J[7] = 1. - z*z*inv1mX; J[8] =     -y*z*inv1mX;
334*2bb1671aSToby Isaac           }
335*2bb1671aSToby Isaac         }
336*2bb1671aSToby Isaac         break;
337*2bb1671aSToby Isaac       case 2:
338*2bb1671aSToby Isaac         {
339*2bb1671aSToby Isaac           PetscReal inorm;
340*2bb1671aSToby Isaac 
341*2bb1671aSToby Isaac           J[2] = J[3] * J[7] - J[6] * J[4];
342*2bb1671aSToby Isaac           J[5] = J[6] * J[1] - J[0] * J[7];
343*2bb1671aSToby Isaac           J[8] = J[0] * J[4] - J[3] * J[1];
344*2bb1671aSToby Isaac 
345*2bb1671aSToby Isaac           inorm = 1./ PetscSqrtReal(J[2]*J[2] + J[5]*J[5] + J[8]*J[8]);
346*2bb1671aSToby Isaac 
347*2bb1671aSToby Isaac           J[2] *= inorm;
348*2bb1671aSToby Isaac           J[5] *= inorm;
349*2bb1671aSToby Isaac           J[8] *= inorm;
350*2bb1671aSToby Isaac         }
351*2bb1671aSToby Isaac         break;
352*2bb1671aSToby Isaac       }
353*2bb1671aSToby Isaac       for (j = 0; j < dE*dE; j++) {
354*2bb1671aSToby Isaac         g->J[i*dE*dE + j] = J[j];
355*2bb1671aSToby Isaac       }
356*2bb1671aSToby Isaac     }
357*2bb1671aSToby Isaac   }
358*2bb1671aSToby Isaac   ierr = PetscFEGeomComplete(g);CHKERRQ(ierr);
359*2bb1671aSToby Isaac   ierr = DMFieldGetFEInvariance(field,pointIS,NULL,&g->isAffine,NULL);CHKERRQ(ierr);
360*2bb1671aSToby Isaac   *geom = g;
361*2bb1671aSToby Isaac   PetscFunctionReturn(0);
362*2bb1671aSToby Isaac }
363