xref: /petsc/src/dm/impls/plex/plexproject.c (revision be1504a2a14612606f10c91ce91b25d1ba71133c)
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);
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     PetscInt subp;
1001ba34526SMatthew G. Knepley 
101559a1558SMatthew G. Knepley     ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr);
1021ba34526SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
1038c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
104496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
105496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
1068c6c5593SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
1078c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
1088c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
10911d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
110496733e2SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
1111ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
1128c6c5593SMatthew G. Knepley   }
1138c6c5593SMatthew G. Knepley   /* Get values for closure */
114c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
1158c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
116c330f8ffSToby Isaac     PetscQuadrature   allPoints;
117c330f8ffSToby Isaac     PetscInt          q, dim, numPoints;
118c330f8ffSToby Isaac     const PetscReal   *points;
119c330f8ffSToby Isaac     PetscScalar       *pointEval;
120c330f8ffSToby Isaac     DM                dm;
121c330f8ffSToby Isaac 
1228c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
1238c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
124c330f8ffSToby Isaac     if (!funcs[f]) {
125*be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
126c330f8ffSToby Isaac       continue;
127c330f8ffSToby Isaac     }
128c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
129c330f8ffSToby Isaac     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
130c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
131c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
1320c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
133c330f8ffSToby Isaac       const PetscReal *v0;
134c330f8ffSToby Isaac       const PetscReal *invJ;
135c330f8ffSToby Isaac 
136c330f8ffSToby Isaac       if (isAffine) {
137*be1504a2SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
138c330f8ffSToby Isaac         v0 = x;
139c330f8ffSToby Isaac         invJ = fegeom->invJ;
1401c531cf8SMatthew G. Knepley       } else {
141c330f8ffSToby Isaac         v0 = &fegeom->v[tp*dE];
142c330f8ffSToby Isaac         invJ = &fegeom->invJ[tp*dE*dE];
1438c6c5593SMatthew G. Knepley       }
144c330f8ffSToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, invJ, coefficients, coefficients_t, u, u_x, u_t);
145c330f8ffSToby Isaac       if (probAux) {EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);}
146*be1504a2SMatthew 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]);
1471c531cf8SMatthew G. Knepley     }
148c330f8ffSToby Isaac     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
149c330f8ffSToby Isaac     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
150c330f8ffSToby Isaac     v += spDim;
1518c6c5593SMatthew G. Knepley   }
1528c6c5593SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1538c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
1548c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1558c6c5593SMatthew G. Knepley }
1568c6c5593SMatthew G. Knepley 
157c330f8ffSToby 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[],
1581c531cf8SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
1591c531cf8SMatthew G. Knepley                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
1608c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
1618c6c5593SMatthew G. Knepley {
1628c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
1638c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
1648c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
1658c6c5593SMatthew G. Knepley 
1668c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1678c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1688c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1698c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
1708c6c5593SMatthew G. Knepley   switch (type) {
1718c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
1728c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
173c330f8ffSToby 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;
1748c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
1758c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
176c330f8ffSToby Isaac     ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
1770c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
1780c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
1790c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1800c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
181191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
1828c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
1838c6c5593SMatthew G. Knepley   }
1848c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1858c6c5593SMatthew G. Knepley }
1868c6c5593SMatthew G. Knepley 
187df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
188df5c1128SToby Isaac {
189df5c1128SToby Isaac   PetscReal      *points;
190df5c1128SToby Isaac   PetscInt       f, numPoints;
191df5c1128SToby Isaac   PetscErrorCode ierr;
192df5c1128SToby Isaac 
193df5c1128SToby Isaac   PetscFunctionBegin;
194df5c1128SToby Isaac   numPoints = 0;
195df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
196df5c1128SToby Isaac     if (funcs[f]) {
197df5c1128SToby Isaac       PetscQuadrature fAllPoints;
198df5c1128SToby Isaac       PetscInt        fNumPoints;
199df5c1128SToby Isaac 
200df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
201df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
202df5c1128SToby Isaac       numPoints += fNumPoints;
203df5c1128SToby Isaac     }
204df5c1128SToby Isaac   }
205df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
206df5c1128SToby Isaac   numPoints = 0;
207df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
208df5c1128SToby Isaac     PetscInt spDim;
209df5c1128SToby Isaac 
210df5c1128SToby Isaac     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
211df5c1128SToby Isaac     if (funcs[f]) {
212df5c1128SToby Isaac       PetscQuadrature fAllPoints;
213df5c1128SToby Isaac       PetscInt        fNumPoints, q;
214df5c1128SToby Isaac       const PetscReal *fPoints;
215df5c1128SToby Isaac 
216df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
217df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
218df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
219df5c1128SToby Isaac       numPoints += fNumPoints;
220df5c1128SToby Isaac     }
221df5c1128SToby Isaac   }
222df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
223df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
224df5c1128SToby Isaac   PetscFunctionReturn(0);
225df5c1128SToby Isaac }
226df5c1128SToby Isaac 
2270de51b56SMatthew G. Knepley /*
2280de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
2290de51b56SMatthew G. Knepley 
2300de51b56SMatthew G. Knepley   There are several different scenarios:
2310de51b56SMatthew G. Knepley 
2320de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
2330de51b56SMatthew G. Knepley 
2340de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
2350de51b56SMatthew G. Knepley 
2360de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
2370de51b56SMatthew G. Knepley 
2380de51b56SMatthew 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.
2390de51b56SMatthew G. Knepley 
2400de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
2410de51b56SMatthew G. Knepley 
2420de51b56SMatthew 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.
2430de51b56SMatthew G. Knepley 
2440de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
2450de51b56SMatthew G. Knepley 
2460de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
2470de51b56SMatthew G. Knepley 
2480de51b56SMatthew 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.
2490de51b56SMatthew 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
2500de51b56SMatthew 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
2510de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
2520de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
2530de51b56SMatthew G. Knepley 
2540de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
2550de51b56SMatthew G. Knepley 
2560de51b56SMatthew 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.
2570de51b56SMatthew G. Knepley */
25846fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
2591c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
2608c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
2618c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
2628c6c5593SMatthew G. Knepley {
2638c6c5593SMatthew G. Knepley   DM              dmAux = NULL;
2640c898c18SMatthew G. Knepley   PetscDS         prob, probAux = NULL;
2658c6c5593SMatthew G. Knepley   Vec             localA = NULL;
26647923291SMatthew G. Knepley   PetscSection    section;
2678c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
2680c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
2698c6c5593SMatthew G. Knepley   PetscInt       *Nc;
270c330f8ffSToby Isaac   PetscInt        dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f;
27188aca1feSMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE;
272c330f8ffSToby Isaac   DMField         coordField;
273c330f8ffSToby Isaac   DMLabel         depthLabel;
274c330f8ffSToby Isaac   PetscQuadrature allPoints = NULL;
27547923291SMatthew G. Knepley   PetscErrorCode  ierr;
27647923291SMatthew G. Knepley 
27747923291SMatthew G. Knepley   PetscFunctionBegin;
27888aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
27988aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
2808c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2812af58022SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
2820de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
2830de51b56SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
28488aca1feSMatthew G. Knepley     if (!minHeight && dmAux) {
28588aca1feSMatthew G. Knepley       DMLabel spmap;
28688aca1feSMatthew G. Knepley 
28788aca1feSMatthew G. Knepley       /* If dmAux is a surface, then force the projection to take place over a surface */
28888aca1feSMatthew G. Knepley       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
28988aca1feSMatthew G. Knepley       if (spmap) {minHeight = 1; auxBd = PETSC_TRUE;}
29088aca1feSMatthew G. Knepley     }
2910de51b56SMatthew G. Knepley   }
2928c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
2932af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
2948c6c5593SMatthew 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);}
2950c898c18SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
2968c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
29747923291SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
29847923291SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
2990c898c18SMatthew G. Knepley   if (dmAux) {
3000c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
3010c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
3020c898c18SMatthew G. Knepley   }
303496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
304496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
3058c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
3068c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
3078c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
3088c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
30947923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
31047923291SMatthew G. Knepley     PetscObject  obj;
31147923291SMatthew G. Knepley     PetscClassId id;
31247923291SMatthew G. Knepley 
31347923291SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
31447923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
31547923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
31647923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
31747923291SMatthew G. Knepley 
31847923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
31947923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
32047923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
32147923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
32247923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
32347923291SMatthew G. Knepley 
32447923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
32547923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
32647923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
32747923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
32847923291SMatthew G. Knepley   }
329c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
3300c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
33174fc349aSMatthew G. Knepley     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
33274fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
3334a3e9fdbSToby Isaac     const PetscReal *points;
3344a3e9fdbSToby Isaac     PetscInt         numPoints;
3350c898c18SMatthew G. Knepley 
3362af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
337c330f8ffSToby Isaac     ierr = PetscDualSpaceGetAllPointsUnion(Nf,cellsp,dim,funcs,&allPoints);CHKERRQ(ierr);
338c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
3390c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
3400c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
3410c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
3420c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
34374fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
34474fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
34574fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
3460c898c18SMatthew G. Knepley     }
3470c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
3480c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
34974fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
35074fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
35174fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
3520c898c18SMatthew G. Knepley     }
3530c898c18SMatthew G. Knepley   }
354c330f8ffSToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
355c330f8ffSToby Isaac   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
35647923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
3572af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
35888aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
3590de51b56SMatthew G. Knepley     PetscDS      probEff         = prob;
3608c6c5593SMatthew G. Knepley     PetscScalar *values;
361c330f8ffSToby Isaac     PetscBool   *fieldActive, isAffine;
36284f7fa7aSMatthew G. Knepley     PetscInt     pStart, pEnd, p, spDim, totDim, numValues;
363c330f8ffSToby Isaac     IS           heightIS;
3648c6c5593SMatthew G. Knepley 
3650de51b56SMatthew G. Knepley     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);}
36647923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
367c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel,depth - h,&heightIS);CHKERRQ(ierr);
3688c6c5593SMatthew G. Knepley     if (!h) {
3698c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
3708c6c5593SMatthew G. Knepley 
3718c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
3728c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
3738c6c5593SMatthew G. Knepley     }
374c330f8ffSToby Isaac     if (pEnd <= pStart) {
375c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
376c330f8ffSToby Isaac       continue;
377c330f8ffSToby Isaac     }
37847923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
37947923291SMatthew G. Knepley     totDim = 0;
38047923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
381bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
38247923291SMatthew G. Knepley         sp[f] = cellsp[f];
38347923291SMatthew G. Knepley       } else {
384bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
38547923291SMatthew G. Knepley         if (!sp[f]) continue;
38647923291SMatthew G. Knepley       }
38747923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
3889c3cf19fSMatthew G. Knepley       totDim += spDim;
38947923291SMatthew G. Knepley     }
39047923291SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
3918c6c5593SMatthew 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);
392c330f8ffSToby Isaac     if (!totDim) {
393c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
394c330f8ffSToby Isaac       continue;
395c330f8ffSToby Isaac     }
39647923291SMatthew G. Knepley     /* Loop over points at this height */
39769291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
39869291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
3998c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
4008c6c5593SMatthew G. Knepley     if (label) {
4018c6c5593SMatthew G. Knepley       PetscInt i;
40247923291SMatthew G. Knepley 
40347923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
404c330f8ffSToby Isaac         IS              pointIS, isectIS;
40547923291SMatthew G. Knepley         const PetscInt *points;
4068c6c5593SMatthew G. Knepley         PetscInt        n;
407c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
408c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
40947923291SMatthew G. Knepley 
41047923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
41147923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
412c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
413c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
414c330f8ffSToby Isaac         if (!isectIS) continue;
415c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
416c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
417c330f8ffSToby Isaac         ierr = DMFieldGetFEInvariance(coordField,isectIS,NULL,&isAffine,NULL);CHKERRQ(ierr);
418c330f8ffSToby Isaac         if (isAffine) {
419c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
420c330f8ffSToby Isaac         }
421c330f8ffSToby Isaac         if (!quad) {
422c330f8ffSToby Isaac           if (!h && allPoints) {
423c330f8ffSToby Isaac             quad = allPoints;
424c330f8ffSToby Isaac             allPoints = NULL;
425c330f8ffSToby Isaac           } else {
426c330f8ffSToby Isaac             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
427c330f8ffSToby Isaac           }
428c330f8ffSToby Isaac         }
429c330f8ffSToby Isaac         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
43047923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
43147923291SMatthew G. Knepley           const PetscInt  point = points[p];
43247923291SMatthew G. Knepley 
433b420b6ccSMatthew G. Knepley           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
434c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
435c330f8ffSToby 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);
43647923291SMatthew G. Knepley           if (ierr) {
43747923291SMatthew G. Knepley             PetscErrorCode ierr2;
43869291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
43969291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
44047923291SMatthew G. Knepley             CHKERRQ(ierr);
44147923291SMatthew G. Knepley           }
442ba322698SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
44347923291SMatthew G. Knepley         }
444c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
445c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
446c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
447c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
448c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
44947923291SMatthew G. Knepley       }
4508c6c5593SMatthew G. Knepley     } else {
451c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
452c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
453c330f8ffSToby Isaac       IS              pointIS;
454c330f8ffSToby Isaac 
455c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
456c330f8ffSToby Isaac       ierr = DMFieldGetFEInvariance(coordField,pointIS,NULL,&isAffine,NULL);CHKERRQ(ierr);
457c330f8ffSToby Isaac       if (isAffine) {
458c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
459c330f8ffSToby Isaac       }
460c330f8ffSToby Isaac       if (!quad) {
461c330f8ffSToby Isaac         if (!h && allPoints) {
462c330f8ffSToby Isaac           quad = allPoints;
463c330f8ffSToby Isaac           allPoints = NULL;
464c330f8ffSToby Isaac         } else {
465c330f8ffSToby Isaac           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
466c330f8ffSToby Isaac         }
467c330f8ffSToby Isaac       }
468c330f8ffSToby Isaac       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
4698c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
470b420b6ccSMatthew G. Knepley         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
471c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
472c330f8ffSToby 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);
4738c6c5593SMatthew G. Knepley         if (ierr) {
4748c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
47569291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
47669291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
4778c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
4788c6c5593SMatthew G. Knepley         }
479ba322698SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
4808c6c5593SMatthew G. Knepley       }
481c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
482c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
483c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
484c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
4858c6c5593SMatthew G. Knepley     }
486c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
48769291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
48869291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
48947923291SMatthew G. Knepley   }
4908c6c5593SMatthew G. Knepley   /* Cleanup */
4910c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
49274fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
49374fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
4940c898c18SMatthew G. Knepley 
4950c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
4960c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
4970c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
49874fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
49974fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
50074fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
5010c898c18SMatthew G. Knepley     }
5020c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
5030c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
50474fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
50574fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
50674fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
5070c898c18SMatthew G. Knepley     }
5080c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
5090c898c18SMatthew G. Knepley   }
510c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
511496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
5128c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
5138c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
51447923291SMatthew G. Knepley }
5158c6c5593SMatthew G. Knepley 
5168c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
5178c6c5593SMatthew G. Knepley {
5188c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
5198c6c5593SMatthew G. Knepley 
5208c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5211c531cf8SMatthew 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);
5228c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
5238c6c5593SMatthew G. Knepley }
5248c6c5593SMatthew G. Knepley 
5251c531cf8SMatthew 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)
5268c6c5593SMatthew G. Knepley {
5278c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
5288c6c5593SMatthew G. Knepley 
5298c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5301c531cf8SMatthew 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);
53147923291SMatthew G. Knepley   PetscFunctionReturn(0);
53247923291SMatthew G. Knepley }
53347923291SMatthew G. Knepley 
5348c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
53547923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
53647923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
53747923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
538191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
53947923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
54047923291SMatthew G. Knepley {
54147923291SMatthew G. Knepley   PetscErrorCode ierr;
54247923291SMatthew G. Knepley 
54347923291SMatthew G. Knepley   PetscFunctionBegin;
5441c531cf8SMatthew 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);
5458c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
54647923291SMatthew G. Knepley }
54747923291SMatthew G. Knepley 
5481c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
5498c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
5508c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
5518c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
552191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
5538c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
5548c6c5593SMatthew G. Knepley {
5558c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
55647923291SMatthew G. Knepley 
5578c6c5593SMatthew G. Knepley   PetscFunctionBegin;
5581c531cf8SMatthew 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);
55947923291SMatthew G. Knepley   PetscFunctionReturn(0);
56047923291SMatthew G. Knepley }
561