xref: /petsc/src/dm/impls/plex/plexproject.c (revision c330f8ffacf2fcaab5f14247658e3a46611746d8)
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*c330f8ffSToby 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 {
9*c330f8ffSToby Isaac   PetscInt       coordDim, Nf, *Nc, f, totDim, spDim, d, v, tp;
10*c330f8ffSToby Isaac   PetscBool      isAffine;
118c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
128c6c5593SMatthew G. Knepley 
138c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
14*c330f8ffSToby Isaac   ierr = DMGetCoordinateDim(dm,&coordDim);CHKERRQ(ierr);
158c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
16496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
179c3cf19fSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
188c6c5593SMatthew G. Knepley   /* Get values for closure */
19*c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
20*c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
218c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
228c6c5593SMatthew G. Knepley 
238c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
248c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
258c6c5593SMatthew G. Knepley     if (funcs[f]) {
26*c330f8ffSToby Isaac       if (isFE[f]) {
27*c330f8ffSToby Isaac         PetscQuadrature   allPoints;
28*c330f8ffSToby Isaac         PetscInt          q, dim, numPoints;
29*c330f8ffSToby Isaac         const PetscReal   *points;
30*c330f8ffSToby Isaac         PetscScalar       *pointEval;
31*c330f8ffSToby Isaac         PetscReal         *x;
32*c330f8ffSToby Isaac         DM                dm;
33*c330f8ffSToby Isaac 
34*c330f8ffSToby Isaac         ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
35*c330f8ffSToby Isaac         ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
36*c330f8ffSToby Isaac         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
37*c330f8ffSToby Isaac         ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
38*c330f8ffSToby Isaac         ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
39*c330f8ffSToby Isaac         ierr = DMGetWorkArray(dm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
40*c330f8ffSToby Isaac         for (q = 0; q < numPoints; q++, tp++) {
41*c330f8ffSToby Isaac           const PetscReal *v0;
42*c330f8ffSToby Isaac 
43*c330f8ffSToby Isaac           if (isAffine) {
44*c330f8ffSToby Isaac             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
45*c330f8ffSToby Isaac             v0 = x;
468c6c5593SMatthew G. Knepley           } else {
47*c330f8ffSToby Isaac             v0 = &fegeom->v[tp*coordDim];
488c6c5593SMatthew G. Knepley           }
49*c330f8ffSToby Isaac           ierr = (*funcs[f])(coordDim,time,v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
50*c330f8ffSToby Isaac         }
51*c330f8ffSToby Isaac         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
52*c330f8ffSToby Isaac         ierr = DMRestoreWorkArray(dm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
53*c330f8ffSToby Isaac         ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
54*c330f8ffSToby Isaac         v += spDim;
55*c330f8ffSToby Isaac       } else {
56*c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
57*c330f8ffSToby Isaac           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
58*c330f8ffSToby Isaac         }
59*c330f8ffSToby Isaac       }
60*c330f8ffSToby Isaac     } else {
61*c330f8ffSToby Isaac       for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;}
628c6c5593SMatthew G. Knepley     }
639c3cf19fSMatthew G. Knepley   }
648c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
658c6c5593SMatthew G. Knepley }
668c6c5593SMatthew G. Knepley 
67*c330f8ffSToby Isaac static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fegeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
680c898c18SMatthew G. Knepley                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
698c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
708c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
718c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
72191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
738c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
748c6c5593SMatthew G. Knepley {
758c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
769c3cf19fSMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL, *bc;
778c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
788c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
79191494d9SMatthew G. Knepley   const PetscScalar *constants;
808c6c5593SMatthew G. Knepley   PetscReal         *x;
81496733e2SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
82*c330f8ffSToby Isaac   const PetscInt     dE = fegeom->dimEmbed;
83*c330f8ffSToby Isaac   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
84*c330f8ffSToby Isaac   PetscBool          isAffine;
858c6c5593SMatthew G. Knepley   PetscErrorCode     ierr;
868c6c5593SMatthew G. Knepley 
878c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
888c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
89496733e2SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
90496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
918c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
928c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
939c3cf19fSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
94496733e2SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
95191494d9SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
968c6c5593SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
97496733e2SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
988c6c5593SMatthew G. Knepley   if (dmAux) {
991ba34526SMatthew G. Knepley     DMLabel  spmap;
1001ba34526SMatthew G. Knepley     PetscInt subp;
1011ba34526SMatthew G. Knepley 
1021ba34526SMatthew G. Knepley     ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
1031ba34526SMatthew G. Knepley     if (spmap) {
1041ba34526SMatthew G. Knepley       IS              subpointIS;
1051ba34526SMatthew G. Knepley       const PetscInt *subpoints;
1061ba34526SMatthew G. Knepley       PetscInt        numSubpoints;
1071ba34526SMatthew G. Knepley 
1081ba34526SMatthew G. Knepley       ierr = DMPlexCreateSubpointIS(dmAux, &subpointIS);CHKERRQ(ierr);
1091ba34526SMatthew G. Knepley       ierr = ISGetLocalSize(subpointIS, &numSubpoints);CHKERRQ(ierr);
1101ba34526SMatthew G. Knepley       ierr = ISGetIndices(subpointIS, &subpoints);CHKERRQ(ierr);
1111ba34526SMatthew G. Knepley       ierr = PetscFindInt(p, numSubpoints, subpoints, &subp);CHKERRQ(ierr);
1121ba34526SMatthew G. Knepley       if (subp < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d not found in submesh", p);
1131ba34526SMatthew G. Knepley       ierr = ISRestoreIndices(subpointIS, &subpoints);CHKERRQ(ierr);
1141ba34526SMatthew G. Knepley       ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
1151ba34526SMatthew G. Knepley     } else subp = p;
1161ba34526SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
1178c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
118496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
119496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
1208c6c5593SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1218c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
1228c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
12311d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
124496733e2SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
1251ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
1268c6c5593SMatthew G. Knepley   }
1278c6c5593SMatthew G. Knepley   /* Get values for closure */
128*c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
1298c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
130*c330f8ffSToby Isaac     PetscQuadrature   allPoints;
131*c330f8ffSToby Isaac     PetscInt          q, dim, numPoints;
132*c330f8ffSToby Isaac     const PetscReal   *points;
133*c330f8ffSToby Isaac     PetscScalar       *pointEval;
134*c330f8ffSToby Isaac     DM                dm;
135*c330f8ffSToby Isaac 
1368c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
1378c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
138*c330f8ffSToby Isaac     if (!funcs[f]) {
139*c330f8ffSToby Isaac       for (d = 0; d < spDim; d++, v++) {
140*c330f8ffSToby Isaac         values[v] = 0.;
141*c330f8ffSToby Isaac       }
142*c330f8ffSToby Isaac       continue;
143*c330f8ffSToby Isaac     }
144*c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
145*c330f8ffSToby Isaac     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
146*c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
147*c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
1480c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
149*c330f8ffSToby Isaac       const PetscReal *v0;
150*c330f8ffSToby Isaac       const PetscReal *invJ;
151*c330f8ffSToby Isaac 
152*c330f8ffSToby Isaac       if (isAffine) {
153*c330f8ffSToby Isaac         CoordinatesRefToReal(dim, dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
154*c330f8ffSToby Isaac         v0 = x;
155*c330f8ffSToby Isaac         invJ = fegeom->invJ;
1561c531cf8SMatthew G. Knepley       } else {
157*c330f8ffSToby Isaac         v0 = &fegeom->v[tp*dE];
158*c330f8ffSToby Isaac         invJ = &fegeom->invJ[tp*dE*dE];
1598c6c5593SMatthew G. Knepley       }
160*c330f8ffSToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, invJ, coefficients, coefficients_t, u, u_x, u_t);
161*c330f8ffSToby Isaac       if (probAux) {EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);}
162*c330f8ffSToby Isaac       (*funcs[f])(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, v0, numConstants, constants, &pointEval[Nc[f]*q]);
1631c531cf8SMatthew G. Knepley     }
164*c330f8ffSToby Isaac     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
165*c330f8ffSToby Isaac     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
166*c330f8ffSToby Isaac     v += spDim;
1678c6c5593SMatthew G. Knepley   }
1688c6c5593SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1698c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
1708c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1718c6c5593SMatthew G. Knepley }
1728c6c5593SMatthew G. Knepley 
173*c330f8ffSToby 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[],
1741c531cf8SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
1751c531cf8SMatthew G. Knepley                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
1768c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
1778c6c5593SMatthew G. Knepley {
1788c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
1798c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
1808c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
1818c6c5593SMatthew G. Knepley 
1828c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1838c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1848c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1858c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
1868c6c5593SMatthew G. Knepley   switch (type) {
1878c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
1888c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
189*c330f8ffSToby 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;
1908c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
1918c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
192*c330f8ffSToby Isaac     ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
1930c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
1940c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
1950c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1960c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
197191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
1988c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
1998c6c5593SMatthew G. Knepley   }
2008c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2018c6c5593SMatthew G. Knepley }
2028c6c5593SMatthew G. Knepley 
203df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
204df5c1128SToby Isaac {
205df5c1128SToby Isaac   PetscReal      *points;
206df5c1128SToby Isaac   PetscInt       f, numPoints;
207df5c1128SToby Isaac   PetscErrorCode ierr;
208df5c1128SToby Isaac 
209df5c1128SToby Isaac   PetscFunctionBegin;
210df5c1128SToby Isaac   numPoints = 0;
211df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
212df5c1128SToby Isaac     if (funcs[f]) {
213df5c1128SToby Isaac       PetscQuadrature fAllPoints;
214df5c1128SToby Isaac       PetscInt        fNumPoints;
215df5c1128SToby Isaac 
216df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
217df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
218df5c1128SToby Isaac       numPoints += fNumPoints;
219df5c1128SToby Isaac     }
220df5c1128SToby Isaac   }
221df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
222df5c1128SToby Isaac   numPoints = 0;
223df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
224df5c1128SToby Isaac     PetscInt spDim;
225df5c1128SToby Isaac 
226df5c1128SToby Isaac     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
227df5c1128SToby Isaac     if (funcs[f]) {
228df5c1128SToby Isaac       PetscQuadrature fAllPoints;
229df5c1128SToby Isaac       PetscInt        fNumPoints, q;
230df5c1128SToby Isaac       const PetscReal *fPoints;
231df5c1128SToby Isaac 
232df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
233df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
234df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
235df5c1128SToby Isaac       numPoints += fNumPoints;
236df5c1128SToby Isaac     }
237df5c1128SToby Isaac   }
238df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
239df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
240df5c1128SToby Isaac   PetscFunctionReturn(0);
241df5c1128SToby Isaac }
242df5c1128SToby Isaac 
2430de51b56SMatthew G. Knepley /*
2440de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
2450de51b56SMatthew G. Knepley 
2460de51b56SMatthew G. Knepley   There are several different scenarios:
2470de51b56SMatthew G. Knepley 
2480de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
2490de51b56SMatthew G. Knepley 
2500de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
2510de51b56SMatthew G. Knepley 
2520de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
2530de51b56SMatthew G. Knepley 
2540de51b56SMatthew 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.
2550de51b56SMatthew G. Knepley 
2560de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
2570de51b56SMatthew G. Knepley 
2580de51b56SMatthew 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.
2590de51b56SMatthew G. Knepley 
2600de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
2610de51b56SMatthew G. Knepley 
2620de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
2630de51b56SMatthew G. Knepley 
2640de51b56SMatthew 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.
2650de51b56SMatthew 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
2660de51b56SMatthew 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
2670de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
2680de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
2690de51b56SMatthew G. Knepley 
2700de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
2710de51b56SMatthew G. Knepley 
2720de51b56SMatthew 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.
2730de51b56SMatthew G. Knepley */
27446fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
2751c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
2768c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
2778c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
2788c6c5593SMatthew G. Knepley {
2798c6c5593SMatthew G. Knepley   DM              dmAux = NULL;
2800c898c18SMatthew G. Knepley   PetscDS         prob, probAux = NULL;
2818c6c5593SMatthew G. Knepley   Vec             localA = NULL;
28247923291SMatthew G. Knepley   PetscSection    section;
2838c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
2840c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
2858c6c5593SMatthew G. Knepley   PetscInt       *Nc;
286*c330f8ffSToby Isaac   PetscInt        dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f;
28788aca1feSMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE;
288*c330f8ffSToby Isaac   DMField         coordField;
289*c330f8ffSToby Isaac   DMLabel         depthLabel;
290*c330f8ffSToby Isaac   PetscQuadrature allPoints = NULL;
29147923291SMatthew G. Knepley   PetscErrorCode  ierr;
29247923291SMatthew G. Knepley 
29347923291SMatthew G. Knepley   PetscFunctionBegin;
29488aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
29588aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
2968c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2972af58022SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
2980de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
2990de51b56SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
30088aca1feSMatthew G. Knepley     if (!minHeight && dmAux) {
30188aca1feSMatthew G. Knepley       DMLabel spmap;
30288aca1feSMatthew G. Knepley 
30388aca1feSMatthew G. Knepley       /* If dmAux is a surface, then force the projection to take place over a surface */
30488aca1feSMatthew G. Knepley       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
30588aca1feSMatthew G. Knepley       if (spmap) {minHeight = 1; auxBd = PETSC_TRUE;}
30688aca1feSMatthew G. Knepley     }
3070de51b56SMatthew G. Knepley   }
3088c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
3092af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
3108c6c5593SMatthew 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);}
3110c898c18SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
3128c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
31347923291SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
31447923291SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3150c898c18SMatthew G. Knepley   if (dmAux) {
3160c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
3170c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
3180c898c18SMatthew G. Knepley   }
319496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
320496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
3218c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
3228c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
3238c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
3248c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
32547923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
32647923291SMatthew G. Knepley     PetscObject  obj;
32747923291SMatthew G. Knepley     PetscClassId id;
32847923291SMatthew G. Knepley 
32947923291SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
33047923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
33147923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
33247923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
33347923291SMatthew G. Knepley 
33447923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
33547923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
33647923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
33747923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
33847923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
33947923291SMatthew G. Knepley 
34047923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
34147923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
34247923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
34347923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
34447923291SMatthew G. Knepley   }
345*c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
3460c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
34774fc349aSMatthew G. Knepley     PetscInt   effectiveHeight = auxBd ? minHeight : 0;
34874fc349aSMatthew G. Knepley     PetscFE    fem, subfem;
3490c898c18SMatthew G. Knepley     PetscReal *points;
35074fc349aSMatthew G. Knepley     PetscInt   numPoints, spDim, qdim = 0, d;
3510c898c18SMatthew G. Knepley 
3522af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
353*c330f8ffSToby Isaac     ierr = PetscDualSpaceGetAllPointsUnion(Nf,cellsp,dim,funcs,&allPoints);CHKERRQ(ierr);
354*c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
3550c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
3560c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
3570c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
3580c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
35974fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
36074fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
36174fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
3620c898c18SMatthew G. Knepley     }
3630c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
3640c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
36574fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
36674fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
36774fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
3680c898c18SMatthew G. Knepley     }
3690c898c18SMatthew G. Knepley   }
370*c330f8ffSToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
371*c330f8ffSToby Isaac   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
37247923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
3732af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
37488aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
3750de51b56SMatthew G. Knepley     PetscDS      probEff         = prob;
3768c6c5593SMatthew G. Knepley     PetscScalar *values;
377*c330f8ffSToby Isaac     PetscBool   *fieldActive, isAffine;
37884f7fa7aSMatthew G. Knepley     PetscInt     pStart, pEnd, p, spDim, totDim, numValues;
379*c330f8ffSToby Isaac     IS           heightIS;
3808c6c5593SMatthew G. Knepley 
3810de51b56SMatthew G. Knepley     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);}
38247923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
383*c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel,depth - h,&heightIS);CHKERRQ(ierr);
3848c6c5593SMatthew G. Knepley     if (!h) {
3858c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
3868c6c5593SMatthew G. Knepley 
3878c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3888c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
3898c6c5593SMatthew G. Knepley     }
390*c330f8ffSToby Isaac     if (pEnd <= pStart) {
391*c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
392*c330f8ffSToby Isaac       continue;
393*c330f8ffSToby Isaac     }
39447923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
39547923291SMatthew G. Knepley     totDim = 0;
39647923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
397bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
39847923291SMatthew G. Knepley         sp[f] = cellsp[f];
39947923291SMatthew G. Knepley       } else {
400bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
40147923291SMatthew G. Knepley         if (!sp[f]) continue;
40247923291SMatthew G. Knepley       }
40347923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
4049c3cf19fSMatthew G. Knepley       totDim += spDim;
40547923291SMatthew G. Knepley     }
40647923291SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
4078c6c5593SMatthew 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);
408*c330f8ffSToby Isaac     if (!totDim) {
409*c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
410*c330f8ffSToby Isaac       continue;
411*c330f8ffSToby Isaac     }
41247923291SMatthew G. Knepley     /* Loop over points at this height */
41369291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
41469291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
4158c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
4168c6c5593SMatthew G. Knepley     if (label) {
4178c6c5593SMatthew G. Knepley       PetscInt i;
41847923291SMatthew G. Knepley 
41947923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
420*c330f8ffSToby Isaac         IS              pointIS, isectIS;
42147923291SMatthew G. Knepley         const PetscInt *points;
4228c6c5593SMatthew G. Knepley         PetscInt        n;
423*c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
424*c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
42547923291SMatthew G. Knepley 
42647923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
42747923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
428*c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
429*c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
430*c330f8ffSToby Isaac         if (!isectIS) continue;
431*c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
432*c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
433*c330f8ffSToby Isaac         ierr = DMFieldGetFEInvariance(coordField,isectIS,NULL,&isAffine,NULL);CHKERRQ(ierr);
434*c330f8ffSToby Isaac         if (isAffine) {
435*c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
436*c330f8ffSToby Isaac         }
437*c330f8ffSToby Isaac         if (!quad) {
438*c330f8ffSToby Isaac           if (!h && allPoints) {
439*c330f8ffSToby Isaac             quad = allPoints;
440*c330f8ffSToby Isaac             allPoints = NULL;
441*c330f8ffSToby Isaac           } else {
442*c330f8ffSToby Isaac             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
443*c330f8ffSToby Isaac           }
444*c330f8ffSToby Isaac         }
445*c330f8ffSToby Isaac         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
44647923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
44747923291SMatthew G. Knepley           const PetscInt  point = points[p];
44847923291SMatthew G. Knepley 
449b420b6ccSMatthew G. Knepley           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
450*c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
451*c330f8ffSToby 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);
45247923291SMatthew G. Knepley           if (ierr) {
45347923291SMatthew G. Knepley             PetscErrorCode ierr2;
45469291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
45569291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
45647923291SMatthew G. Knepley             CHKERRQ(ierr);
45747923291SMatthew G. Knepley           }
458ba322698SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
45947923291SMatthew G. Knepley         }
460*c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
461*c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
462*c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
463*c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
464*c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
46547923291SMatthew G. Knepley       }
4668c6c5593SMatthew G. Knepley     } else {
467*c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
468*c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
469*c330f8ffSToby Isaac       IS              pointIS;
470*c330f8ffSToby Isaac 
471*c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
472*c330f8ffSToby Isaac       ierr = DMFieldGetFEInvariance(coordField,pointIS,NULL,&isAffine,NULL);CHKERRQ(ierr);
473*c330f8ffSToby Isaac       if (isAffine) {
474*c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
475*c330f8ffSToby Isaac       }
476*c330f8ffSToby Isaac       if (!quad) {
477*c330f8ffSToby Isaac         if (!h && allPoints) {
478*c330f8ffSToby Isaac           quad = allPoints;
479*c330f8ffSToby Isaac           allPoints = NULL;
480*c330f8ffSToby Isaac         } else {
481*c330f8ffSToby Isaac           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
482*c330f8ffSToby Isaac         }
483*c330f8ffSToby Isaac       }
484*c330f8ffSToby Isaac       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
4858c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
486b420b6ccSMatthew G. Knepley         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
487*c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
488*c330f8ffSToby 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);
4898c6c5593SMatthew G. Knepley         if (ierr) {
4908c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
49169291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
49269291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
4938c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
4948c6c5593SMatthew G. Knepley         }
495ba322698SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
4968c6c5593SMatthew G. Knepley       }
497*c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
498*c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
499*c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
500*c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
5018c6c5593SMatthew G. Knepley     }
502*c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
50369291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
50469291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
50547923291SMatthew G. Knepley   }
5068c6c5593SMatthew G. Knepley   /* Cleanup */
5070c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
50874fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
50974fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
5100c898c18SMatthew G. Knepley 
5110c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
5120c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
5130c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
51474fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
51574fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
51674fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
5170c898c18SMatthew G. Knepley     }
5180c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
5190c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
52074fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
52174fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
52274fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
5230c898c18SMatthew G. Knepley     }
5240c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
5250c898c18SMatthew G. Knepley   }
526*c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
527496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
5288c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
5298c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
53047923291SMatthew G. Knepley }
5318c6c5593SMatthew G. Knepley 
5328c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
5338c6c5593SMatthew G. Knepley {
5348c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
5358c6c5593SMatthew G. Knepley 
5368c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5371c531cf8SMatthew 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);
5388c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
5398c6c5593SMatthew G. Knepley }
5408c6c5593SMatthew G. Knepley 
5411c531cf8SMatthew 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)
5428c6c5593SMatthew G. Knepley {
5438c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
5448c6c5593SMatthew G. Knepley 
5458c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5461c531cf8SMatthew 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);
54747923291SMatthew G. Knepley   PetscFunctionReturn(0);
54847923291SMatthew G. Knepley }
54947923291SMatthew G. Knepley 
5508c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
55147923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
55247923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
55347923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
554191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
55547923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
55647923291SMatthew G. Knepley {
55747923291SMatthew G. Knepley   PetscErrorCode ierr;
55847923291SMatthew G. Knepley 
55947923291SMatthew G. Knepley   PetscFunctionBegin;
5601c531cf8SMatthew 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);
5618c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
56247923291SMatthew G. Knepley }
56347923291SMatthew G. Knepley 
5641c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
5658c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
5668c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
5678c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
568191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
5698c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
5708c6c5593SMatthew G. Knepley {
5718c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
57247923291SMatthew G. Knepley 
5738c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5741c531cf8SMatthew 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);
57547923291SMatthew G. Knepley   PetscFunctionReturn(0);
57647923291SMatthew G. Knepley }
577