xref: /petsc/src/dm/impls/plex/plexproject.c (revision e39fcb4e08cf2f5d08228d3d9a1b955ec2b97d02)
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;
10c330f8ffSToby Isaac   PetscBool      isAffine;
118c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
128c6c5593SMatthew G. Knepley 
138c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
14c330f8ffSToby 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 */
19c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
20c330f8ffSToby 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]) {
26c330f8ffSToby Isaac       if (isFE[f]) {
27c330f8ffSToby Isaac         PetscQuadrature   allPoints;
28c330f8ffSToby Isaac         PetscInt          q, dim, numPoints;
29c330f8ffSToby Isaac         const PetscReal   *points;
30c330f8ffSToby Isaac         PetscScalar       *pointEval;
31c330f8ffSToby Isaac         PetscReal         *x;
32c330f8ffSToby Isaac         DM                dm;
33c330f8ffSToby Isaac 
34c330f8ffSToby Isaac         ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
35c330f8ffSToby Isaac         ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
36c330f8ffSToby Isaac         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
37c330f8ffSToby Isaac         ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
38c330f8ffSToby Isaac         ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
39c330f8ffSToby Isaac         ierr = DMGetWorkArray(dm,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           }
49c330f8ffSToby Isaac           ierr = (*funcs[f])(coordDim,time,v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
50c330f8ffSToby Isaac         }
51c330f8ffSToby Isaac         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
52c330f8ffSToby Isaac         ierr = DMRestoreWorkArray(dm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
53c330f8ffSToby Isaac         ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
54c330f8ffSToby Isaac         v += spDim;
55c330f8ffSToby Isaac       } else {
56c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
57c330f8ffSToby Isaac           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
58c330f8ffSToby Isaac         }
59c330f8ffSToby Isaac       }
60c330f8ffSToby Isaac     } else {
61c330f8ffSToby 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 
67c330f8ffSToby 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;
82c330f8ffSToby Isaac   const PetscInt     dE = fegeom->dimEmbed;
83c330f8ffSToby Isaac   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
84c330f8ffSToby 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);
96e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
97496733e2SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
988c6c5593SMatthew G. Knepley   if (dmAux) {
99dcc01f9aSMatthew G. Knepley     DMLabel  spmap;
100dcc01f9aSMatthew G. Knepley     PetscInt subp = p;
1011ba34526SMatthew G. Knepley 
102dcc01f9aSMatthew G. Knepley     /* If dm is a submesh, do not get subpoint */
103dcc01f9aSMatthew G. Knepley     ierr = DMPlexGetSubpointMap(dm, &spmap);CHKERRQ(ierr);
104*e39fcb4eSStefano Zampini     if (!spmap) {
105*e39fcb4eSStefano Zampini       PetscInt h;
106*e39fcb4eSStefano Zampini       ierr = DMPlexGetVTKCellHeight(dmAux, &h);CHKERRQ(ierr);
107*e39fcb4eSStefano Zampini       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
108*e39fcb4eSStefano Zampini       if (spmap && !h) { ierr = DMLabelGetValue(spmap, p, &subp);CHKERRQ(ierr); }
109*e39fcb4eSStefano Zampini       else             { ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr); }
110*e39fcb4eSStefano Zampini     }
1111ba34526SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
1128c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
113496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
114496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
115e87a4003SBarry Smith     ierr = DMGetSection(dmAux, &sectionAux);CHKERRQ(ierr);
1168c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
1178c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
11811d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
119496733e2SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
1201ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
1218c6c5593SMatthew G. Knepley   }
1228c6c5593SMatthew G. Knepley   /* Get values for closure */
123c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
1248c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
125c330f8ffSToby Isaac     PetscQuadrature   allPoints;
126c330f8ffSToby Isaac     PetscInt          q, dim, numPoints;
127c330f8ffSToby Isaac     const PetscReal   *points;
128c330f8ffSToby Isaac     PetscScalar       *pointEval;
129c330f8ffSToby Isaac     DM                dm;
130c330f8ffSToby Isaac 
1318c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
1328c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
133c330f8ffSToby Isaac     if (!funcs[f]) {
134be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
135c330f8ffSToby Isaac       continue;
136c330f8ffSToby Isaac     }
137c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
138c330f8ffSToby Isaac     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
139c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
140c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
1410c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
142c330f8ffSToby Isaac       const PetscReal *v0;
143c330f8ffSToby Isaac       const PetscReal *invJ;
144c330f8ffSToby Isaac 
145c330f8ffSToby Isaac       if (isAffine) {
146be1504a2SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
147c330f8ffSToby Isaac         v0 = x;
148c330f8ffSToby Isaac         invJ = fegeom->invJ;
1491c531cf8SMatthew G. Knepley       } else {
150c330f8ffSToby Isaac         v0 = &fegeom->v[tp*dE];
151c330f8ffSToby Isaac         invJ = &fegeom->invJ[tp*dE*dE];
1528c6c5593SMatthew G. Knepley       }
153c330f8ffSToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, invJ, coefficients, coefficients_t, u, u_x, u_t);
154*e39fcb4eSStefano Zampini       if (probAux) EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);
155be1504a2SMatthew G. Knepley       (*funcs[f])(dE, 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]);
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 
166c330f8ffSToby 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[],
1671c531cf8SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
1681c531cf8SMatthew G. Knepley                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
1698c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
1708c6c5593SMatthew G. Knepley {
1718c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
1728c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
1738c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
1748c6c5593SMatthew G. Knepley 
1758c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1768c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1778c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1788c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
1798c6c5593SMatthew G. Knepley   switch (type) {
1808c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
1818c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
182c330f8ffSToby 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;
1838c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
1848c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
185c330f8ffSToby Isaac     ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
1860c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
1870c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
1880c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1890c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
190191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
1918c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
1928c6c5593SMatthew G. Knepley   }
1938c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1948c6c5593SMatthew G. Knepley }
1958c6c5593SMatthew G. Knepley 
196df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
197df5c1128SToby Isaac {
198df5c1128SToby Isaac   PetscReal      *points;
199df5c1128SToby Isaac   PetscInt       f, numPoints;
200df5c1128SToby Isaac   PetscErrorCode ierr;
201df5c1128SToby Isaac 
202df5c1128SToby Isaac   PetscFunctionBegin;
203df5c1128SToby Isaac   numPoints = 0;
204df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
205df5c1128SToby Isaac     if (funcs[f]) {
206df5c1128SToby Isaac       PetscQuadrature fAllPoints;
207df5c1128SToby Isaac       PetscInt        fNumPoints;
208df5c1128SToby Isaac 
209df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
210df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
211df5c1128SToby Isaac       numPoints += fNumPoints;
212df5c1128SToby Isaac     }
213df5c1128SToby Isaac   }
214df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
215df5c1128SToby Isaac   numPoints = 0;
216df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
217df5c1128SToby Isaac     PetscInt spDim;
218df5c1128SToby Isaac 
219df5c1128SToby Isaac     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
220df5c1128SToby Isaac     if (funcs[f]) {
221df5c1128SToby Isaac       PetscQuadrature fAllPoints;
22254ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
223df5c1128SToby Isaac       const PetscReal *fPoints;
224df5c1128SToby Isaac 
225df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
22654ee0cceSMatthew G. Knepley       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
22754ee0cceSMatthew 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);
228df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
229df5c1128SToby Isaac       numPoints += fNumPoints;
230df5c1128SToby Isaac     }
231df5c1128SToby Isaac   }
232df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
233df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
234df5c1128SToby Isaac   PetscFunctionReturn(0);
235df5c1128SToby Isaac }
236df5c1128SToby Isaac 
2370de51b56SMatthew G. Knepley /*
2380de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
2390de51b56SMatthew G. Knepley 
2400de51b56SMatthew G. Knepley   There are several different scenarios:
2410de51b56SMatthew G. Knepley 
2420de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
2430de51b56SMatthew G. Knepley 
2440de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
2450de51b56SMatthew G. Knepley 
2460de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
2470de51b56SMatthew G. Knepley 
2480de51b56SMatthew 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.
2490de51b56SMatthew G. Knepley 
2500de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
2510de51b56SMatthew G. Knepley 
2520de51b56SMatthew 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.
2530de51b56SMatthew G. Knepley 
2540de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
2550de51b56SMatthew G. Knepley 
2560de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
2570de51b56SMatthew G. Knepley 
2580de51b56SMatthew 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.
2590de51b56SMatthew 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
2600de51b56SMatthew 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
2610de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
2620de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
2630de51b56SMatthew G. Knepley 
2640de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
2650de51b56SMatthew G. Knepley 
2660de51b56SMatthew 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.
2670de51b56SMatthew G. Knepley */
26846fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
2691c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
2708c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
2718c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
2728c6c5593SMatthew G. Knepley {
2738c6c5593SMatthew G. Knepley   DM              dmAux = NULL;
2740c898c18SMatthew G. Knepley   PetscDS         prob, probAux = NULL;
2758c6c5593SMatthew G. Knepley   Vec             localA = NULL;
27647923291SMatthew G. Knepley   PetscSection    section;
2778c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
2780c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
2798c6c5593SMatthew G. Knepley   PetscInt       *Nc;
280c330f8ffSToby Isaac   PetscInt        dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f;
28188aca1feSMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE;
282c330f8ffSToby Isaac   DMField         coordField;
283c330f8ffSToby Isaac   DMLabel         depthLabel;
284c330f8ffSToby Isaac   PetscQuadrature allPoints = NULL;
28547923291SMatthew G. Knepley   PetscErrorCode  ierr;
28647923291SMatthew G. Knepley 
28747923291SMatthew G. Knepley   PetscFunctionBegin;
28888aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
28988aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
2908c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2912af58022SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
2920de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
2930de51b56SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
294*e39fcb4eSStefano Zampini     if (dmAux && !localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
29588aca1feSMatthew G. Knepley     if (!minHeight && dmAux) {
29688aca1feSMatthew G. Knepley       DMLabel spmap;
29788aca1feSMatthew G. Knepley 
29888aca1feSMatthew G. Knepley       /* If dmAux is a surface, then force the projection to take place over a surface */
29988aca1feSMatthew G. Knepley       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
300*e39fcb4eSStefano Zampini       if (spmap) {
301*e39fcb4eSStefano Zampini         ierr = DMPlexGetVTKCellHeight(dmAux, &minHeight);CHKERRQ(ierr);
302*e39fcb4eSStefano Zampini         auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
303*e39fcb4eSStefano Zampini       }
30488aca1feSMatthew G. Knepley     }
3050de51b56SMatthew G. Knepley   }
3068c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
3072af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
308*e39fcb4eSStefano Zampini   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
3090c898c18SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
3108c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
311e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
31247923291SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3130c898c18SMatthew G. Knepley   if (dmAux) {
3140c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
3150c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
3160c898c18SMatthew G. Knepley   }
317496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
318496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
3198c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
3208c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
3218c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
3228c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
32347923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
32447923291SMatthew G. Knepley     PetscObject  obj;
32547923291SMatthew G. Knepley     PetscClassId id;
32647923291SMatthew G. Knepley 
32747923291SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
32847923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
32947923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
33047923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
33147923291SMatthew G. Knepley 
33247923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
33347923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
33447923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
33547923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
33647923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
33747923291SMatthew G. Knepley 
33847923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
33947923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
34047923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
34147923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
34247923291SMatthew G. Knepley   }
343c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
3440c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
34574fc349aSMatthew G. Knepley     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
34674fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
3474a3e9fdbSToby Isaac     const PetscReal *points;
3484a3e9fdbSToby Isaac     PetscInt         numPoints;
3490c898c18SMatthew G. Knepley 
3502af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
35154ee0cceSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
35254ee0cceSMatthew G. Knepley       if (!effectiveHeight) {sp[f] = cellsp[f];}
35354ee0cceSMatthew G. Knepley       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
35454ee0cceSMatthew G. Knepley     }
35554ee0cceSMatthew G. Knepley     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
356c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
3570c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
3580c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
3590c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
3600c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
36174fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
36274fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
36374fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
3640c898c18SMatthew G. Knepley     }
3650c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
3660c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
36774fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
36874fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
36974fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
3700c898c18SMatthew G. Knepley     }
3710c898c18SMatthew G. Knepley   }
372c330f8ffSToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
373c330f8ffSToby Isaac   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
37447923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
3752af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
37688aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
3770de51b56SMatthew G. Knepley     PetscDS      probEff         = prob;
3788c6c5593SMatthew G. Knepley     PetscScalar *values;
379b7260050SToby Isaac     PetscBool   *fieldActive;
380b7260050SToby Isaac     PetscInt     maxDegree;
38184f7fa7aSMatthew G. Knepley     PetscInt     pStart, pEnd, p, spDim, totDim, numValues;
382c330f8ffSToby Isaac     IS           heightIS;
3838c6c5593SMatthew G. Knepley 
3840de51b56SMatthew G. Knepley     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);}
38547923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
386c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel,depth - h,&heightIS);CHKERRQ(ierr);
3878c6c5593SMatthew G. Knepley     if (!h) {
3888c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
3898c6c5593SMatthew G. Knepley 
3908c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3918c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
3928c6c5593SMatthew G. Knepley     }
393c330f8ffSToby Isaac     if (pEnd <= pStart) {
394c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
395c330f8ffSToby Isaac       continue;
396c330f8ffSToby Isaac     }
39747923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
39847923291SMatthew G. Knepley     totDim = 0;
39947923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
400bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
40147923291SMatthew G. Knepley         sp[f] = cellsp[f];
40247923291SMatthew G. Knepley       } else {
403bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
40447923291SMatthew G. Knepley         if (!sp[f]) continue;
40547923291SMatthew G. Knepley       }
40647923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
4079c3cf19fSMatthew G. Knepley       totDim += spDim;
40847923291SMatthew G. Knepley     }
40947923291SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
4108c6c5593SMatthew 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);
411c330f8ffSToby Isaac     if (!totDim) {
412c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
413c330f8ffSToby Isaac       continue;
414c330f8ffSToby Isaac     }
41547923291SMatthew G. Knepley     /* Loop over points at this height */
41669291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
41769291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
4188c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
4198c6c5593SMatthew G. Knepley     if (label) {
4208c6c5593SMatthew G. Knepley       PetscInt i;
42147923291SMatthew G. Knepley 
42247923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
423c330f8ffSToby Isaac         IS              pointIS, isectIS;
42447923291SMatthew G. Knepley         const PetscInt *points;
4258c6c5593SMatthew G. Knepley         PetscInt        n;
426c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
427c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
42847923291SMatthew G. Knepley 
42947923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
43047923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
431c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
432c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
433c330f8ffSToby Isaac         if (!isectIS) continue;
434c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
435c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
436b7260050SToby Isaac         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
437b7260050SToby Isaac         if (maxDegree <= 1) {
438c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
439c330f8ffSToby Isaac         }
440c330f8ffSToby Isaac         if (!quad) {
441c330f8ffSToby Isaac           if (!h && allPoints) {
442c330f8ffSToby Isaac             quad = allPoints;
443c330f8ffSToby Isaac             allPoints = NULL;
444c330f8ffSToby Isaac           } else {
445c330f8ffSToby Isaac             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
446c330f8ffSToby Isaac           }
447c330f8ffSToby Isaac         }
448c330f8ffSToby Isaac         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
44947923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
45047923291SMatthew G. Knepley           const PetscInt  point = points[p];
45147923291SMatthew G. Knepley 
452b420b6ccSMatthew G. Knepley           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
453c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
454c330f8ffSToby 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);
45547923291SMatthew G. Knepley           if (ierr) {
45647923291SMatthew G. Knepley             PetscErrorCode ierr2;
45769291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
45869291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
45947923291SMatthew G. Knepley             CHKERRQ(ierr);
46047923291SMatthew G. Knepley           }
461ba322698SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
46247923291SMatthew G. Knepley         }
463c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
464c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
465c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
466c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
467c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
46847923291SMatthew G. Knepley       }
4698c6c5593SMatthew G. Knepley     } else {
470c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
471c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
472c330f8ffSToby Isaac       IS              pointIS;
473c330f8ffSToby Isaac 
474c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
475b7260050SToby Isaac       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
476b7260050SToby Isaac       if (maxDegree <= 1) {
477c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
478c330f8ffSToby Isaac       }
479c330f8ffSToby Isaac       if (!quad) {
480c330f8ffSToby Isaac         if (!h && allPoints) {
481c330f8ffSToby Isaac           quad = allPoints;
482c330f8ffSToby Isaac           allPoints = NULL;
483c330f8ffSToby Isaac         } else {
484c330f8ffSToby Isaac           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
485c330f8ffSToby Isaac         }
486c330f8ffSToby Isaac       }
487c330f8ffSToby Isaac       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
4888c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
489b420b6ccSMatthew G. Knepley         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
490c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
491c330f8ffSToby 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);
4928c6c5593SMatthew G. Knepley         if (ierr) {
4938c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
49469291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
49569291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
4968c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
4978c6c5593SMatthew G. Knepley         }
498ba322698SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
4998c6c5593SMatthew G. Knepley       }
500c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
501c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
502c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
503c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
5048c6c5593SMatthew G. Knepley     }
505c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
50669291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
50769291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
50847923291SMatthew G. Knepley   }
5098c6c5593SMatthew G. Knepley   /* Cleanup */
5100c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
51174fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
51274fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
5130c898c18SMatthew G. Knepley 
5140c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
5150c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
5160c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
51774fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
51874fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
51974fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
5200c898c18SMatthew G. Knepley     }
5210c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
5220c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
52374fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
52474fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
52574fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
5260c898c18SMatthew G. Knepley     }
5270c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
5280c898c18SMatthew G. Knepley   }
529c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
530496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
5318c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
5328c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
53347923291SMatthew G. Knepley }
5348c6c5593SMatthew G. Knepley 
5358c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
5368c6c5593SMatthew G. Knepley {
5378c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
5388c6c5593SMatthew G. Knepley 
5398c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5401c531cf8SMatthew 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);
5418c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
5428c6c5593SMatthew G. Knepley }
5438c6c5593SMatthew G. Knepley 
5441c531cf8SMatthew 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)
5458c6c5593SMatthew G. Knepley {
5468c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
5478c6c5593SMatthew G. Knepley 
5488c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5491c531cf8SMatthew 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);
55047923291SMatthew G. Knepley   PetscFunctionReturn(0);
55147923291SMatthew G. Knepley }
55247923291SMatthew G. Knepley 
5538c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
55447923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
55547923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
55647923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
557191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
55847923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
55947923291SMatthew G. Knepley {
56047923291SMatthew G. Knepley   PetscErrorCode ierr;
56147923291SMatthew G. Knepley 
56247923291SMatthew G. Knepley   PetscFunctionBegin;
5631c531cf8SMatthew 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);
5648c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
56547923291SMatthew G. Knepley }
56647923291SMatthew G. Knepley 
5671c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
5688c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
5698c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
5708c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
571191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
5728c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
5738c6c5593SMatthew G. Knepley {
5748c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
57547923291SMatthew G. Knepley 
5768c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5771c531cf8SMatthew 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);
57847923291SMatthew G. Knepley   PetscFunctionReturn(0);
57947923291SMatthew G. Knepley }
580