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 5c330f8ffSToby Isaac static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS prob, PetscReal time, PetscFEGeom *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 { 9c330f8ffSToby Isaac PetscInt coordDim, Nf, *Nc, f, totDim, spDim, d, v, tp; 10ca3d3a14SMatthew G. Knepley PetscBool isAffine, transform; 118c6c5593SMatthew G. Knepley PetscErrorCode ierr; 128c6c5593SMatthew G. Knepley 138c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 14c330f8ffSToby Isaac ierr = DMGetCoordinateDim(dm,&coordDim);CHKERRQ(ierr); 15ca3d3a14SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 168c6c5593SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 17496733e2SMatthew G. Knepley ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 189c3cf19fSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 198c6c5593SMatthew G. Knepley /* Get values for closure */ 20c330f8ffSToby Isaac isAffine = fegeom->isAffine; 21c330f8ffSToby Isaac for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 228c6c5593SMatthew G. Knepley void * const ctx = ctxs ? ctxs[f] : NULL; 238c6c5593SMatthew G. Knepley 248c6c5593SMatthew G. Knepley if (!sp[f]) continue; 258c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 268c6c5593SMatthew G. Knepley if (funcs[f]) { 27c330f8ffSToby Isaac if (isFE[f]) { 28c330f8ffSToby Isaac PetscQuadrature allPoints; 29c330f8ffSToby Isaac PetscInt q, dim, numPoints; 30c330f8ffSToby Isaac const PetscReal *points; 31c330f8ffSToby Isaac PetscScalar *pointEval; 32c330f8ffSToby Isaac PetscReal *x; 33ca3d3a14SMatthew G. Knepley DM rdm; 34c330f8ffSToby Isaac 35ca3d3a14SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr); 36c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 37c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 38ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 39ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 40c330f8ffSToby Isaac for (q = 0; q < numPoints; q++, tp++) { 41c330f8ffSToby Isaac const PetscReal *v0; 42c330f8ffSToby Isaac 43c330f8ffSToby Isaac if (isAffine) { 44c330f8ffSToby Isaac CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x); 45c330f8ffSToby Isaac v0 = x; 468c6c5593SMatthew G. Knepley } else { 47c330f8ffSToby Isaac v0 = &fegeom->v[tp*coordDim]; 488c6c5593SMatthew G. Knepley } 49ec277c0fSMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;} 50c330f8ffSToby Isaac ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr); 51c330f8ffSToby Isaac } 524bee2e38SMatthew G. Knepley /* Transform point evaluations pointEval[q,c] */ 534bee2e38SMatthew G. Knepley ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr); 54c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 55ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 56ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 57c330f8ffSToby Isaac v += spDim; 58c330f8ffSToby Isaac } else { 59c330f8ffSToby Isaac for (d = 0; d < spDim; ++d, ++v) { 60c330f8ffSToby Isaac ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 61c330f8ffSToby Isaac } 62c330f8ffSToby Isaac } 63c330f8ffSToby Isaac } else { 64c330f8ffSToby Isaac for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;} 658c6c5593SMatthew G. Knepley } 669c3cf19fSMatthew G. Knepley } 678c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 688c6c5593SMatthew G. Knepley } 698c6c5593SMatthew G. Knepley 704bee2e38SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 710c898c18SMatthew G. Knepley PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 728c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 738c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 748c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 75191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 768c6c5593SMatthew G. Knepley PetscScalar values[]) 778c6c5593SMatthew G. Knepley { 788c6c5593SMatthew G. Knepley PetscSection section, sectionAux = NULL; 794bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 808c6c5593SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 818c6c5593SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 82191494d9SMatthew G. Knepley const PetscScalar *constants; 838c6c5593SMatthew G. Knepley PetscReal *x; 84496733e2SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 854bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 864bee2e38SMatthew G. Knepley const PetscInt dE = cgeom->dimEmbed; 87c330f8ffSToby Isaac PetscInt dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 88ca3d3a14SMatthew G. Knepley PetscBool isAffine, transform; 898c6c5593SMatthew G. Knepley PetscErrorCode ierr; 908c6c5593SMatthew G. Knepley 918c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 928c6c5593SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 93496733e2SMatthew G. Knepley ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 94496733e2SMatthew G. Knepley ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 958c6c5593SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 968c6c5593SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 979c3cf19fSMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 984bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 99191494d9SMatthew G. Knepley ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 100ca3d3a14SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 101e87a4003SBarry Smith ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 102496733e2SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 1038c6c5593SMatthew G. Knepley if (dmAux) { 10444171101SMatthew G. Knepley PetscInt subp; 1051ba34526SMatthew G. Knepley 10644171101SMatthew G. Knepley ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, p, &subp);CHKERRQ(ierr); 1071ba34526SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 1088c6c5593SMatthew G. Knepley ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 109496733e2SMatthew G. Knepley ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 110496733e2SMatthew G. Knepley ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 111e87a4003SBarry Smith ierr = DMGetSection(dmAux, §ionAux);CHKERRQ(ierr); 1128c6c5593SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 1138c6c5593SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 11411d189daSMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 1151ba34526SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 1168c6c5593SMatthew G. Knepley } 1178c6c5593SMatthew G. Knepley /* Get values for closure */ 1184bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 1194bee2e38SMatthew G. Knepley if (isAffine) { 1204bee2e38SMatthew G. Knepley fegeom.v = x; 1214bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 1224bee2e38SMatthew G. Knepley fegeom.J = cgeom->J; 1234bee2e38SMatthew G. Knepley fegeom.invJ = cgeom->invJ; 1244bee2e38SMatthew G. Knepley fegeom.detJ = cgeom->detJ; 1254bee2e38SMatthew G. Knepley } 1268c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 127c330f8ffSToby Isaac PetscQuadrature allPoints; 128c330f8ffSToby Isaac PetscInt q, dim, numPoints; 129c330f8ffSToby Isaac const PetscReal *points; 130c330f8ffSToby Isaac PetscScalar *pointEval; 131c330f8ffSToby Isaac DM dm; 132c330f8ffSToby Isaac 1338c6c5593SMatthew G. Knepley if (!sp[f]) continue; 1348c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 135c330f8ffSToby Isaac if (!funcs[f]) { 136be1504a2SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 137c330f8ffSToby Isaac continue; 138c330f8ffSToby Isaac } 139c330f8ffSToby Isaac ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 140c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 141c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 142c330f8ffSToby Isaac ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 1430c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 144c330f8ffSToby Isaac if (isAffine) { 1454bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x); 1461c531cf8SMatthew G. Knepley } else { 1474bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[tp*dE]; 1484bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[tp*dE*dE]; 1494bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[tp*dE*dE]; 1504bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[tp]; 1518c6c5593SMatthew G. Knepley } 152a8f1f9e5SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 153a8f1f9e5SMatthew G. Knepley if (probAux) {ierr = PetscFEEvaluateFieldJets_Internal(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 1544bee2e38SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);} 1554bee2e38SMatthew G. Knepley (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, numConstants, constants, &pointEval[Nc[f]*q]); 1561c531cf8SMatthew G. Knepley } 157c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 158c330f8ffSToby Isaac ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 159c330f8ffSToby Isaac v += spDim; 1608c6c5593SMatthew G. Knepley } 1618c6c5593SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 1628c6c5593SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 1638c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 1648c6c5593SMatthew G. Knepley } 1658c6c5593SMatthew G. Knepley 1664bee2e38SMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 167b18799e0SMatthew G. Knepley PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 168b18799e0SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 169b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 170b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 171b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 172b18799e0SMatthew G. Knepley PetscScalar values[]) 173b18799e0SMatthew G. Knepley { 174b18799e0SMatthew G. Knepley PetscSection section, sectionAux = NULL; 1754bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 176b18799e0SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 177b18799e0SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 178b18799e0SMatthew G. Knepley const PetscScalar *constants; 179b18799e0SMatthew G. Knepley PetscReal *x; 180b18799e0SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 1819f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 1824bee2e38SMatthew G. Knepley const PetscInt dE = fgeom->dimEmbed; 183b18799e0SMatthew G. Knepley PetscInt dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 184b18799e0SMatthew G. Knepley PetscBool isAffine; 185b18799e0SMatthew G. Knepley PetscErrorCode ierr; 186b18799e0SMatthew G. Knepley 187b18799e0SMatthew G. Knepley PetscFunctionBeginHot; 188b18799e0SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 189b18799e0SMatthew G. Knepley ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 190b18799e0SMatthew G. Knepley ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 191b18799e0SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 192b18799e0SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 193b18799e0SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 1944bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 195b18799e0SMatthew G. Knepley ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 196b18799e0SMatthew G. Knepley ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 197b18799e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 198b18799e0SMatthew G. Knepley if (dmAux) { 199b18799e0SMatthew G. Knepley DMLabel spmap; 200b18799e0SMatthew G. Knepley PetscInt subp = p; 201b18799e0SMatthew G. Knepley 202b18799e0SMatthew G. Knepley /* If dm is a submesh, do not get subpoint */ 203b18799e0SMatthew G. Knepley ierr = DMPlexGetSubpointMap(dm, &spmap);CHKERRQ(ierr); 204b18799e0SMatthew G. Knepley if (!spmap) {ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr);} 205b18799e0SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 206b18799e0SMatthew G. Knepley ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 207b18799e0SMatthew G. Knepley ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 208b18799e0SMatthew G. Knepley ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 209b18799e0SMatthew G. Knepley ierr = DMGetSection(dmAux, §ionAux);CHKERRQ(ierr); 210b18799e0SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 211b18799e0SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 212b18799e0SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 213b18799e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 214b18799e0SMatthew G. Knepley } 215b18799e0SMatthew G. Knepley /* Get values for closure */ 2164bee2e38SMatthew G. Knepley isAffine = fgeom->isAffine; 217*4dcf62a8SSatish Balay fegeom.n = 0; 218*4dcf62a8SSatish Balay fegeom.J = 0; 219*4dcf62a8SSatish Balay fegeom.v = 0; 220*4dcf62a8SSatish Balay fegeom.xi = 0; 2214bee2e38SMatthew G. Knepley if (isAffine) { 2224bee2e38SMatthew G. Knepley fegeom.v = x; 2234bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 2244bee2e38SMatthew G. Knepley fegeom.J = fgeom->J; 2254bee2e38SMatthew G. Knepley fegeom.invJ = fgeom->invJ; 2264bee2e38SMatthew G. Knepley fegeom.detJ = fgeom->detJ; 2274bee2e38SMatthew G. Knepley fegeom.n = fgeom->n; 2289f209ee4SMatthew G. Knepley 2299f209ee4SMatthew G. Knepley cgeom.J = fgeom->suppJ[0]; 2309f209ee4SMatthew G. Knepley cgeom.invJ = fgeom->suppInvJ[0]; 2319f209ee4SMatthew G. Knepley cgeom.detJ = fgeom->suppDetJ[0]; 2324bee2e38SMatthew G. Knepley } 233b18799e0SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 234b18799e0SMatthew G. Knepley PetscQuadrature allPoints; 235b18799e0SMatthew G. Knepley PetscInt q, dim, numPoints; 236b18799e0SMatthew G. Knepley const PetscReal *points; 237b18799e0SMatthew G. Knepley PetscScalar *pointEval; 238b18799e0SMatthew G. Knepley DM dm; 239b18799e0SMatthew G. Knepley 240b18799e0SMatthew G. Knepley if (!sp[f]) continue; 241b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 242b18799e0SMatthew G. Knepley if (!funcs[f]) { 243b18799e0SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 244b18799e0SMatthew G. Knepley continue; 245b18799e0SMatthew G. Knepley } 246b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 247b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 248b18799e0SMatthew G. Knepley ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 249b18799e0SMatthew G. Knepley ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 250b18799e0SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 251b18799e0SMatthew G. Knepley if (isAffine) { 2524bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x); 253b18799e0SMatthew G. Knepley } else { 2543fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[tp*dE]; 2559f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[tp*dE*dE]; 2569f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[tp*dE*dE]; 2574bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[tp]; 2584bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[tp*dE]; 2599f209ee4SMatthew G. Knepley 2609f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][tp*dE*dE]; 2619f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][tp*dE*dE]; 2629f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][tp]; 263b18799e0SMatthew G. Knepley } 2649f209ee4SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 2659f209ee4SMatthew G. Knepley if (probAux) {ierr = PetscFEEvaluateFieldJets_Internal(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 2664bee2e38SMatthew G. Knepley (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, fegeom.n, numConstants, constants, &pointEval[Nc[f]*q]); 267b18799e0SMatthew G. Knepley } 268b18799e0SMatthew G. Knepley ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 269b18799e0SMatthew G. Knepley ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 270b18799e0SMatthew G. Knepley v += spDim; 271b18799e0SMatthew G. Knepley } 272b18799e0SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 273b18799e0SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 274b18799e0SMatthew G. Knepley PetscFunctionReturn(0); 275b18799e0SMatthew G. Knepley } 276b18799e0SMatthew G. Knepley 277c330f8ffSToby Isaac static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS prob, PetscFEGeom *fegeom, DM dmAux, PetscDS probAux, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[], 2781c531cf8SMatthew G. Knepley PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[], 2791c531cf8SMatthew G. Knepley PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 2808c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 2818c6c5593SMatthew G. Knepley { 2828c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom; 2838c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed; 2848c6c5593SMatthew G. Knepley PetscErrorCode ierr; 2858c6c5593SMatthew G. Knepley 2868c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 2878c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2888c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 2898c6c5593SMatthew G. Knepley if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 2908c6c5593SMatthew G. Knepley switch (type) { 2918c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL: 2928c6c5593SMatthew G. Knepley case DM_BC_NATURAL: 293c330f8ffSToby Isaac 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; 2948c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 2958c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD: 296c330f8ffSToby Isaac ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps, 2970c898c18SMatthew G. Knepley basisTab, basisDerTab, basisTabAux, basisDerTabAux, 2980c898c18SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 2990c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 3000c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 301191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 302b18799e0SMatthew G. Knepley case DM_BC_ESSENTIAL_BD_FIELD: 303b18799e0SMatthew G. Knepley ierr = DMProjectPoint_BdField_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps, 304b18799e0SMatthew G. Knepley basisTab, basisDerTab, basisTabAux, basisDerTabAux, 305b18799e0SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 306b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 307b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 308b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 3098c6c5593SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 3108c6c5593SMatthew G. Knepley } 3118c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 3128c6c5593SMatthew G. Knepley } 3138c6c5593SMatthew G. Knepley 314df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 315df5c1128SToby Isaac { 316df5c1128SToby Isaac PetscReal *points; 317df5c1128SToby Isaac PetscInt f, numPoints; 318df5c1128SToby Isaac PetscErrorCode ierr; 319df5c1128SToby Isaac 320df5c1128SToby Isaac PetscFunctionBegin; 321df5c1128SToby Isaac numPoints = 0; 322df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 323df5c1128SToby Isaac if (funcs[f]) { 324df5c1128SToby Isaac PetscQuadrature fAllPoints; 325df5c1128SToby Isaac PetscInt fNumPoints; 326df5c1128SToby Isaac 327df5c1128SToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 328df5c1128SToby Isaac ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr); 329df5c1128SToby Isaac numPoints += fNumPoints; 330df5c1128SToby Isaac } 331df5c1128SToby Isaac } 332df5c1128SToby Isaac ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 333df5c1128SToby Isaac numPoints = 0; 334df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 335df5c1128SToby Isaac PetscInt spDim; 336df5c1128SToby Isaac 337df5c1128SToby Isaac ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 338df5c1128SToby Isaac if (funcs[f]) { 339df5c1128SToby Isaac PetscQuadrature fAllPoints; 34054ee0cceSMatthew G. Knepley PetscInt qdim, fNumPoints, q; 341df5c1128SToby Isaac const PetscReal *fPoints; 342df5c1128SToby Isaac 343df5c1128SToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 34454ee0cceSMatthew G. Knepley ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr); 34554ee0cceSMatthew G. Knepley if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %D for dual basis does not match input dimension %D", qdim, dim); 346df5c1128SToby Isaac for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q]; 347df5c1128SToby Isaac numPoints += fNumPoints; 348df5c1128SToby Isaac } 349df5c1128SToby Isaac } 350df5c1128SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 351df5c1128SToby Isaac ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 352df5c1128SToby Isaac PetscFunctionReturn(0); 353df5c1128SToby Isaac } 354df5c1128SToby Isaac 355e5e52638SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *prob) 356e5e52638SMatthew G. Knepley { 357e5e52638SMatthew G. Knepley DMLabel depthLabel; 358e5e52638SMatthew G. Knepley PetscInt dim, cdepth, ls = -1, i; 359e5e52638SMatthew G. Knepley PetscErrorCode ierr; 360e5e52638SMatthew G. Knepley 361e5e52638SMatthew G. Knepley PetscFunctionBegin; 362e5e52638SMatthew G. Knepley if (lStart) *lStart = -1; 363e5e52638SMatthew G. Knepley if (!label) PetscFunctionReturn(0); 364e5e52638SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 365e5e52638SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr); 366e5e52638SMatthew G. Knepley cdepth = dim - height; 367e5e52638SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 368e5e52638SMatthew G. Knepley IS pointIS; 369e5e52638SMatthew G. Knepley const PetscInt *points; 370e5e52638SMatthew G. Knepley PetscInt pdepth; 371e5e52638SMatthew G. Knepley 372e5e52638SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 373e5e52638SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 374e5e52638SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 375e5e52638SMatthew G. Knepley ierr = DMLabelGetValue(depthLabel, points[0], &pdepth);CHKERRQ(ierr); 376e5e52638SMatthew G. Knepley if (pdepth == cdepth) { 377e5e52638SMatthew G. Knepley ls = points[0]; 378e5e52638SMatthew G. Knepley if (prob) {ierr = DMGetCellDS(dm, ls, prob);CHKERRQ(ierr);} 379e5e52638SMatthew G. Knepley } 380e5e52638SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 381e5e52638SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 382e5e52638SMatthew G. Knepley if (ls >= 0) break; 383e5e52638SMatthew G. Knepley } 384e5e52638SMatthew G. Knepley if (lStart) *lStart = ls; 385e5e52638SMatthew G. Knepley PetscFunctionReturn(0); 386e5e52638SMatthew G. Knepley } 387e5e52638SMatthew G. Knepley 3880de51b56SMatthew G. Knepley /* 3890de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 3900de51b56SMatthew G. Knepley 3910de51b56SMatthew G. Knepley There are several different scenarios: 3920de51b56SMatthew G. Knepley 3930de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data 3940de51b56SMatthew G. Knepley 3950de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells. 3960de51b56SMatthew G. Knepley 3970de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data 3980de51b56SMatthew G. Knepley 3990de51b56SMatthew 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. 4000de51b56SMatthew G. Knepley 4010de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data 4020de51b56SMatthew G. Knepley 4030de51b56SMatthew 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. 4040de51b56SMatthew G. Knepley 4050de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest. 4060de51b56SMatthew G. Knepley 4070de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 4080de51b56SMatthew G. Knepley 4090de51b56SMatthew 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. 4100de51b56SMatthew 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 4110de51b56SMatthew 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 4120de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces. 4130de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 4140de51b56SMatthew G. Knepley 4150de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 4160de51b56SMatthew G. Knepley 4170de51b56SMatthew 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. 4180de51b56SMatthew G. Knepley */ 41946fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 4201c531cf8SMatthew G. Knepley PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], 4218c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 4228c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 4238c6c5593SMatthew G. Knepley { 424ca3d3a14SMatthew G. Knepley DM dmAux = NULL, tdm; 425e5e52638SMatthew G. Knepley PetscDS prob = NULL, probAux = NULL; 426ca3d3a14SMatthew G. Knepley Vec localA = NULL, tv; 42747923291SMatthew G. Knepley PetscSection section; 4288c6c5593SMatthew G. Knepley PetscDualSpace *sp, *cellsp; 4290c898c18SMatthew G. Knepley PetscReal **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL; 4308c6c5593SMatthew G. Knepley PetscInt *Nc; 431c330f8ffSToby Isaac PetscInt dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f; 432ca3d3a14SMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, transform; 433c330f8ffSToby Isaac DMField coordField; 434c330f8ffSToby Isaac DMLabel depthLabel; 435c330f8ffSToby Isaac PetscQuadrature allPoints = NULL; 43647923291SMatthew G. Knepley PetscErrorCode ierr; 43747923291SMatthew G. Knepley 43847923291SMatthew G. Knepley PetscFunctionBegin; 43988aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 44088aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr); 4418c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 4422af58022SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr); 443ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr); 444ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr); 445ca3d3a14SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 4460de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */ 4470de51b56SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 448e39fcb4eSStefano Zampini if (dmAux && !localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector"); 44988aca1feSMatthew G. Knepley if (!minHeight && dmAux) { 45088aca1feSMatthew G. Knepley DMLabel spmap; 45188aca1feSMatthew G. Knepley 45288aca1feSMatthew G. Knepley /* If dmAux is a surface, then force the projection to take place over a surface */ 45388aca1feSMatthew G. Knepley ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr); 454e39fcb4eSStefano Zampini if (spmap) { 455e39fcb4eSStefano Zampini ierr = DMPlexGetVTKCellHeight(dmAux, &minHeight);CHKERRQ(ierr); 456e39fcb4eSStefano Zampini auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE; 457e39fcb4eSStefano Zampini } 45888aca1feSMatthew G. Knepley } 4590de51b56SMatthew G. Knepley } 460e5e52638SMatthew G. Knepley ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 461e5e52638SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr); 4628c6c5593SMatthew G. Knepley ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr); 4632af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight); 464e39fcb4eSStefano Zampini if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim); 465e5e52638SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, 0, NULL, &prob);CHKERRQ(ierr); 466e5e52638SMatthew G. Knepley if (!prob) {ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);} 467e5e52638SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 4688c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 469e87a4003SBarry Smith ierr = DMGetSection(dm, §ion);CHKERRQ(ierr); 4700c898c18SMatthew G. Knepley if (dmAux) { 4710c898c18SMatthew G. Knepley ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 4720c898c18SMatthew G. Knepley ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 4730c898c18SMatthew G. Knepley } 474496733e2SMatthew G. Knepley ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 475496733e2SMatthew G. Knepley ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr); 4768c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);} 4778c6c5593SMatthew G. Knepley else {cellsp = sp;} 4788c6c5593SMatthew G. Knepley if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 4798c6c5593SMatthew G. Knepley /* Get cell dual spaces */ 48047923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 48147923291SMatthew G. Knepley PetscObject obj; 48247923291SMatthew G. Knepley PetscClassId id; 48347923291SMatthew G. Knepley 484e5e52638SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 48547923291SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 48647923291SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 48747923291SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 48847923291SMatthew G. Knepley 48947923291SMatthew G. Knepley hasFE = PETSC_TRUE; 49047923291SMatthew G. Knepley isFE[f] = PETSC_TRUE; 49147923291SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 49247923291SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 49347923291SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 49447923291SMatthew G. Knepley 49547923291SMatthew G. Knepley hasFV = PETSC_TRUE; 49647923291SMatthew G. Knepley isFE[f] = PETSC_FALSE; 49747923291SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 49847923291SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 49947923291SMatthew G. Knepley } 500c330f8ffSToby Isaac ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr); 5010c898c18SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 50274fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 50374fc349aSMatthew G. Knepley PetscFE fem, subfem; 5044a3e9fdbSToby Isaac const PetscReal *points; 5054a3e9fdbSToby Isaac PetscInt numPoints; 5060c898c18SMatthew G. Knepley 5072af58022SMatthew G. Knepley if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 50854ee0cceSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 50954ee0cceSMatthew G. Knepley if (!effectiveHeight) {sp[f] = cellsp[f];} 51054ee0cceSMatthew G. Knepley else {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);} 51154ee0cceSMatthew G. Knepley } 51254ee0cceSMatthew G. Knepley ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr); 513c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 5140c898c18SMatthew G. Knepley ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr); 5150c898c18SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 5160c898c18SMatthew G. Knepley if (!isFE[f]) continue; 5170c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 51874fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 51974fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 52074fc349aSMatthew G. Knepley ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 5210c898c18SMatthew G. Knepley } 5220c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 5230c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 52474fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 52574fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 52674fc349aSMatthew G. Knepley ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 5270c898c18SMatthew G. Knepley } 5280c898c18SMatthew G. Knepley } 52947923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 5302af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) { 53188aca1feSMatthew G. Knepley PetscInt effectiveHeight = h - (auxBd ? 0 : minHeight); 5320de51b56SMatthew G. Knepley PetscDS probEff = prob; 5338c6c5593SMatthew G. Knepley PetscScalar *values; 534b7260050SToby Isaac PetscBool *fieldActive; 535b7260050SToby Isaac PetscInt maxDegree; 536e5e52638SMatthew G. Knepley PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 537c330f8ffSToby Isaac IS heightIS; 5388c6c5593SMatthew G. Knepley 53947923291SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr); 540e5e52638SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr); 541c330f8ffSToby Isaac ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr); 5428c6c5593SMatthew G. Knepley if (!h) { 5438c6c5593SMatthew G. Knepley PetscInt cEndInterior; 5448c6c5593SMatthew G. Knepley 5458c6c5593SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 5468c6c5593SMatthew G. Knepley pEnd = cEndInterior < 0 ? pEnd : cEndInterior; 5478c6c5593SMatthew G. Knepley } 548c330f8ffSToby Isaac if (pEnd <= pStart) { 549c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 550c330f8ffSToby Isaac continue; 551c330f8ffSToby Isaac } 55247923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */ 55347923291SMatthew G. Knepley totDim = 0; 55447923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 555bec1f4aaSMatthew G. Knepley if (!effectiveHeight) { 55647923291SMatthew G. Knepley sp[f] = cellsp[f]; 55747923291SMatthew G. Knepley } else { 558bec1f4aaSMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr); 55947923291SMatthew G. Knepley if (!sp[f]) continue; 56047923291SMatthew G. Knepley } 56147923291SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 5629c3cf19fSMatthew G. Knepley totDim += spDim; 56347923291SMatthew G. Knepley } 564e5e52638SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr); 5658c6c5593SMatthew 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); 566c330f8ffSToby Isaac if (!totDim) { 567c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 568c330f8ffSToby Isaac continue; 569c330f8ffSToby Isaac } 570c9731667SStefano Zampini if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);} 57147923291SMatthew G. Knepley /* Loop over points at this height */ 57269291d52SBarry Smith ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 57369291d52SBarry Smith ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 5748c6c5593SMatthew G. Knepley for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 5758c6c5593SMatthew G. Knepley if (label) { 5768c6c5593SMatthew G. Knepley PetscInt i; 57747923291SMatthew G. Knepley 57847923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 579c330f8ffSToby Isaac IS pointIS, isectIS; 58047923291SMatthew G. Knepley const PetscInt *points; 5818c6c5593SMatthew G. Knepley PetscInt n; 582c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 583c330f8ffSToby Isaac PetscQuadrature quad = NULL; 58447923291SMatthew G. Knepley 58547923291SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 58647923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 587c330f8ffSToby Isaac ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr); 588c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 589c330f8ffSToby Isaac if (!isectIS) continue; 590c330f8ffSToby Isaac ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr); 591c330f8ffSToby Isaac ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr); 592b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr); 593b7260050SToby Isaac if (maxDegree <= 1) { 594c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr); 595c330f8ffSToby Isaac } 596c330f8ffSToby Isaac if (!quad) { 597c330f8ffSToby Isaac if (!h && allPoints) { 598c330f8ffSToby Isaac quad = allPoints; 599c330f8ffSToby Isaac allPoints = NULL; 600c330f8ffSToby Isaac } else { 601c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 602c330f8ffSToby Isaac } 603c330f8ffSToby Isaac } 604c330f8ffSToby Isaac ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr); 60547923291SMatthew G. Knepley for (p = 0; p < n; ++p) { 60647923291SMatthew G. Knepley const PetscInt point = points[p]; 60747923291SMatthew G. Knepley 608b420b6ccSMatthew G. Knepley ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 609c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 610c330f8ffSToby Isaac ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 61147923291SMatthew G. Knepley if (ierr) { 61247923291SMatthew G. Knepley PetscErrorCode ierr2; 61369291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 61469291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 61547923291SMatthew G. Knepley CHKERRQ(ierr); 61647923291SMatthew G. Knepley } 617ca3d3a14SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 618ba322698SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr); 61947923291SMatthew G. Knepley } 620c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 621c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 622c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 623c330f8ffSToby Isaac ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr); 624c330f8ffSToby Isaac ierr = ISDestroy(&isectIS);CHKERRQ(ierr); 62547923291SMatthew G. Knepley } 6268c6c5593SMatthew G. Knepley } else { 627c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 628c330f8ffSToby Isaac PetscQuadrature quad = NULL; 629c330f8ffSToby Isaac IS pointIS; 630c330f8ffSToby Isaac 631c330f8ffSToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr); 632b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 633b7260050SToby Isaac if (maxDegree <= 1) { 634c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr); 635c330f8ffSToby Isaac } 636c330f8ffSToby Isaac if (!quad) { 637c330f8ffSToby Isaac if (!h && allPoints) { 638c330f8ffSToby Isaac quad = allPoints; 639c330f8ffSToby Isaac allPoints = NULL; 640c330f8ffSToby Isaac } else { 641c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 642c330f8ffSToby Isaac } 643c330f8ffSToby Isaac } 644c330f8ffSToby Isaac ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr); 6458c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 646b420b6ccSMatthew G. Knepley ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr); 647c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr); 648c330f8ffSToby Isaac ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 6498c6c5593SMatthew G. Knepley if (ierr) { 6508c6c5593SMatthew G. Knepley PetscErrorCode ierr2; 65169291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 65269291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 6538c6c5593SMatthew G. Knepley CHKERRQ(ierr); 6548c6c5593SMatthew G. Knepley } 655ca3d3a14SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 656ba322698SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr); 6578c6c5593SMatthew G. Knepley } 658c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr); 659c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 660c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 661c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 6628c6c5593SMatthew G. Knepley } 663c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 66469291d52SBarry Smith ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 66569291d52SBarry Smith ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 66647923291SMatthew G. Knepley } 6678c6c5593SMatthew G. Knepley /* Cleanup */ 6680c898c18SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 66974fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 67074fc349aSMatthew G. Knepley PetscFE fem, subfem; 6710c898c18SMatthew G. Knepley 6720c898c18SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 6730c898c18SMatthew G. Knepley if (!isFE[f]) continue; 6740c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr); 67574fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 67674fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 67774fc349aSMatthew G. Knepley ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 6780c898c18SMatthew G. Knepley } 6790c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 6800c898c18SMatthew G. Knepley ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 68174fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 68274fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 68374fc349aSMatthew G. Knepley ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 6840c898c18SMatthew G. Knepley } 6850c898c18SMatthew G. Knepley ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr); 6860c898c18SMatthew G. Knepley } 687c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr); 688496733e2SMatthew G. Knepley ierr = PetscFree2(isFE, sp);CHKERRQ(ierr); 6898c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);} 6908c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 69147923291SMatthew G. Knepley } 6928c6c5593SMatthew G. Knepley 6938c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 6948c6c5593SMatthew G. Knepley { 6958c6c5593SMatthew G. Knepley PetscErrorCode ierr; 6968c6c5593SMatthew G. Knepley 6978c6c5593SMatthew G. Knepley PetscFunctionBegin; 6981c531cf8SMatthew 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); 6998c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 7008c6c5593SMatthew G. Knepley } 7018c6c5593SMatthew G. Knepley 7021c531cf8SMatthew 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) 7038c6c5593SMatthew G. Knepley { 7048c6c5593SMatthew G. Knepley PetscErrorCode ierr; 7058c6c5593SMatthew G. Knepley 7068c6c5593SMatthew G. Knepley PetscFunctionBegin; 7071c531cf8SMatthew 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); 70847923291SMatthew G. Knepley PetscFunctionReturn(0); 70947923291SMatthew G. Knepley } 71047923291SMatthew G. Knepley 7118c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 71247923291SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 71347923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 71447923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 715191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 71647923291SMatthew G. Knepley InsertMode mode, Vec localX) 71747923291SMatthew G. Knepley { 71847923291SMatthew G. Knepley PetscErrorCode ierr; 71947923291SMatthew G. Knepley 72047923291SMatthew G. Knepley PetscFunctionBegin; 7211c531cf8SMatthew 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); 7228c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 72347923291SMatthew G. Knepley } 72447923291SMatthew G. Knepley 7251c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 7268c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 7278c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 7288c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 729191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 7308c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 7318c6c5593SMatthew G. Knepley { 7328c6c5593SMatthew G. Knepley PetscErrorCode ierr; 73347923291SMatthew G. Knepley 7348c6c5593SMatthew G. Knepley PetscFunctionBegin; 7351c531cf8SMatthew 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); 73647923291SMatthew G. Knepley PetscFunctionReturn(0); 73747923291SMatthew G. Knepley } 738