xref: /petsc/src/dm/dt/fe/interface/fegeom.c (revision 5fedec97950c19de564efaecd0f125b1a6cb2b20)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
32b99622eSMatthew G. Knepley /*@C
42b99622eSMatthew G. Knepley   PetscFEGeomCreate - Create a PetscFEGeom object to manage geometry for a group of cells
52b99622eSMatthew G. Knepley 
62b99622eSMatthew G. Knepley   Input Parameters:
72b99622eSMatthew G. Knepley + quad     - A PetscQuadrature determining the tabulation
82b99622eSMatthew G. Knepley . numCells - The number of cells in the group
92b99622eSMatthew G. Knepley . dimEmbed - The coordinate dimension
102b99622eSMatthew G. Knepley - faceData - Flag to construct geometry data for the faces
112b99622eSMatthew G. Knepley 
122b99622eSMatthew G. Knepley   Output Parameter:
132b99622eSMatthew G. Knepley . geom     - The PetscFEGeom object
142b99622eSMatthew G. Knepley 
152b99622eSMatthew G. Knepley   Level: beginner
162b99622eSMatthew G. Knepley 
172b99622eSMatthew G. Knepley .seealso: PetscFEGeomDestroy(), PetscFEGeomComplete()
182b99622eSMatthew G. Knepley @*/
1920cf1dd8SToby Isaac PetscErrorCode PetscFEGeomCreate(PetscQuadrature quad, PetscInt numCells, PetscInt dimEmbed, PetscBool faceData, PetscFEGeom **geom)
2020cf1dd8SToby Isaac {
2120cf1dd8SToby Isaac   PetscFEGeom     *g;
227489efa5SMatthew G. Knepley   PetscInt        dim, Nq, N;
2320cf1dd8SToby Isaac   const PetscReal *p;
2420cf1dd8SToby Isaac   PetscErrorCode  ierr;
2520cf1dd8SToby Isaac 
2620cf1dd8SToby Isaac   PetscFunctionBegin;
277489efa5SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad,&dim,NULL,&Nq,&p,NULL);CHKERRQ(ierr);
2820cf1dd8SToby Isaac   ierr = PetscNew(&g);CHKERRQ(ierr);
2920cf1dd8SToby Isaac   g->xi         = p;
3020cf1dd8SToby Isaac   g->numCells   = numCells;
317489efa5SMatthew G. Knepley   g->numPoints  = Nq;
3220cf1dd8SToby Isaac   g->dim        = dim;
3320cf1dd8SToby Isaac   g->dimEmbed   = dimEmbed;
34*5fedec97SMatthew G. Knepley   g->isCohesive = PETSC_FALSE;
357489efa5SMatthew G. Knepley   N = numCells * Nq;
3620cf1dd8SToby Isaac   ierr = PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ);CHKERRQ(ierr);
3720cf1dd8SToby Isaac   if (faceData) {
3841418987SMatthew G. Knepley     ierr = PetscCalloc2(numCells, &g->face, N * dimEmbed, &g->n);CHKERRQ(ierr);
3941418987SMatthew G. Knepley     ierr = PetscCalloc6(N * dimEmbed * dimEmbed, &(g->suppJ[0]),    N * dimEmbed * dimEmbed, &(g->suppJ[1]),
4041418987SMatthew G. Knepley                         N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]),
4141418987SMatthew G. Knepley                         N,                       &(g->suppDetJ[0]), N,                       &(g->suppDetJ[1]));CHKERRQ(ierr);
4220cf1dd8SToby Isaac   }
434e2f966aSMatthew G. Knepley   ierr = PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ);CHKERRQ(ierr);
4420cf1dd8SToby Isaac   *geom = g;
4520cf1dd8SToby Isaac   PetscFunctionReturn(0);
4620cf1dd8SToby Isaac }
4720cf1dd8SToby Isaac 
482b99622eSMatthew G. Knepley /*@C
492b99622eSMatthew G. Knepley   PetscFEGeomDestroy - Destroy a PetscFEGeom object
502b99622eSMatthew G. Knepley 
512b99622eSMatthew G. Knepley   Input Parameter:
522b99622eSMatthew G. Knepley . geom - PetscFEGeom object
532b99622eSMatthew G. Knepley 
542b99622eSMatthew G. Knepley   Level: beginner
552b99622eSMatthew G. Knepley 
562b99622eSMatthew G. Knepley .seealso: PetscFEGeomCreate()
572b99622eSMatthew G. Knepley @*/
5820cf1dd8SToby Isaac PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom)
5920cf1dd8SToby Isaac {
6020cf1dd8SToby Isaac   PetscErrorCode ierr;
6120cf1dd8SToby Isaac 
6220cf1dd8SToby Isaac   PetscFunctionBegin;
6320cf1dd8SToby Isaac   if (!*geom) PetscFunctionReturn(0);
6420cf1dd8SToby Isaac   ierr = PetscFree3((*geom)->v,(*geom)->J,(*geom)->detJ);CHKERRQ(ierr);
6520cf1dd8SToby Isaac   ierr = PetscFree((*geom)->invJ);CHKERRQ(ierr);
6641418987SMatthew G. Knepley   ierr = PetscFree2((*geom)->face,(*geom)->n);CHKERRQ(ierr);
6741418987SMatthew G. Knepley   ierr = PetscFree6((*geom)->suppJ[0],(*geom)->suppJ[1],(*geom)->suppInvJ[0],(*geom)->suppInvJ[1],(*geom)->suppDetJ[0],(*geom)->suppDetJ[1]);CHKERRQ(ierr);
6820cf1dd8SToby Isaac   ierr = PetscFree(*geom);CHKERRQ(ierr);
6920cf1dd8SToby Isaac   PetscFunctionReturn(0);
7020cf1dd8SToby Isaac }
7120cf1dd8SToby Isaac 
722b99622eSMatthew G. Knepley /*@C
732b99622eSMatthew G. Knepley   PetscFEGeomGetChunk - Get a chunk of cells in the group as a PetscFEGeom
742b99622eSMatthew G. Knepley 
752b99622eSMatthew G. Knepley   Input Parameters:
762b99622eSMatthew G. Knepley + geom   - PetscFEGeom object
772b99622eSMatthew G. Knepley . cStart - The first cell in the chunk
782b99622eSMatthew G. Knepley - cEnd   - The first cell not in the chunk
792b99622eSMatthew G. Knepley 
802b99622eSMatthew G. Knepley   Output Parameter:
812b99622eSMatthew G. Knepley . chunkGeom - The chunk of cells
822b99622eSMatthew G. Knepley 
832b99622eSMatthew G. Knepley   Level: intermediate
842b99622eSMatthew G. Knepley 
852b99622eSMatthew G. Knepley .seealso: PetscFEGeomRestoreChunk(), PetscFEGeomCreate()
862b99622eSMatthew G. Knepley @*/
8720cf1dd8SToby Isaac PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
8820cf1dd8SToby Isaac {
8920cf1dd8SToby Isaac   PetscInt       Nq;
9020cf1dd8SToby Isaac   PetscInt       dE;
9120cf1dd8SToby Isaac   PetscErrorCode ierr;
9220cf1dd8SToby Isaac 
9320cf1dd8SToby Isaac   PetscFunctionBegin;
9420cf1dd8SToby Isaac   PetscValidPointer(geom,1);
95064a246eSJacob Faibussowitsch   PetscValidPointer(chunkGeom,4);
9620cf1dd8SToby Isaac   if (!(*chunkGeom)) {
9720cf1dd8SToby Isaac     ierr = PetscNew(chunkGeom);CHKERRQ(ierr);
9820cf1dd8SToby Isaac   }
9920cf1dd8SToby Isaac   Nq = geom->numPoints;
10020cf1dd8SToby Isaac   dE= geom->dimEmbed;
10120cf1dd8SToby Isaac   (*chunkGeom)->dim = geom->dim;
10220cf1dd8SToby Isaac   (*chunkGeom)->dimEmbed = geom->dimEmbed;
10320cf1dd8SToby Isaac   (*chunkGeom)->numPoints = geom->numPoints;
10420cf1dd8SToby Isaac   (*chunkGeom)->numCells = cEnd - cStart;
10520cf1dd8SToby Isaac   (*chunkGeom)->xi = geom->xi;
10620cf1dd8SToby Isaac   (*chunkGeom)->v = &geom->v[Nq*dE*cStart];
10720cf1dd8SToby Isaac   (*chunkGeom)->J = &geom->J[Nq*dE*dE*cStart];
10820cf1dd8SToby Isaac   (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq*dE*dE*cStart] : NULL;
10920cf1dd8SToby Isaac   (*chunkGeom)->detJ = &geom->detJ[Nq*cStart];
11020cf1dd8SToby Isaac   (*chunkGeom)->n = geom->n ? &geom->n[Nq*dE*cStart] : NULL;
11120cf1dd8SToby Isaac   (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL;
11241418987SMatthew G. Knepley   (*chunkGeom)->suppJ[0]    = geom->suppJ[0]    ? &geom->suppJ[0][Nq*dE*dE*cStart]    : NULL;
11341418987SMatthew G. Knepley   (*chunkGeom)->suppJ[1]    = geom->suppJ[1]    ? &geom->suppJ[1][Nq*dE*dE*cStart]    : NULL;
11420cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq*dE*dE*cStart] : NULL;
11520cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq*dE*dE*cStart] : NULL;
11641418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[0] = geom->suppDetJ[0] ? &geom->suppDetJ[0][Nq*cStart]       : NULL;
11741418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[1] = geom->suppDetJ[1] ? &geom->suppDetJ[1][Nq*cStart]       : NULL;
11820cf1dd8SToby Isaac   (*chunkGeom)->isAffine = geom->isAffine;
11920cf1dd8SToby Isaac   PetscFunctionReturn(0);
12020cf1dd8SToby Isaac }
12120cf1dd8SToby Isaac 
1222b99622eSMatthew G. Knepley /*@C
1232b99622eSMatthew G. Knepley   PetscFEGeomRestoreChunk - Restore the chunk
1242b99622eSMatthew G. Knepley 
1252b99622eSMatthew G. Knepley   Input Parameters:
1262b99622eSMatthew G. Knepley + geom      - PetscFEGeom object
1272b99622eSMatthew G. Knepley . cStart    - The first cell in the chunk
1282b99622eSMatthew G. Knepley . cEnd      - The first cell not in the chunk
1292b99622eSMatthew G. Knepley - chunkGeom - The chunk of cells
1302b99622eSMatthew G. Knepley 
1312b99622eSMatthew G. Knepley   Level: intermediate
1322b99622eSMatthew G. Knepley 
1332b99622eSMatthew G. Knepley .seealso: PetscFEGeomGetChunk(), PetscFEGeomCreate()
1342b99622eSMatthew G. Knepley @*/
13520cf1dd8SToby Isaac PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
13620cf1dd8SToby Isaac {
13720cf1dd8SToby Isaac   PetscErrorCode ierr;
13820cf1dd8SToby Isaac 
13920cf1dd8SToby Isaac   PetscFunctionBegin;
14020cf1dd8SToby Isaac   ierr = PetscFree(*chunkGeom);CHKERRQ(ierr);
14120cf1dd8SToby Isaac   PetscFunctionReturn(0);
14220cf1dd8SToby Isaac }
14320cf1dd8SToby Isaac 
1446587ee25SMatthew G. Knepley /*@C
1456587ee25SMatthew G. Knepley   PetscFEGeomGetPoint - Get the geometry for cell c at point p as a PetscFEGeom
1466587ee25SMatthew G. Knepley 
1476587ee25SMatthew G. Knepley   Input Parameters:
1486587ee25SMatthew G. Knepley + geom    - PetscFEGeom object
1496587ee25SMatthew G. Knepley . c       - The cell
1506587ee25SMatthew G. Knepley . p       - The point
1516587ee25SMatthew G. Knepley - pcoords - The reference coordinates of point p, or NULL
1526587ee25SMatthew G. Knepley 
1536587ee25SMatthew G. Knepley   Output Parameter:
1546587ee25SMatthew G. Knepley . pgeom - The geometry of cell c at point p
1556587ee25SMatthew G. Knepley 
1566587ee25SMatthew G. Knepley   Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom,
1576587ee25SMatthew G. Knepley   nothing needs to be done with it afterwards.
1586587ee25SMatthew G. Knepley 
1596587ee25SMatthew G. Knepley   In the affine case, pgeom must have storage for the integration point coordinates in pgeom->v if pcoords is passed in.
1606587ee25SMatthew G. Knepley 
1616587ee25SMatthew G. Knepley   Level: intermediate
1626587ee25SMatthew G. Knepley 
1636587ee25SMatthew G. Knepley .seealso: PetscFEGeomRestoreChunk(), PetscFEGeomCreate()
1646587ee25SMatthew G. Knepley @*/
1656587ee25SMatthew G. Knepley PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, const PetscReal pcoords[], PetscFEGeom *pgeom)
1666587ee25SMatthew G. Knepley {
1676587ee25SMatthew G. Knepley   const PetscInt dim = geom->dim;
1686587ee25SMatthew G. Knepley   const PetscInt dE  = geom->dimEmbed;
1696587ee25SMatthew G. Knepley   const PetscInt Np  = geom->numPoints;
1706587ee25SMatthew G. Knepley 
1716587ee25SMatthew G. Knepley   PetscFunctionBeginHot;
1726587ee25SMatthew G. Knepley   pgeom->dim      = dim;
1736587ee25SMatthew G. Knepley   pgeom->dimEmbed = dE;
1746587ee25SMatthew G. Knepley   //pgeom->isAffine = geom->isAffine;
1756587ee25SMatthew G. Knepley   if (geom->isAffine) {
1766587ee25SMatthew G. Knepley     if (!p) {
1776587ee25SMatthew G. Knepley       pgeom->xi   = geom->xi;
1786587ee25SMatthew G. Knepley       pgeom->J    = &geom->J[c*Np*dE*dE];
1796587ee25SMatthew G. Knepley       pgeom->invJ = &geom->invJ[c*Np*dE*dE];
1806587ee25SMatthew G. Knepley       pgeom->detJ = &geom->detJ[c*Np];
1816587ee25SMatthew G. Knepley       pgeom->n    = geom->n ? &geom->n[c*Np*dE] : NULL;
1826587ee25SMatthew G. Knepley     }
1836587ee25SMatthew G. Knepley     if (pcoords) {CoordinatesRefToReal(dE, dim, pgeom->xi, &geom->v[c*Np*dE], pgeom->J, pcoords, pgeom->v);}
1846587ee25SMatthew G. Knepley   } else {
1856587ee25SMatthew G. Knepley     pgeom->v    = &geom->v[(c*Np+p)*dE];
1866587ee25SMatthew G. Knepley     pgeom->J    = &geom->J[(c*Np+p)*dE*dE];
1876587ee25SMatthew G. Knepley     pgeom->invJ = &geom->invJ[(c*Np+p)*dE*dE];
1886587ee25SMatthew G. Knepley     pgeom->detJ = &geom->detJ[c*Np+p];
1896587ee25SMatthew G. Knepley     pgeom->n    = geom->n ? &geom->n[(c*Np+p)*dE] : NULL;
1906587ee25SMatthew G. Knepley   }
1916587ee25SMatthew G. Knepley   PetscFunctionReturn(0);
1926587ee25SMatthew G. Knepley }
1936587ee25SMatthew G. Knepley 
1946587ee25SMatthew G. Knepley /*@C
1956587ee25SMatthew G. Knepley   PetscFEGeomGetCellPoint - Get the cell geometry for face f at point p as a PetscFEGeom
1966587ee25SMatthew G. Knepley 
1976587ee25SMatthew G. Knepley   Input Parameters:
1986587ee25SMatthew G. Knepley + geom    - PetscFEGeom object
1996587ee25SMatthew G. Knepley . f       - The face
2006587ee25SMatthew G. Knepley - p       - The point
2016587ee25SMatthew G. Knepley 
2026587ee25SMatthew G. Knepley   Output Parameter:
2036587ee25SMatthew G. Knepley . pgeom - The cell geometry of face f at point p
2046587ee25SMatthew G. Knepley 
2056587ee25SMatthew G. Knepley   Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom,
2066587ee25SMatthew G. Knepley   nothing needs to be done with it afterwards.
2076587ee25SMatthew G. Knepley 
2086587ee25SMatthew G. Knepley   Level: intermediate
2096587ee25SMatthew G. Knepley 
2106587ee25SMatthew G. Knepley .seealso: PetscFEGeomRestoreChunk(), PetscFEGeomCreate()
2116587ee25SMatthew G. Knepley @*/
2126587ee25SMatthew G. Knepley PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, PetscFEGeom *pgeom)
2136587ee25SMatthew G. Knepley {
214*5fedec97SMatthew G. Knepley   const PetscBool bd  = geom->dimEmbed > geom->dim && !geom->isCohesive ? PETSC_TRUE : PETSC_FALSE;
2156587ee25SMatthew G. Knepley   const PetscInt  dim = bd ? geom->dimEmbed : geom->dim;
2166587ee25SMatthew G. Knepley   const PetscInt  dE  = geom->dimEmbed;
2176587ee25SMatthew G. Knepley   const PetscInt  Np  = geom->numPoints;
2186587ee25SMatthew G. Knepley 
2196587ee25SMatthew G. Knepley   PetscFunctionBeginHot;
2206587ee25SMatthew G. Knepley   pgeom->dim      = dim;
2216587ee25SMatthew G. Knepley   pgeom->dimEmbed = dE;
2226587ee25SMatthew G. Knepley   //pgeom->isAffine = geom->isAffine;
2236587ee25SMatthew G. Knepley   if (geom->isAffine) {
2246587ee25SMatthew G. Knepley     if (!p) {
2256587ee25SMatthew G. Knepley       if (bd) {
2266587ee25SMatthew G. Knepley         pgeom->J     = &geom->suppJ[0][c*Np*dE*dE];
2276587ee25SMatthew G. Knepley         pgeom->invJ  = &geom->suppInvJ[0][c*Np*dE*dE];
2286587ee25SMatthew G. Knepley         pgeom->detJ  = &geom->suppDetJ[0][c*Np];
2296587ee25SMatthew G. Knepley       } else {
2306587ee25SMatthew G. Knepley         pgeom->J    = &geom->J[c*Np*dE*dE];
2316587ee25SMatthew G. Knepley         pgeom->invJ = &geom->invJ[c*Np*dE*dE];
2326587ee25SMatthew G. Knepley         pgeom->detJ = &geom->detJ[c*Np];
2336587ee25SMatthew G. Knepley       }
2346587ee25SMatthew G. Knepley     }
2356587ee25SMatthew G. Knepley   } else {
2366587ee25SMatthew G. Knepley     if (bd) {
2376587ee25SMatthew G. Knepley       pgeom->J     = &geom->suppJ[0][(c*Np+p)*dE*dE];
2386587ee25SMatthew G. Knepley       pgeom->invJ  = &geom->suppInvJ[0][(c*Np+p)*dE*dE];
2396587ee25SMatthew G. Knepley       pgeom->detJ  = &geom->suppDetJ[0][c*Np+p];
2406587ee25SMatthew G. Knepley     } else {
2416587ee25SMatthew G. Knepley       pgeom->J    = &geom->J[(c*Np+p)*dE*dE];
2426587ee25SMatthew G. Knepley       pgeom->invJ = &geom->invJ[(c*Np+p)*dE*dE];
2436587ee25SMatthew G. Knepley       pgeom->detJ = &geom->detJ[c*Np+p];
2446587ee25SMatthew G. Knepley     }
2456587ee25SMatthew G. Knepley   }
2466587ee25SMatthew G. Knepley   PetscFunctionReturn(0);
2476587ee25SMatthew G. Knepley }
2486587ee25SMatthew G. Knepley 
2492b99622eSMatthew G. Knepley /*@
2502b99622eSMatthew G. Knepley   PetscFEGeomComplete - Calculate derived quntites from base geometry specification
2512b99622eSMatthew G. Knepley 
2522b99622eSMatthew G. Knepley   Input Parameter:
2532b99622eSMatthew G. Knepley . geom - PetscFEGeom object
2542b99622eSMatthew G. Knepley 
2552b99622eSMatthew G. Knepley   Level: intermediate
2562b99622eSMatthew G. Knepley 
2572b99622eSMatthew G. Knepley .seealso: PetscFEGeomCreate()
2582b99622eSMatthew G. Knepley @*/
25920cf1dd8SToby Isaac PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom)
26020cf1dd8SToby Isaac {
26120cf1dd8SToby Isaac   PetscInt i, j, N, dE;
26220cf1dd8SToby Isaac 
26320cf1dd8SToby Isaac   PetscFunctionBeginHot;
26420cf1dd8SToby Isaac   N = geom->numPoints * geom->numCells;
26520cf1dd8SToby Isaac   dE = geom->dimEmbed;
26620cf1dd8SToby Isaac   switch (dE) {
26720cf1dd8SToby Isaac   case 3:
26820cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
26920cf1dd8SToby Isaac       DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
27020cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
27120cf1dd8SToby Isaac     }
27220cf1dd8SToby Isaac     break;
27320cf1dd8SToby Isaac   case 2:
27420cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
27520cf1dd8SToby Isaac       DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
27620cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
27720cf1dd8SToby Isaac     }
27820cf1dd8SToby Isaac     break;
27920cf1dd8SToby Isaac   case 1:
28020cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
28120cf1dd8SToby Isaac       geom->detJ[i] = PetscAbsReal(geom->J[i]);
28220cf1dd8SToby Isaac       if (geom->invJ) geom->invJ[i] = 1. / geom->J[i];
28320cf1dd8SToby Isaac     }
28420cf1dd8SToby Isaac     break;
28520cf1dd8SToby Isaac   }
28620cf1dd8SToby Isaac   if (geom->n) {
28720cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
28820cf1dd8SToby Isaac       for (j = 0; j < dE; j++) {
28920cf1dd8SToby Isaac         geom->n[dE*i + j] = geom->J[dE*dE*i + dE*j + dE-1] * ((dE == 2) ? -1. : 1.);
29020cf1dd8SToby Isaac       }
29120cf1dd8SToby Isaac     }
29220cf1dd8SToby Isaac   }
29320cf1dd8SToby Isaac   PetscFunctionReturn(0);
29420cf1dd8SToby Isaac }
295