xref: /petsc/src/dm/impls/plex/plexproject.c (revision 559a1558ba4e4cd6f39f4336c756d8eb59ac12ff)
147923291SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
247923291SMatthew G. Knepley 
38c6c5593SMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
48c6c5593SMatthew G. Knepley 
50de51b56SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS prob, PetscReal time, PetscFECellGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
68c6c5593SMatthew G. Knepley                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
78c6c5593SMatthew G. Knepley                                                   PetscScalar values[])
847923291SMatthew G. Knepley {
99c3cf19fSMatthew G. Knepley   PetscInt       Nf, *Nc, f, totDim, spDim, d, v;
108c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
118c6c5593SMatthew G. Knepley 
128c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
138c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
14496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
159c3cf19fSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
168c6c5593SMatthew G. Knepley   /* Get values for closure */
178c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
188c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
198c6c5593SMatthew G. Knepley 
208c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
218c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
229c3cf19fSMatthew G. Knepley     for (d = 0; d < spDim; ++d, ++v) {
238c6c5593SMatthew G. Knepley       if (funcs[f]) {
248c6c5593SMatthew G. Knepley         if (isFE[f]) {ierr = PetscDualSpaceApply(sp[f], d, time, fegeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
258c6c5593SMatthew G. Knepley         else         {ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
268c6c5593SMatthew G. Knepley       } else {
279c3cf19fSMatthew G. Knepley         values[v] = 0.0;
288c6c5593SMatthew G. Knepley       }
298c6c5593SMatthew G. Knepley     }
309c3cf19fSMatthew G. Knepley   }
318c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
328c6c5593SMatthew G. Knepley }
338c6c5593SMatthew G. Knepley 
340de51b56SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFECellGeom *fegeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
350c898c18SMatthew G. Knepley                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
368c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
378c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
388c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
39191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
408c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
418c6c5593SMatthew G. Knepley {
428c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
439c3cf19fSMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL, *bc;
448c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
458c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
46191494d9SMatthew G. Knepley   const PetscScalar *constants;
478c6c5593SMatthew G. Knepley   PetscReal         *x;
48496733e2SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
499b168fc5SMatthew G. Knepley   const PetscInt     dim = fegeom->dim, dimEmbed = fegeom->dimEmbed;
509b168fc5SMatthew G. Knepley   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, c, v, tp = 0;
518c6c5593SMatthew G. Knepley   PetscErrorCode     ierr;
528c6c5593SMatthew G. Knepley 
538c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
548c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
55496733e2SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
56496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
578c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
588c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
599c3cf19fSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
60496733e2SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
61191494d9SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
628c6c5593SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
63496733e2SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
648c6c5593SMatthew G. Knepley   if (dmAux) {
651ba34526SMatthew G. Knepley     PetscInt subp;
661ba34526SMatthew G. Knepley 
67*559a1558SMatthew G. Knepley     ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr);
681ba34526SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
698c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
70496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
71496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
728c6c5593SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
738c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
748c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
7511d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
76496733e2SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
771ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
788c6c5593SMatthew G. Knepley   }
798c6c5593SMatthew G. Knepley   /* Get values for closure */
808c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
818c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
828c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
839c3cf19fSMatthew G. Knepley     for (d = 0; d < spDim; ++d, ++v) {
848c6c5593SMatthew G. Knepley       if (funcs[f]) {
858c6c5593SMatthew G. Knepley         PetscQuadrature  quad;
869c3cf19fSMatthew G. Knepley         const PetscReal *points, *weights;
878c6c5593SMatthew G. Knepley         PetscInt         numPoints, q;
888c6c5593SMatthew G. Knepley 
898c6c5593SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
909c3cf19fSMatthew G. Knepley         ierr = PetscQuadratureGetData(quad, NULL, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
910c898c18SMatthew G. Knepley         for (q = 0; q < numPoints; ++q, ++tp) {
929b168fc5SMatthew G. Knepley           CoordinatesRefToReal(dimEmbed, dim, fegeom->v0, fegeom->J, &points[q*dim], x);
93496733e2SMatthew G. Knepley           EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, fegeom->invJ, coefficients, coefficients_t, u, u_x, u_t);
941ba34526SMatthew G. Knepley           if (probAux) {EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, fegeom->invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);}
95191494d9SMatthew G. Knepley           (*funcs[f])(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, x, numConstants, constants, bc);
961c531cf8SMatthew G. Knepley           if (Ncc) {
971c531cf8SMatthew G. Knepley             for (c = 0; c < Ncc; ++c) values[v] += bc[comps[c]]*weights[comps[c]];
981c531cf8SMatthew G. Knepley           } else {
999c3cf19fSMatthew G. Knepley             for (c = 0; c < Nc[f]; ++c) values[v] += bc[c]*weights[c];
1008c6c5593SMatthew G. Knepley           }
1011c531cf8SMatthew G. Knepley         }
1028c6c5593SMatthew G. Knepley       } else {
1039c3cf19fSMatthew G. Knepley         values[v] = 0.0;
1048c6c5593SMatthew G. Knepley       }
1058c6c5593SMatthew G. Knepley     }
1068c6c5593SMatthew G. Knepley   }
1078c6c5593SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1088c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
1098c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1108c6c5593SMatthew G. Knepley }
1118c6c5593SMatthew G. Knepley 
1120de51b56SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
1131c531cf8SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
1141c531cf8SMatthew G. Knepley                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
1158c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
1168c6c5593SMatthew G. Knepley {
1178c6c5593SMatthew G. Knepley   PetscFECellGeom fegeom;
1188c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
1198c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
1208c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
1218c6c5593SMatthew G. Knepley 
1228c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1238c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1248c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1258c6c5593SMatthew G. Knepley   if (hasFE) {
12643048b12SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, fegeom.v0, fegeom.J, fegeom.invJ, &fegeom.detJ);CHKERRQ(ierr);
1270de51b56SMatthew G. Knepley     fegeom.dim      = dim - effectiveHeight;
1288c6c5593SMatthew G. Knepley     fegeom.dimEmbed = dimEmbed;
1298c6c5593SMatthew G. Knepley   }
1308c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
1318c6c5593SMatthew G. Knepley   switch (type) {
1328c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
1338c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
1340de51b56SMatthew G. Knepley     ierr = DMProjectPoint_Func_Private(dm, prob, time, &fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
1358c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
1368c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
1370de51b56SMatthew G. Knepley     ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, &fegeom, sp, p, Ncc, comps,
1380c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
1390c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
1400c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1410c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
142191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
1438c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
1448c6c5593SMatthew G. Knepley   }
1458c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1468c6c5593SMatthew G. Knepley }
1478c6c5593SMatthew G. Knepley 
1480de51b56SMatthew G. Knepley /*
1490de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
1500de51b56SMatthew G. Knepley 
1510de51b56SMatthew G. Knepley   There are several different scenarios:
1520de51b56SMatthew G. Knepley 
1530de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
1540de51b56SMatthew G. Knepley 
1550de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
1560de51b56SMatthew G. Knepley 
1570de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
1580de51b56SMatthew G. Knepley 
1590de51b56SMatthew G. Knepley      Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation.
1600de51b56SMatthew G. Knepley 
1610de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
1620de51b56SMatthew G. Knepley 
1630de51b56SMatthew G. Knepley      Here minHeight=1 and auxbd=PETSC_TRUE since we loop over faces and use data only supported on those faces. This is common when imposing Dirichlet boundary conditions.
1640de51b56SMatthew G. Knepley 
1650de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
1660de51b56SMatthew G. Knepley 
1670de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
1680de51b56SMatthew G. Knepley 
1690de51b56SMatthew G. Knepley   If we are using an input field (DM_BC_ESSENTIAL_FIELD or DM_BC_NATURAL_FIELD), we need to evaluate it at all the quadrature points of the dual basis functionals.
1700de51b56SMatthew G. Knepley     - We use effectiveHeight to mean the height above our incoming DS. For example, if the DS is for a submesh then the effective height is zero, whereas if the DS
1710de51b56SMatthew G. Knepley       is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When th effective height is nonzero, we need to extract
1720de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
1730de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
1740de51b56SMatthew G. Knepley 
1750de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
1760de51b56SMatthew G. Knepley 
1770de51b56SMatthew G. Knepley   If we have a label, we iterate over those points. This will probably break the maxHeight functionality since we do not check the height of those points.
1780de51b56SMatthew G. Knepley */
17946fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
1801c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
1818c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
1828c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
1838c6c5593SMatthew G. Knepley {
1848c6c5593SMatthew G. Knepley   DM              dmAux = NULL;
1850c898c18SMatthew G. Knepley   PetscDS         prob, probAux = NULL;
1868c6c5593SMatthew G. Knepley   Vec             localA = NULL;
18747923291SMatthew G. Knepley   PetscSection    section;
1888c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
1890c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
1908c6c5593SMatthew G. Knepley   PetscInt       *Nc;
1912af58022SMatthew G. Knepley   PetscInt        dim, dimEmbed, minHeight, maxHeight, h, Nf, NfAux = 0, f;
19288aca1feSMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE;
19347923291SMatthew G. Knepley   PetscErrorCode  ierr;
19447923291SMatthew G. Knepley 
19547923291SMatthew G. Knepley   PetscFunctionBegin;
19688aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
19788aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
1988c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1992af58022SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
2000de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
2010de51b56SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
20288aca1feSMatthew G. Knepley     if (!minHeight && dmAux) {
20388aca1feSMatthew G. Knepley       DMLabel spmap;
20488aca1feSMatthew G. Knepley 
20588aca1feSMatthew G. Knepley       /* If dmAux is a surface, then force the projection to take place over a surface */
20688aca1feSMatthew G. Knepley       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
20788aca1feSMatthew G. Knepley       if (spmap) {minHeight = 1; auxBd = PETSC_TRUE;}
20888aca1feSMatthew G. Knepley     }
2090de51b56SMatthew G. Knepley   }
2108c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
2112af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
2128c6c5593SMatthew G. Knepley   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %d not in [0, %d)\n", maxHeight, dim);}
2130c898c18SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
2148c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
21547923291SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
21647923291SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
2170c898c18SMatthew G. Knepley   if (dmAux) {
2180c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
2190c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
2200c898c18SMatthew G. Knepley   }
221496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
222496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
2238c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
2248c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
2258c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
2268c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
22747923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
22847923291SMatthew G. Knepley     PetscObject  obj;
22947923291SMatthew G. Knepley     PetscClassId id;
23047923291SMatthew G. Knepley 
23147923291SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
23247923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
23347923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
23447923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
23547923291SMatthew G. Knepley 
23647923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
23747923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
23847923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
23947923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
24047923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
24147923291SMatthew G. Knepley 
24247923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
24347923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
24447923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
24547923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
24647923291SMatthew G. Knepley   }
2470c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
24874fc349aSMatthew G. Knepley     PetscInt   effectiveHeight = auxBd ? minHeight : 0;
24974fc349aSMatthew G. Knepley     PetscFE    fem, subfem;
2500c898c18SMatthew G. Knepley     PetscReal *points;
25174fc349aSMatthew G. Knepley     PetscInt   numPoints, spDim, qdim = 0, d;
2520c898c18SMatthew G. Knepley 
2532af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
2540c898c18SMatthew G. Knepley     numPoints = 0;
2550c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
25674fc349aSMatthew G. Knepley       if (!effectiveHeight) {sp[f] = cellsp[f];}
25774fc349aSMatthew G. Knepley       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
25874fc349aSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
2590c898c18SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
2600c898c18SMatthew G. Knepley         if (funcs[f]) {
2610c898c18SMatthew G. Knepley           PetscQuadrature quad;
2620c898c18SMatthew G. Knepley           PetscInt        Nq;
2630c898c18SMatthew G. Knepley 
26474fc349aSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
26574fc349aSMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
2660c898c18SMatthew G. Knepley           numPoints += Nq;
2670c898c18SMatthew G. Knepley         }
2680c898c18SMatthew G. Knepley       }
2690c898c18SMatthew G. Knepley     }
27074fc349aSMatthew G. Knepley     ierr = PetscMalloc1(numPoints*qdim, &points);CHKERRQ(ierr);
2710c898c18SMatthew G. Knepley     numPoints = 0;
2720c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
27374fc349aSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
2740c898c18SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
2750c898c18SMatthew G. Knepley         if (funcs[f]) {
2760c898c18SMatthew G. Knepley           PetscQuadrature  quad;
2770c898c18SMatthew G. Knepley           const PetscReal *qpoints;
2780c898c18SMatthew G. Knepley           PetscInt         Nq, q;
2790c898c18SMatthew G. Knepley 
28074fc349aSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
2819c3cf19fSMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, NULL);CHKERRQ(ierr);
28274fc349aSMatthew G. Knepley           for (q = 0; q < Nq*qdim; ++q) points[numPoints*qdim+q] = qpoints[q];
2830c898c18SMatthew G. Knepley           numPoints += Nq;
2840c898c18SMatthew G. Knepley         }
2850c898c18SMatthew G. Knepley       }
2860c898c18SMatthew G. Knepley     }
2870c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
2880c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
2890c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
2900c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
29174fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
29274fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
29374fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
2940c898c18SMatthew G. Knepley     }
2950c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
2960c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
29774fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
29874fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
29974fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
3000c898c18SMatthew G. Knepley     }
3010c898c18SMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
3020c898c18SMatthew G. Knepley   }
30347923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
3042af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
30588aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
3060de51b56SMatthew G. Knepley     PetscDS      probEff         = prob;
3078c6c5593SMatthew G. Knepley     PetscScalar *values;
3088c6c5593SMatthew G. Knepley     PetscBool   *fieldActive;
30984f7fa7aSMatthew G. Knepley     PetscInt     pStart, pEnd, p, spDim, totDim, numValues;
3108c6c5593SMatthew G. Knepley 
3110de51b56SMatthew G. Knepley     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);}
31247923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
3138c6c5593SMatthew G. Knepley     if (!h) {
3148c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
3158c6c5593SMatthew G. Knepley 
3168c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3178c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
3188c6c5593SMatthew G. Knepley     }
31947923291SMatthew G. Knepley     if (pEnd <= pStart) continue;
32047923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
32147923291SMatthew G. Knepley     totDim = 0;
32247923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
323bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
32447923291SMatthew G. Knepley         sp[f] = cellsp[f];
32547923291SMatthew G. Knepley       } else {
326bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
32747923291SMatthew G. Knepley         if (!sp[f]) continue;
32847923291SMatthew G. Knepley       }
32947923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
3309c3cf19fSMatthew G. Knepley       totDim += spDim;
33147923291SMatthew G. Knepley     }
33247923291SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
3338c6c5593SMatthew G. Knepley     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim);
33447923291SMatthew G. Knepley     if (!totDim) continue;
33547923291SMatthew G. Knepley     /* Loop over points at this height */
33669291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
33769291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
3388c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
3398c6c5593SMatthew G. Knepley     if (label) {
3408c6c5593SMatthew G. Knepley       PetscInt i;
34147923291SMatthew G. Knepley 
34247923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
34347923291SMatthew G. Knepley         IS              pointIS;
34447923291SMatthew G. Knepley         const PetscInt *points;
3458c6c5593SMatthew G. Knepley         PetscInt        n;
34647923291SMatthew G. Knepley 
34747923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
34847923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
34947923291SMatthew G. Knepley         ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
35047923291SMatthew G. Knepley         ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
35147923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
35247923291SMatthew G. Knepley           const PetscInt  point = points[p];
35347923291SMatthew G. Knepley 
35447923291SMatthew G. Knepley           if ((point < pStart) || (point >= pEnd)) continue;
355b420b6ccSMatthew G. Knepley           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
3560de51b56SMatthew G. Knepley           ierr = DMProjectPoint_Private(dm, probEff, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
35747923291SMatthew G. Knepley           if (ierr) {
35847923291SMatthew G. Knepley             PetscErrorCode ierr2;
35969291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
36069291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
36147923291SMatthew G. Knepley             CHKERRQ(ierr);
36247923291SMatthew G. Knepley           }
363ba322698SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
36447923291SMatthew G. Knepley         }
36547923291SMatthew G. Knepley         ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
36647923291SMatthew G. Knepley         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
36747923291SMatthew G. Knepley       }
3688c6c5593SMatthew G. Knepley     } else {
3698c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
370b420b6ccSMatthew G. Knepley         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
3710de51b56SMatthew G. Knepley         ierr = DMProjectPoint_Private(dm, probEff, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
3728c6c5593SMatthew G. Knepley         if (ierr) {
3738c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
37469291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
37569291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
3768c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
3778c6c5593SMatthew G. Knepley         }
378ba322698SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
3798c6c5593SMatthew G. Knepley       }
3808c6c5593SMatthew G. Knepley     }
38169291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
38269291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
38347923291SMatthew G. Knepley   }
3848c6c5593SMatthew G. Knepley   /* Cleanup */
3850c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
38674fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
38774fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
3880c898c18SMatthew G. Knepley 
3890c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
3900c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
3910c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
39274fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
39374fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
39474fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
3950c898c18SMatthew G. Knepley     }
3960c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
3970c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
39874fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
39974fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
40074fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
4010c898c18SMatthew G. Knepley     }
4020c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
4030c898c18SMatthew G. Knepley   }
404496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
4058c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
4068c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
40747923291SMatthew G. Knepley }
4088c6c5593SMatthew G. Knepley 
4098c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
4108c6c5593SMatthew G. Knepley {
4118c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
4128c6c5593SMatthew G. Knepley 
4138c6c5593SMatthew G. Knepley   PetscFunctionBegin;
4141c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
4158c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
4168c6c5593SMatthew G. Knepley }
4178c6c5593SMatthew G. Knepley 
4181c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
4198c6c5593SMatthew G. Knepley {
4208c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
4218c6c5593SMatthew G. Knepley 
4228c6c5593SMatthew G. Knepley   PetscFunctionBegin;
4231c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
42447923291SMatthew G. Knepley   PetscFunctionReturn(0);
42547923291SMatthew G. Knepley }
42647923291SMatthew G. Knepley 
4278c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
42847923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
42947923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
43047923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
431191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
43247923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
43347923291SMatthew G. Knepley {
43447923291SMatthew G. Knepley   PetscErrorCode ierr;
43547923291SMatthew G. Knepley 
43647923291SMatthew G. Knepley   PetscFunctionBegin;
4371c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
4388c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
43947923291SMatthew G. Knepley }
44047923291SMatthew G. Knepley 
4411c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
4428c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
4438c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
4448c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
445191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
4468c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
4478c6c5593SMatthew G. Knepley {
4488c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
44947923291SMatthew G. Knepley 
4508c6c5593SMatthew G. Knepley   PetscFunctionBegin;
4511c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
45247923291SMatthew G. Knepley   PetscFunctionReturn(0);
45347923291SMatthew G. Knepley }
454