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 5*0de51b56SMatthew 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 34*0de51b56SMatthew 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, §ion);CHKERRQ(ierr); 63496733e2SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 648c6c5593SMatthew G. Knepley if (dmAux) { 651ba34526SMatthew G. Knepley DMLabel spmap; 661ba34526SMatthew G. Knepley PetscInt subp; 671ba34526SMatthew G. Knepley 681ba34526SMatthew G. Knepley ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr); 691ba34526SMatthew G. Knepley if (spmap) { 701ba34526SMatthew G. Knepley IS subpointIS; 711ba34526SMatthew G. Knepley const PetscInt *subpoints; 721ba34526SMatthew G. Knepley PetscInt numSubpoints; 731ba34526SMatthew G. Knepley 741ba34526SMatthew G. Knepley ierr = DMPlexCreateSubpointIS(dmAux, &subpointIS);CHKERRQ(ierr); 751ba34526SMatthew G. Knepley ierr = ISGetLocalSize(subpointIS, &numSubpoints);CHKERRQ(ierr); 761ba34526SMatthew G. Knepley ierr = ISGetIndices(subpointIS, &subpoints);CHKERRQ(ierr); 771ba34526SMatthew G. Knepley ierr = PetscFindInt(p, numSubpoints, subpoints, &subp);CHKERRQ(ierr); 781ba34526SMatthew G. Knepley if (subp < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d not found in submesh", p); 791ba34526SMatthew G. Knepley ierr = ISRestoreIndices(subpointIS, &subpoints);CHKERRQ(ierr); 801ba34526SMatthew G. Knepley ierr = ISDestroy(&subpointIS);CHKERRQ(ierr); 811ba34526SMatthew G. Knepley } else subp = p; 821ba34526SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 838c6c5593SMatthew G. Knepley ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 84496733e2SMatthew G. Knepley ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 85496733e2SMatthew G. Knepley ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 868c6c5593SMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 878c6c5593SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 888c6c5593SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 8911d189daSMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 90496733e2SMatthew G. Knepley ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 911ba34526SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 928c6c5593SMatthew G. Knepley } 938c6c5593SMatthew G. Knepley /* Get values for closure */ 948c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 958c6c5593SMatthew G. Knepley if (!sp[f]) continue; 968c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 979c3cf19fSMatthew G. Knepley for (d = 0; d < spDim; ++d, ++v) { 988c6c5593SMatthew G. Knepley if (funcs[f]) { 998c6c5593SMatthew G. Knepley PetscQuadrature quad; 1009c3cf19fSMatthew G. Knepley const PetscReal *points, *weights; 1018c6c5593SMatthew G. Knepley PetscInt numPoints, q; 1028c6c5593SMatthew G. Knepley 1038c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr); 1049c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, NULL, &numPoints, &points, &weights);CHKERRQ(ierr); 1050c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 1069b168fc5SMatthew G. Knepley CoordinatesRefToReal(dimEmbed, dim, fegeom->v0, fegeom->J, &points[q*dim], x); 107496733e2SMatthew G. Knepley EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, fegeom->invJ, coefficients, coefficients_t, u, u_x, u_t); 1081ba34526SMatthew G. Knepley if (probAux) {EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, fegeom->invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);} 109191494d9SMatthew 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); 1101c531cf8SMatthew G. Knepley if (Ncc) { 1111c531cf8SMatthew G. Knepley for (c = 0; c < Ncc; ++c) values[v] += bc[comps[c]]*weights[comps[c]]; 1121c531cf8SMatthew G. Knepley } else { 1139c3cf19fSMatthew G. Knepley for (c = 0; c < Nc[f]; ++c) values[v] += bc[c]*weights[c]; 1148c6c5593SMatthew G. Knepley } 1151c531cf8SMatthew G. Knepley } 1168c6c5593SMatthew G. Knepley } else { 1179c3cf19fSMatthew G. Knepley values[v] = 0.0; 1188c6c5593SMatthew G. Knepley } 1198c6c5593SMatthew G. Knepley } 1208c6c5593SMatthew G. Knepley } 1218c6c5593SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 1228c6c5593SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 1238c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 1248c6c5593SMatthew G. Knepley } 1258c6c5593SMatthew G. Knepley 126*0de51b56SMatthew 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[], 1271c531cf8SMatthew G. Knepley PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 1281c531cf8SMatthew G. Knepley PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 1298c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 1308c6c5593SMatthew G. Knepley { 1318c6c5593SMatthew G. Knepley PetscFECellGeom fegeom; 1328c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom; 1338c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed; 1348c6c5593SMatthew G. Knepley PetscErrorCode ierr; 1358c6c5593SMatthew G. Knepley 1368c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 1378c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1388c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 1398c6c5593SMatthew G. Knepley if (hasFE) { 14043048b12SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, fegeom.v0, fegeom.J, fegeom.invJ, &fegeom.detJ);CHKERRQ(ierr); 141*0de51b56SMatthew G. Knepley fegeom.dim = dim - effectiveHeight; 1428c6c5593SMatthew G. Knepley fegeom.dimEmbed = dimEmbed; 1438c6c5593SMatthew G. Knepley } 1448c6c5593SMatthew G. Knepley if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 1458c6c5593SMatthew G. Knepley switch (type) { 1468c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL: 1478c6c5593SMatthew G. Knepley case DM_BC_NATURAL: 148*0de51b56SMatthew 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; 1498c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 1508c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD: 151*0de51b56SMatthew G. Knepley ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, &fegeom, sp, p, Ncc, comps, 1520c898c18SMatthew G. Knepley basisTab, basisDerTab, basisTabAux, basisDerTabAux, 1530c898c18SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 1540c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1550c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 156191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 1578c6c5593SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 1588c6c5593SMatthew G. Knepley } 1598c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 1608c6c5593SMatthew G. Knepley } 1618c6c5593SMatthew G. Knepley 162*0de51b56SMatthew G. Knepley /* 163*0de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 164*0de51b56SMatthew G. Knepley 165*0de51b56SMatthew G. Knepley There are several different scenarios: 166*0de51b56SMatthew G. Knepley 167*0de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data 168*0de51b56SMatthew G. Knepley 169*0de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells. 170*0de51b56SMatthew G. Knepley 171*0de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data 172*0de51b56SMatthew G. Knepley 173*0de51b56SMatthew 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. 174*0de51b56SMatthew G. Knepley 175*0de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data 176*0de51b56SMatthew G. Knepley 177*0de51b56SMatthew 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. 178*0de51b56SMatthew G. Knepley 179*0de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest. 180*0de51b56SMatthew G. Knepley 181*0de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 182*0de51b56SMatthew G. Knepley 183*0de51b56SMatthew 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. 184*0de51b56SMatthew 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 185*0de51b56SMatthew 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 186*0de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces. 187*0de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 188*0de51b56SMatthew G. Knepley 189*0de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 190*0de51b56SMatthew G. Knepley 191*0de51b56SMatthew 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. 192*0de51b56SMatthew G. Knepley */ 19346fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 1941c531cf8SMatthew G. Knepley PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], 1958c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 1968c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 1978c6c5593SMatthew G. Knepley { 1988c6c5593SMatthew G. Knepley DM dmAux = NULL; 1990c898c18SMatthew G. Knepley PetscDS prob, probAux = NULL; 2008c6c5593SMatthew G. Knepley Vec localA = NULL; 20147923291SMatthew G. Knepley PetscSection section; 2028c6c5593SMatthew G. Knepley PetscDualSpace *sp, *cellsp; 2030c898c18SMatthew G. Knepley PetscReal **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL; 2048c6c5593SMatthew G. Knepley PetscInt *Nc; 2052af58022SMatthew G. Knepley PetscInt dim, dimEmbed, minHeight, maxHeight, h, Nf, NfAux = 0, f; 20688aca1feSMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE; 20747923291SMatthew G. Knepley PetscErrorCode ierr; 20847923291SMatthew G. Knepley 20947923291SMatthew G. Knepley PetscFunctionBegin; 21088aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 21188aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr); 2128c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2132af58022SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr); 214*0de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */ 215*0de51b56SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 21688aca1feSMatthew G. Knepley if (!minHeight && dmAux) { 21788aca1feSMatthew G. Knepley DMLabel spmap; 21888aca1feSMatthew G. Knepley 21988aca1feSMatthew G. Knepley /* If dmAux is a surface, then force the projection to take place over a surface */ 22088aca1feSMatthew G. Knepley ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr); 22188aca1feSMatthew G. Knepley if (spmap) {minHeight = 1; auxBd = PETSC_TRUE;} 22288aca1feSMatthew G. Knepley } 223*0de51b56SMatthew G. Knepley } 2248c6c5593SMatthew G. Knepley ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr); 2252af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight); 2268c6c5593SMatthew 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);} 2270c898c18SMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 2288c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 22947923291SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 23047923291SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 2310c898c18SMatthew G. Knepley if (dmAux) { 2320c898c18SMatthew G. Knepley ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 2330c898c18SMatthew G. Knepley ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 2340c898c18SMatthew G. Knepley } 235496733e2SMatthew G. Knepley ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 236496733e2SMatthew G. Knepley ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr); 2378c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);} 2388c6c5593SMatthew G. Knepley else {cellsp = sp;} 2398c6c5593SMatthew G. Knepley if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 2408c6c5593SMatthew G. Knepley /* Get cell dual spaces */ 24147923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 24247923291SMatthew G. Knepley PetscObject obj; 24347923291SMatthew G. Knepley PetscClassId id; 24447923291SMatthew G. Knepley 24547923291SMatthew G. Knepley ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 24647923291SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 24747923291SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 24847923291SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 24947923291SMatthew G. Knepley 25047923291SMatthew G. Knepley hasFE = PETSC_TRUE; 25147923291SMatthew G. Knepley isFE[f] = PETSC_TRUE; 25247923291SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 25347923291SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 25447923291SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 25547923291SMatthew G. Knepley 25647923291SMatthew G. Knepley hasFV = PETSC_TRUE; 25747923291SMatthew G. Knepley isFE[f] = PETSC_FALSE; 25847923291SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 25947923291SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 26047923291SMatthew G. Knepley } 2610c898c18SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 26274fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 26374fc349aSMatthew G. Knepley PetscFE fem, subfem; 2640c898c18SMatthew G. Knepley PetscReal *points; 26574fc349aSMatthew G. Knepley PetscInt numPoints, spDim, qdim = 0, d; 2660c898c18SMatthew G. Knepley 2672af58022SMatthew G. Knepley if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 2680c898c18SMatthew G. Knepley numPoints = 0; 2690c898c18SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 27074fc349aSMatthew G. Knepley if (!effectiveHeight) {sp[f] = cellsp[f];} 27174fc349aSMatthew G. Knepley else {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);} 27274fc349aSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 2730c898c18SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 2740c898c18SMatthew G. Knepley if (funcs[f]) { 2750c898c18SMatthew G. Knepley PetscQuadrature quad; 2760c898c18SMatthew G. Knepley PetscInt Nq; 2770c898c18SMatthew G. Knepley 27874fc349aSMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr); 27974fc349aSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 2800c898c18SMatthew G. Knepley numPoints += Nq; 2810c898c18SMatthew G. Knepley } 2820c898c18SMatthew G. Knepley } 2830c898c18SMatthew G. Knepley } 28474fc349aSMatthew G. Knepley ierr = PetscMalloc1(numPoints*qdim, &points);CHKERRQ(ierr); 2850c898c18SMatthew G. Knepley numPoints = 0; 2860c898c18SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 28774fc349aSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 2880c898c18SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 2890c898c18SMatthew G. Knepley if (funcs[f]) { 2900c898c18SMatthew G. Knepley PetscQuadrature quad; 2910c898c18SMatthew G. Knepley const PetscReal *qpoints; 2920c898c18SMatthew G. Knepley PetscInt Nq, q; 2930c898c18SMatthew G. Knepley 29474fc349aSMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr); 2959c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, NULL);CHKERRQ(ierr); 29674fc349aSMatthew G. Knepley for (q = 0; q < Nq*qdim; ++q) points[numPoints*qdim+q] = qpoints[q]; 2970c898c18SMatthew G. Knepley numPoints += Nq; 2980c898c18SMatthew G. Knepley } 2990c898c18SMatthew G. Knepley } 3000c898c18SMatthew G. Knepley } 3010c898c18SMatthew G. Knepley ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr); 3020c898c18SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 3030c898c18SMatthew G. Knepley if (!isFE[f]) continue; 3040c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 30574fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 30674fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 30774fc349aSMatthew G. Knepley ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 3080c898c18SMatthew G. Knepley } 3090c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 3100c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 31174fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 31274fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 31374fc349aSMatthew G. Knepley ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 3140c898c18SMatthew G. Knepley } 3150c898c18SMatthew G. Knepley ierr = PetscFree(points);CHKERRQ(ierr); 3160c898c18SMatthew G. Knepley } 31747923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 3182af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) { 31988aca1feSMatthew G. Knepley PetscInt effectiveHeight = h - (auxBd ? 0 : minHeight); 320*0de51b56SMatthew G. Knepley PetscDS probEff = prob; 3218c6c5593SMatthew G. Knepley PetscScalar *values; 3228c6c5593SMatthew G. Knepley PetscBool *fieldActive; 32384f7fa7aSMatthew G. Knepley PetscInt pStart, pEnd, p, spDim, totDim, numValues; 3248c6c5593SMatthew G. Knepley 325*0de51b56SMatthew G. Knepley if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);} 32647923291SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 3278c6c5593SMatthew G. Knepley if (!h) { 3288c6c5593SMatthew G. Knepley PetscInt cEndInterior; 3298c6c5593SMatthew G. Knepley 3308c6c5593SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 3318c6c5593SMatthew G. Knepley pEnd = cEndInterior < 0 ? pEnd : cEndInterior; 3328c6c5593SMatthew G. Knepley } 33347923291SMatthew G. Knepley if (pEnd <= pStart) continue; 33447923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */ 33547923291SMatthew G. Knepley totDim = 0; 33647923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 337bec1f4aaSMatthew G. Knepley if (!effectiveHeight) { 33847923291SMatthew G. Knepley sp[f] = cellsp[f]; 33947923291SMatthew G. Knepley } else { 340bec1f4aaSMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr); 34147923291SMatthew G. Knepley if (!sp[f]) continue; 34247923291SMatthew G. Knepley } 34347923291SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 3449c3cf19fSMatthew G. Knepley totDim += spDim; 34547923291SMatthew G. Knepley } 34647923291SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr); 3478c6c5593SMatthew 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); 34847923291SMatthew G. Knepley if (!totDim) continue; 34947923291SMatthew G. Knepley /* Loop over points at this height */ 35069291d52SBarry Smith ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 35169291d52SBarry Smith ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 3528c6c5593SMatthew G. Knepley for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 3538c6c5593SMatthew G. Knepley if (label) { 3548c6c5593SMatthew G. Knepley PetscInt i; 35547923291SMatthew G. Knepley 35647923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 35747923291SMatthew G. Knepley IS pointIS; 35847923291SMatthew G. Knepley const PetscInt *points; 3598c6c5593SMatthew G. Knepley PetscInt n; 36047923291SMatthew G. Knepley 36147923291SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 36247923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 36347923291SMatthew G. Knepley ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 36447923291SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 36547923291SMatthew G. Knepley for (p = 0; p < n; ++p) { 36647923291SMatthew G. Knepley const PetscInt point = points[p]; 36747923291SMatthew G. Knepley 36847923291SMatthew G. Knepley if ((point < pStart) || (point >= pEnd)) continue; 369b420b6ccSMatthew G. Knepley ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 370*0de51b56SMatthew 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); 37147923291SMatthew G. Knepley if (ierr) { 37247923291SMatthew G. Knepley PetscErrorCode ierr2; 37369291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 37469291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 37547923291SMatthew G. Knepley CHKERRQ(ierr); 37647923291SMatthew G. Knepley } 37747923291SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr); 37847923291SMatthew G. Knepley } 37947923291SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 38047923291SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 38147923291SMatthew G. Knepley } 3828c6c5593SMatthew G. Knepley } else { 3838c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 384b420b6ccSMatthew G. Knepley ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 385*0de51b56SMatthew 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); 3868c6c5593SMatthew G. Knepley if (ierr) { 3878c6c5593SMatthew G. Knepley PetscErrorCode ierr2; 38869291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 38969291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 3908c6c5593SMatthew G. Knepley CHKERRQ(ierr); 3918c6c5593SMatthew G. Knepley } 3928c6c5593SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, values, mode);CHKERRQ(ierr); 3938c6c5593SMatthew G. Knepley } 3948c6c5593SMatthew G. Knepley } 39569291d52SBarry Smith ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 39669291d52SBarry Smith ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 39747923291SMatthew G. Knepley } 3988c6c5593SMatthew G. Knepley /* Cleanup */ 3990c898c18SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 40074fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 40174fc349aSMatthew G. Knepley PetscFE fem, subfem; 4020c898c18SMatthew G. Knepley 4030c898c18SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 4040c898c18SMatthew G. Knepley if (!isFE[f]) continue; 4050c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 40674fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 40774fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 40874fc349aSMatthew G. Knepley ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 4090c898c18SMatthew G. Knepley } 4100c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 4110c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 41274fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 41374fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 41474fc349aSMatthew G. Knepley ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 4150c898c18SMatthew G. Knepley } 4160c898c18SMatthew G. Knepley ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr); 4170c898c18SMatthew G. Knepley } 418496733e2SMatthew G. Knepley ierr = PetscFree2(isFE, sp);CHKERRQ(ierr); 4198c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);} 4208c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 42147923291SMatthew G. Knepley } 4228c6c5593SMatthew G. Knepley 4238c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 4248c6c5593SMatthew G. Knepley { 4258c6c5593SMatthew G. Knepley PetscErrorCode ierr; 4268c6c5593SMatthew G. Knepley 4278c6c5593SMatthew G. Knepley PetscFunctionBegin; 4281c531cf8SMatthew 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); 4298c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 4308c6c5593SMatthew G. Knepley } 4318c6c5593SMatthew G. Knepley 4321c531cf8SMatthew 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) 4338c6c5593SMatthew G. Knepley { 4348c6c5593SMatthew G. Knepley PetscErrorCode ierr; 4358c6c5593SMatthew G. Knepley 4368c6c5593SMatthew G. Knepley PetscFunctionBegin; 4371c531cf8SMatthew 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); 43847923291SMatthew G. Knepley PetscFunctionReturn(0); 43947923291SMatthew G. Knepley } 44047923291SMatthew G. Knepley 4418c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 44247923291SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 44347923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 44447923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 445191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 44647923291SMatthew G. Knepley InsertMode mode, Vec localX) 44747923291SMatthew G. Knepley { 44847923291SMatthew G. Knepley PetscErrorCode ierr; 44947923291SMatthew G. Knepley 45047923291SMatthew G. Knepley PetscFunctionBegin; 4511c531cf8SMatthew 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); 4528c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 45347923291SMatthew G. Knepley } 45447923291SMatthew G. Knepley 4551c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 4568c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 4578c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4588c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 459191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 4608c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 4618c6c5593SMatthew G. Knepley { 4628c6c5593SMatthew G. Knepley PetscErrorCode ierr; 46347923291SMatthew G. Knepley 4648c6c5593SMatthew G. Knepley PetscFunctionBegin; 4651c531cf8SMatthew 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); 46647923291SMatthew G. Knepley PetscFunctionReturn(0); 46747923291SMatthew G. Knepley } 468