xref: /petsc/src/dm/dt/fe/interface/fegeom.c (revision 6587ee25fe32caa17c23170c94394ddbe328f69a)
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;
347489efa5SMatthew G. Knepley   N = numCells * Nq;
3520cf1dd8SToby Isaac   ierr = PetscCalloc3(N * dimEmbed, &g->v, N * dimEmbed * dimEmbed, &g->J, N, &g->detJ);CHKERRQ(ierr);
3620cf1dd8SToby Isaac   if (faceData) {
3741418987SMatthew G. Knepley     ierr = PetscCalloc2(numCells, &g->face, N * dimEmbed, &g->n);CHKERRQ(ierr);
3841418987SMatthew G. Knepley     ierr = PetscCalloc6(N * dimEmbed * dimEmbed, &(g->suppJ[0]),    N * dimEmbed * dimEmbed, &(g->suppJ[1]),
3941418987SMatthew G. Knepley                         N * dimEmbed * dimEmbed, &(g->suppInvJ[0]), N * dimEmbed * dimEmbed, &(g->suppInvJ[1]),
4041418987SMatthew G. Knepley                         N,                       &(g->suppDetJ[0]), N,                       &(g->suppDetJ[1]));CHKERRQ(ierr);
4120cf1dd8SToby Isaac   }
424e2f966aSMatthew G. Knepley   ierr = PetscCalloc1(N * dimEmbed * dimEmbed, &g->invJ);CHKERRQ(ierr);
4320cf1dd8SToby Isaac   *geom = g;
4420cf1dd8SToby Isaac   PetscFunctionReturn(0);
4520cf1dd8SToby Isaac }
4620cf1dd8SToby Isaac 
472b99622eSMatthew G. Knepley /*@C
482b99622eSMatthew G. Knepley   PetscFEGeomDestroy - Destroy a PetscFEGeom object
492b99622eSMatthew G. Knepley 
502b99622eSMatthew G. Knepley   Input Parameter:
512b99622eSMatthew G. Knepley . geom - PetscFEGeom object
522b99622eSMatthew G. Knepley 
532b99622eSMatthew G. Knepley   Level: beginner
542b99622eSMatthew G. Knepley 
552b99622eSMatthew G. Knepley .seealso: PetscFEGeomCreate()
562b99622eSMatthew G. Knepley @*/
5720cf1dd8SToby Isaac PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **geom)
5820cf1dd8SToby Isaac {
5920cf1dd8SToby Isaac   PetscErrorCode ierr;
6020cf1dd8SToby Isaac 
6120cf1dd8SToby Isaac   PetscFunctionBegin;
6220cf1dd8SToby Isaac   if (!*geom) PetscFunctionReturn(0);
6320cf1dd8SToby Isaac   ierr = PetscFree3((*geom)->v,(*geom)->J,(*geom)->detJ);CHKERRQ(ierr);
6420cf1dd8SToby Isaac   ierr = PetscFree((*geom)->invJ);CHKERRQ(ierr);
6541418987SMatthew G. Knepley   ierr = PetscFree2((*geom)->face,(*geom)->n);CHKERRQ(ierr);
6641418987SMatthew G. Knepley   ierr = PetscFree6((*geom)->suppJ[0],(*geom)->suppJ[1],(*geom)->suppInvJ[0],(*geom)->suppInvJ[1],(*geom)->suppDetJ[0],(*geom)->suppDetJ[1]);CHKERRQ(ierr);
6720cf1dd8SToby Isaac   ierr = PetscFree(*geom);CHKERRQ(ierr);
6820cf1dd8SToby Isaac   PetscFunctionReturn(0);
6920cf1dd8SToby Isaac }
7020cf1dd8SToby Isaac 
712b99622eSMatthew G. Knepley /*@C
722b99622eSMatthew G. Knepley   PetscFEGeomGetChunk - Get a chunk of cells in the group as a PetscFEGeom
732b99622eSMatthew G. Knepley 
742b99622eSMatthew G. Knepley   Input Parameters:
752b99622eSMatthew G. Knepley + geom   - PetscFEGeom object
762b99622eSMatthew G. Knepley . cStart - The first cell in the chunk
772b99622eSMatthew G. Knepley - cEnd   - The first cell not in the chunk
782b99622eSMatthew G. Knepley 
792b99622eSMatthew G. Knepley   Output Parameter:
802b99622eSMatthew G. Knepley . chunkGeom - The chunk of cells
812b99622eSMatthew G. Knepley 
822b99622eSMatthew G. Knepley   Level: intermediate
832b99622eSMatthew G. Knepley 
842b99622eSMatthew G. Knepley .seealso: PetscFEGeomRestoreChunk(), PetscFEGeomCreate()
852b99622eSMatthew G. Knepley @*/
8620cf1dd8SToby Isaac PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
8720cf1dd8SToby Isaac {
8820cf1dd8SToby Isaac   PetscInt       Nq;
8920cf1dd8SToby Isaac   PetscInt       dE;
9020cf1dd8SToby Isaac   PetscErrorCode ierr;
9120cf1dd8SToby Isaac 
9220cf1dd8SToby Isaac   PetscFunctionBegin;
9320cf1dd8SToby Isaac   PetscValidPointer(geom,1);
94064a246eSJacob Faibussowitsch   PetscValidPointer(chunkGeom,4);
9520cf1dd8SToby Isaac   if (!(*chunkGeom)) {
9620cf1dd8SToby Isaac     ierr = PetscNew(chunkGeom);CHKERRQ(ierr);
9720cf1dd8SToby Isaac   }
9820cf1dd8SToby Isaac   Nq = geom->numPoints;
9920cf1dd8SToby Isaac   dE= geom->dimEmbed;
10020cf1dd8SToby Isaac   (*chunkGeom)->dim = geom->dim;
10120cf1dd8SToby Isaac   (*chunkGeom)->dimEmbed = geom->dimEmbed;
10220cf1dd8SToby Isaac   (*chunkGeom)->numPoints = geom->numPoints;
10320cf1dd8SToby Isaac   (*chunkGeom)->numCells = cEnd - cStart;
10420cf1dd8SToby Isaac   (*chunkGeom)->xi = geom->xi;
10520cf1dd8SToby Isaac   (*chunkGeom)->v = &geom->v[Nq*dE*cStart];
10620cf1dd8SToby Isaac   (*chunkGeom)->J = &geom->J[Nq*dE*dE*cStart];
10720cf1dd8SToby Isaac   (*chunkGeom)->invJ = (geom->invJ) ? &geom->invJ[Nq*dE*dE*cStart] : NULL;
10820cf1dd8SToby Isaac   (*chunkGeom)->detJ = &geom->detJ[Nq*cStart];
10920cf1dd8SToby Isaac   (*chunkGeom)->n = geom->n ? &geom->n[Nq*dE*cStart] : NULL;
11020cf1dd8SToby Isaac   (*chunkGeom)->face = geom->face ? &geom->face[cStart] : NULL;
11141418987SMatthew G. Knepley   (*chunkGeom)->suppJ[0]    = geom->suppJ[0]    ? &geom->suppJ[0][Nq*dE*dE*cStart]    : NULL;
11241418987SMatthew G. Knepley   (*chunkGeom)->suppJ[1]    = geom->suppJ[1]    ? &geom->suppJ[1][Nq*dE*dE*cStart]    : NULL;
11320cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[0] = geom->suppInvJ[0] ? &geom->suppInvJ[0][Nq*dE*dE*cStart] : NULL;
11420cf1dd8SToby Isaac   (*chunkGeom)->suppInvJ[1] = geom->suppInvJ[1] ? &geom->suppInvJ[1][Nq*dE*dE*cStart] : NULL;
11541418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[0] = geom->suppDetJ[0] ? &geom->suppDetJ[0][Nq*cStart]       : NULL;
11641418987SMatthew G. Knepley   (*chunkGeom)->suppDetJ[1] = geom->suppDetJ[1] ? &geom->suppDetJ[1][Nq*cStart]       : NULL;
11720cf1dd8SToby Isaac   (*chunkGeom)->isAffine = geom->isAffine;
11820cf1dd8SToby Isaac   PetscFunctionReturn(0);
11920cf1dd8SToby Isaac }
12020cf1dd8SToby Isaac 
1212b99622eSMatthew G. Knepley /*@C
1222b99622eSMatthew G. Knepley   PetscFEGeomRestoreChunk - Restore the chunk
1232b99622eSMatthew G. Knepley 
1242b99622eSMatthew G. Knepley   Input Parameters:
1252b99622eSMatthew G. Knepley + geom      - PetscFEGeom object
1262b99622eSMatthew G. Knepley . cStart    - The first cell in the chunk
1272b99622eSMatthew G. Knepley . cEnd      - The first cell not in the chunk
1282b99622eSMatthew G. Knepley - chunkGeom - The chunk of cells
1292b99622eSMatthew G. Knepley 
1302b99622eSMatthew G. Knepley   Level: intermediate
1312b99622eSMatthew G. Knepley 
1322b99622eSMatthew G. Knepley .seealso: PetscFEGeomGetChunk(), PetscFEGeomCreate()
1332b99622eSMatthew G. Knepley @*/
13420cf1dd8SToby Isaac PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *geom, PetscInt cStart, PetscInt cEnd, PetscFEGeom **chunkGeom)
13520cf1dd8SToby Isaac {
13620cf1dd8SToby Isaac   PetscErrorCode ierr;
13720cf1dd8SToby Isaac 
13820cf1dd8SToby Isaac   PetscFunctionBegin;
13920cf1dd8SToby Isaac   ierr = PetscFree(*chunkGeom);CHKERRQ(ierr);
14020cf1dd8SToby Isaac   PetscFunctionReturn(0);
14120cf1dd8SToby Isaac }
14220cf1dd8SToby Isaac 
143*6587ee25SMatthew G. Knepley /*@C
144*6587ee25SMatthew G. Knepley   PetscFEGeomGetPoint - Get the geometry for cell c at point p as a PetscFEGeom
145*6587ee25SMatthew G. Knepley 
146*6587ee25SMatthew G. Knepley   Input Parameters:
147*6587ee25SMatthew G. Knepley + geom    - PetscFEGeom object
148*6587ee25SMatthew G. Knepley . c       - The cell
149*6587ee25SMatthew G. Knepley . p       - The point
150*6587ee25SMatthew G. Knepley - pcoords - The reference coordinates of point p, or NULL
151*6587ee25SMatthew G. Knepley 
152*6587ee25SMatthew G. Knepley   Output Parameter:
153*6587ee25SMatthew G. Knepley . pgeom - The geometry of cell c at point p
154*6587ee25SMatthew G. Knepley 
155*6587ee25SMatthew G. Knepley   Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom,
156*6587ee25SMatthew G. Knepley   nothing needs to be done with it afterwards.
157*6587ee25SMatthew G. Knepley 
158*6587ee25SMatthew G. Knepley   In the affine case, pgeom must have storage for the integration point coordinates in pgeom->v if pcoords is passed in.
159*6587ee25SMatthew G. Knepley 
160*6587ee25SMatthew G. Knepley   Level: intermediate
161*6587ee25SMatthew G. Knepley 
162*6587ee25SMatthew G. Knepley .seealso: PetscFEGeomRestoreChunk(), PetscFEGeomCreate()
163*6587ee25SMatthew G. Knepley @*/
164*6587ee25SMatthew G. Knepley PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, const PetscReal pcoords[], PetscFEGeom *pgeom)
165*6587ee25SMatthew G. Knepley {
166*6587ee25SMatthew G. Knepley   const PetscInt dim = geom->dim;
167*6587ee25SMatthew G. Knepley   const PetscInt dE  = geom->dimEmbed;
168*6587ee25SMatthew G. Knepley   const PetscInt Np  = geom->numPoints;
169*6587ee25SMatthew G. Knepley 
170*6587ee25SMatthew G. Knepley   PetscFunctionBeginHot;
171*6587ee25SMatthew G. Knepley   pgeom->dim      = dim;
172*6587ee25SMatthew G. Knepley   pgeom->dimEmbed = dE;
173*6587ee25SMatthew G. Knepley   //pgeom->isAffine = geom->isAffine;
174*6587ee25SMatthew G. Knepley   if (geom->isAffine) {
175*6587ee25SMatthew G. Knepley     if (!p) {
176*6587ee25SMatthew G. Knepley       pgeom->xi   = geom->xi;
177*6587ee25SMatthew G. Knepley       pgeom->J    = &geom->J[c*Np*dE*dE];
178*6587ee25SMatthew G. Knepley       pgeom->invJ = &geom->invJ[c*Np*dE*dE];
179*6587ee25SMatthew G. Knepley       pgeom->detJ = &geom->detJ[c*Np];
180*6587ee25SMatthew G. Knepley       pgeom->n    = geom->n ? &geom->n[c*Np*dE] : NULL;
181*6587ee25SMatthew G. Knepley     }
182*6587ee25SMatthew G. Knepley     if (pcoords) {CoordinatesRefToReal(dE, dim, pgeom->xi, &geom->v[c*Np*dE], pgeom->J, pcoords, pgeom->v);}
183*6587ee25SMatthew G. Knepley   } else {
184*6587ee25SMatthew G. Knepley     pgeom->v    = &geom->v[(c*Np+p)*dE];
185*6587ee25SMatthew G. Knepley     pgeom->J    = &geom->J[(c*Np+p)*dE*dE];
186*6587ee25SMatthew G. Knepley     pgeom->invJ = &geom->invJ[(c*Np+p)*dE*dE];
187*6587ee25SMatthew G. Knepley     pgeom->detJ = &geom->detJ[c*Np+p];
188*6587ee25SMatthew G. Knepley     pgeom->n    = geom->n ? &geom->n[(c*Np+p)*dE] : NULL;
189*6587ee25SMatthew G. Knepley   }
190*6587ee25SMatthew G. Knepley   PetscFunctionReturn(0);
191*6587ee25SMatthew G. Knepley }
192*6587ee25SMatthew G. Knepley 
193*6587ee25SMatthew G. Knepley /*@C
194*6587ee25SMatthew G. Knepley   PetscFEGeomGetCellPoint - Get the cell geometry for face f at point p as a PetscFEGeom
195*6587ee25SMatthew G. Knepley 
196*6587ee25SMatthew G. Knepley   Input Parameters:
197*6587ee25SMatthew G. Knepley + geom    - PetscFEGeom object
198*6587ee25SMatthew G. Knepley . f       - The face
199*6587ee25SMatthew G. Knepley - p       - The point
200*6587ee25SMatthew G. Knepley 
201*6587ee25SMatthew G. Knepley   Output Parameter:
202*6587ee25SMatthew G. Knepley . pgeom - The cell geometry of face f at point p
203*6587ee25SMatthew G. Knepley 
204*6587ee25SMatthew G. Knepley   Note: For affine geometries, this only copies to pgeom at point 0. Since we copy pointers into pgeom,
205*6587ee25SMatthew G. Knepley   nothing needs to be done with it afterwards.
206*6587ee25SMatthew G. Knepley 
207*6587ee25SMatthew G. Knepley   Level: intermediate
208*6587ee25SMatthew G. Knepley 
209*6587ee25SMatthew G. Knepley .seealso: PetscFEGeomRestoreChunk(), PetscFEGeomCreate()
210*6587ee25SMatthew G. Knepley @*/
211*6587ee25SMatthew G. Knepley PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *geom, PetscInt c, PetscInt p, PetscFEGeom *pgeom)
212*6587ee25SMatthew G. Knepley {
213*6587ee25SMatthew G. Knepley   const PetscBool bd  = geom->dimEmbed > geom->dim ? PETSC_TRUE : PETSC_FALSE;
214*6587ee25SMatthew G. Knepley   const PetscInt  dim = bd ? geom->dimEmbed : geom->dim;
215*6587ee25SMatthew G. Knepley   const PetscInt  dE  = geom->dimEmbed;
216*6587ee25SMatthew G. Knepley   const PetscInt  Np  = geom->numPoints;
217*6587ee25SMatthew G. Knepley 
218*6587ee25SMatthew G. Knepley   PetscFunctionBeginHot;
219*6587ee25SMatthew G. Knepley   pgeom->dim      = dim;
220*6587ee25SMatthew G. Knepley   pgeom->dimEmbed = dE;
221*6587ee25SMatthew G. Knepley   //pgeom->isAffine = geom->isAffine;
222*6587ee25SMatthew G. Knepley   if (geom->isAffine) {
223*6587ee25SMatthew G. Knepley     if (!p) {
224*6587ee25SMatthew G. Knepley       if (bd) {
225*6587ee25SMatthew G. Knepley         pgeom->J     = &geom->suppJ[0][c*Np*dE*dE];
226*6587ee25SMatthew G. Knepley         pgeom->invJ  = &geom->suppInvJ[0][c*Np*dE*dE];
227*6587ee25SMatthew G. Knepley         pgeom->detJ  = &geom->suppDetJ[0][c*Np];
228*6587ee25SMatthew G. Knepley       } else {
229*6587ee25SMatthew G. Knepley         pgeom->J    = &geom->J[c*Np*dE*dE];
230*6587ee25SMatthew G. Knepley         pgeom->invJ = &geom->invJ[c*Np*dE*dE];
231*6587ee25SMatthew G. Knepley         pgeom->detJ = &geom->detJ[c*Np];
232*6587ee25SMatthew G. Knepley       }
233*6587ee25SMatthew G. Knepley     }
234*6587ee25SMatthew G. Knepley   } else {
235*6587ee25SMatthew G. Knepley     if (bd) {
236*6587ee25SMatthew G. Knepley       pgeom->J     = &geom->suppJ[0][(c*Np+p)*dE*dE];
237*6587ee25SMatthew G. Knepley       pgeom->invJ  = &geom->suppInvJ[0][(c*Np+p)*dE*dE];
238*6587ee25SMatthew G. Knepley       pgeom->detJ  = &geom->suppDetJ[0][c*Np+p];
239*6587ee25SMatthew G. Knepley     } else {
240*6587ee25SMatthew G. Knepley       pgeom->J    = &geom->J[(c*Np+p)*dE*dE];
241*6587ee25SMatthew G. Knepley       pgeom->invJ = &geom->invJ[(c*Np+p)*dE*dE];
242*6587ee25SMatthew G. Knepley       pgeom->detJ = &geom->detJ[c*Np+p];
243*6587ee25SMatthew G. Knepley     }
244*6587ee25SMatthew G. Knepley   }
245*6587ee25SMatthew G. Knepley   PetscFunctionReturn(0);
246*6587ee25SMatthew G. Knepley }
247*6587ee25SMatthew G. Knepley 
2482b99622eSMatthew G. Knepley /*@
2492b99622eSMatthew G. Knepley   PetscFEGeomComplete - Calculate derived quntites from base geometry specification
2502b99622eSMatthew G. Knepley 
2512b99622eSMatthew G. Knepley   Input Parameter:
2522b99622eSMatthew G. Knepley . geom - PetscFEGeom object
2532b99622eSMatthew G. Knepley 
2542b99622eSMatthew G. Knepley   Level: intermediate
2552b99622eSMatthew G. Knepley 
2562b99622eSMatthew G. Knepley .seealso: PetscFEGeomCreate()
2572b99622eSMatthew G. Knepley @*/
25820cf1dd8SToby Isaac PetscErrorCode PetscFEGeomComplete(PetscFEGeom *geom)
25920cf1dd8SToby Isaac {
26020cf1dd8SToby Isaac   PetscInt i, j, N, dE;
26120cf1dd8SToby Isaac 
26220cf1dd8SToby Isaac   PetscFunctionBeginHot;
26320cf1dd8SToby Isaac   N = geom->numPoints * geom->numCells;
26420cf1dd8SToby Isaac   dE = geom->dimEmbed;
26520cf1dd8SToby Isaac   switch (dE) {
26620cf1dd8SToby Isaac   case 3:
26720cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
26820cf1dd8SToby Isaac       DMPlex_Det3D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
26920cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert3D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
27020cf1dd8SToby Isaac     }
27120cf1dd8SToby Isaac     break;
27220cf1dd8SToby Isaac   case 2:
27320cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
27420cf1dd8SToby Isaac       DMPlex_Det2D_Internal(&geom->detJ[i], &geom->J[dE*dE*i]);
27520cf1dd8SToby Isaac       if (geom->invJ) DMPlex_Invert2D_Internal(&geom->invJ[dE*dE*i], &geom->J[dE*dE*i], geom->detJ[i]);
27620cf1dd8SToby Isaac     }
27720cf1dd8SToby Isaac     break;
27820cf1dd8SToby Isaac   case 1:
27920cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
28020cf1dd8SToby Isaac       geom->detJ[i] = PetscAbsReal(geom->J[i]);
28120cf1dd8SToby Isaac       if (geom->invJ) geom->invJ[i] = 1. / geom->J[i];
28220cf1dd8SToby Isaac     }
28320cf1dd8SToby Isaac     break;
28420cf1dd8SToby Isaac   }
28520cf1dd8SToby Isaac   if (geom->n) {
28620cf1dd8SToby Isaac     for (i = 0; i < N; i++) {
28720cf1dd8SToby Isaac       for (j = 0; j < dE; j++) {
28820cf1dd8SToby Isaac         geom->n[dE*i + j] = geom->J[dE*dE*i + dE*j + dE-1] * ((dE == 2) ? -1. : 1.);
28920cf1dd8SToby Isaac       }
29020cf1dd8SToby Isaac     }
29120cf1dd8SToby Isaac   }
29220cf1dd8SToby Isaac   PetscFunctionReturn(0);
29320cf1dd8SToby Isaac }
294