xref: /petsc/src/dm/impls/plex/plexproject.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #include <petsc/private/petscfeimpl.h>
4 
5 /*
6   DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point
7 
8   Input Parameters:
9 + dm     - The output DM
10 . ds     - The output DS
11 . dmIn   - The input DM
12 . dsIn   - The input DS
13 . time   - The time for this evaluation
14 . fegeom - The FE geometry for this point
15 . fvgeom - The FV geometry for this point
16 . isFE   - Flag indicating whether each output field has an FE discretization
17 . sp     - The output PetscDualSpace for each field
18 . funcs  - The evaluation function for each field
19 - ctxs   - The user context for each field
20 
21   Output Parameter:
22 . values - The value for each dual basis vector in the output dual space
23 
24   Level: developer
25 
26 .seealso: DMProjectPoint_Field_Private()
27 */
28 static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
29                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
30                                                   PetscScalar values[])
31 {
32   PetscInt       coordDim, Nf, *Nc, f, spDim, d, v, tp;
33   PetscBool      isAffine, transform;
34   PetscErrorCode ierr;
35 
36   PetscFunctionBeginHot;
37   ierr = DMGetCoordinateDim(dmIn, &coordDim);CHKERRQ(ierr);
38   ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr);
39   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
40   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
41   /* Get values for closure */
42   isAffine = fegeom->isAffine;
43   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
44     void * const ctx = ctxs ? ctxs[f] : NULL;
45 
46     if (!sp[f]) continue;
47     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
48     if (funcs[f]) {
49       if (isFE[f]) {
50         PetscQuadrature   allPoints;
51         PetscInt          q, dim, numPoints;
52         const PetscReal   *points;
53         PetscScalar       *pointEval;
54         PetscReal         *x;
55         DM                rdm;
56 
57         ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr);
58         ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
59         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
60         ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
61         ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
62         for (q = 0; q < numPoints; q++, tp++) {
63           const PetscReal *v0;
64 
65           if (isAffine) {
66             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
67             v0 = x;
68           } else {
69             v0 = &fegeom->v[tp*coordDim];
70           }
71           if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;}
72           ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
73         }
74         /* Transform point evaluations pointEval[q,c] */
75         ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr);
76         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
77         ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
78         ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
79         v += spDim;
80       } else {
81         for (d = 0; d < spDim; ++d, ++v) {
82           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
83         }
84       }
85     } else {
86       for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;}
87     }
88   }
89   PetscFunctionReturn(0);
90 }
91 
92 /*
93   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
94 
95   Input Parameters:
96 + dm             - The output DM
97 . ds             - The output DS
98 . dmIn           - The input DM
99 . dsIn           - The input DS
100 . dmAux          - The auxiliary DM, which is always for the input space
101 . dsAux          - The auxiliary DS, which is always for the input space
102 . time           - The time for this evaluation
103 . localU         - The local solution
104 . localA         - The local auziliary fields
105 . cgeom          - The FE geometry for this point
106 . sp             - The output PetscDualSpace for each field
107 . p              - The point in the output DM
108 . T              - Input basis and derviatives for each field tabulated on the quadrature points
109 . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
110 . funcs          - The evaluation function for each field
111 - ctxs           - The user context for each field
112 
113   Output Parameter:
114 . values         - The value for each dual basis vector in the output dual space
115 
116   Note: Not supported for FV
117 
118   Level: developer
119 
120 .seealso: DMProjectPoint_Field_Private()
121 */
122 static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p,
123                                                    PetscTabulation *T, PetscTabulation *TAux,
124                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
125                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
126                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
127                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
128                                                    PetscScalar values[])
129 {
130   PetscSection       section, sectionAux = NULL;
131   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
132   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
133   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
134   const PetscScalar *constants;
135   PetscReal         *x;
136   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
137   PetscFEGeom        fegeom;
138   const PetscInt     dE = cgeom->dimEmbed;
139   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
140   PetscBool          isAffine, transform;
141   PetscErrorCode     ierr;
142 
143   PetscFunctionBeginHot;
144   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
145   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
146   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
147   ierr = PetscDSGetComponentOffsets(dsIn, &uOff);CHKERRQ(ierr);
148   ierr = PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x);CHKERRQ(ierr);
149   ierr = PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
150   ierr = PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
151   ierr = PetscDSGetConstants(dsIn, &numConstants, &constants);CHKERRQ(ierr);
152   ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr);
153   ierr = DMGetLocalSection(dmIn, &section);CHKERRQ(ierr);
154   ierr = DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp);CHKERRQ(ierr);
155   ierr = DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
156   if (dmAux) {
157     PetscInt subp;
158 
159     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
160     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
161     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
162     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
163     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
164     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
165     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
166   }
167   /* Get values for closure */
168   isAffine = cgeom->isAffine;
169   if (isAffine) {
170     fegeom.v    = x;
171     fegeom.xi   = cgeom->xi;
172     fegeom.J    = cgeom->J;
173     fegeom.invJ = cgeom->invJ;
174     fegeom.detJ = cgeom->detJ;
175   }
176   for (f = 0, v = 0; f < Nf; ++f) {
177     PetscQuadrature   allPoints;
178     PetscInt          q, dim, numPoints;
179     const PetscReal   *points;
180     PetscScalar       *pointEval;
181     DM                dm;
182 
183     if (!sp[f]) continue;
184     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
185     if (!funcs[f]) {
186       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
187       continue;
188     }
189     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
190     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
191     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
192     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
193     for (q = 0; q < numPoints; ++q, ++tp) {
194       if (isAffine) {
195         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
196       } else {
197         fegeom.v    = &cgeom->v[tp*dE];
198         fegeom.J    = &cgeom->J[tp*dE*dE];
199         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
200         fegeom.detJ = &cgeom->detJ[tp];
201       }
202       ierr = PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
203       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
204       if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);}
205       (*funcs[f])(dE, NfIn, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, numConstants, constants, &pointEval[Nc[f]*q]);
206     }
207     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
208     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
209     v += spDim;
210   }
211   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
212   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
213   PetscFunctionReturn(0);
214 }
215 
216 static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p,
217                                                      PetscTabulation *T, PetscTabulation *TAux,
218                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
219                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
220                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
221                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
222                                                      PetscScalar values[])
223 {
224   PetscSection       section, sectionAux = NULL;
225   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
226   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
227   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
228   const PetscScalar *constants;
229   PetscReal         *x;
230   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
231   PetscFEGeom        fegeom, cgeom;
232   const PetscInt     dE = fgeom->dimEmbed;
233   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
234   PetscBool          isAffine;
235   PetscErrorCode     ierr;
236 
237   PetscFunctionBeginHot;
238   if (dm != dmIn) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
239   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
240   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
241   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
242   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
243   ierr = PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
244   ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
245   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
246   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
247   ierr = DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
248   if (dmAux) {
249     PetscInt subp;
250 
251     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
252     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
253     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
254     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
255     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
256     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
257     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
258   }
259   /* Get values for closure */
260   isAffine = fgeom->isAffine;
261   fegeom.n  = 0;
262   fegeom.J  = 0;
263   fegeom.v  = 0;
264   fegeom.xi = 0;
265   cgeom.dim      = fgeom->dim;
266   cgeom.dimEmbed = fgeom->dimEmbed;
267   if (isAffine) {
268     fegeom.v    = x;
269     fegeom.xi   = fgeom->xi;
270     fegeom.J    = fgeom->J;
271     fegeom.invJ = fgeom->invJ;
272     fegeom.detJ = fgeom->detJ;
273     fegeom.n    = fgeom->n;
274 
275     cgeom.J     = fgeom->suppJ[0];
276     cgeom.invJ  = fgeom->suppInvJ[0];
277     cgeom.detJ  = fgeom->suppDetJ[0];
278   }
279   for (f = 0, v = 0; f < Nf; ++f) {
280     PetscQuadrature   allPoints;
281     PetscInt          q, dim, numPoints;
282     const PetscReal   *points;
283     PetscScalar       *pointEval;
284     DM                dm;
285 
286     if (!sp[f]) continue;
287     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
288     if (!funcs[f]) {
289       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
290       continue;
291     }
292     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
293     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
294     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
295     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
296     for (q = 0; q < numPoints; ++q, ++tp) {
297       if (isAffine) {
298         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
299       } else {
300         fegeom.v    = &fgeom->v[tp*dE];
301         fegeom.J    = &fgeom->J[tp*dE*dE];
302         fegeom.invJ = &fgeom->invJ[tp*dE*dE];
303         fegeom.detJ = &fgeom->detJ[tp];
304         fegeom.n    = &fgeom->n[tp*dE];
305 
306         cgeom.J     = &fgeom->suppJ[0][tp*dE*dE];
307         cgeom.invJ  = &fgeom->suppInvJ[0][tp*dE*dE];
308         cgeom.detJ  = &fgeom->suppDetJ[0][tp];
309       }
310       /* TODO We should use cgeom here, instead of fegeom, however the geometry coming in through fgeom does not have the support cell geometry */
311       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
312       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
313       (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, fegeom.n, numConstants, constants, &pointEval[Nc[f]*q]);
314     }
315     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
316     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
317     v += spDim;
318   }
319   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
320   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
321   PetscFunctionReturn(0);
322 }
323 
324 static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscFEGeom *fegeom, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
325                                              PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux,
326                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
327 {
328   PetscFVCellGeom fvgeom;
329   PetscInt        dim, dimEmbed;
330   PetscErrorCode  ierr;
331 
332   PetscFunctionBeginHot;
333   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
334   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
335   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
336   switch (type) {
337   case DM_BC_ESSENTIAL:
338   case DM_BC_NATURAL:
339     ierr = DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
340   case DM_BC_ESSENTIAL_FIELD:
341   case DM_BC_NATURAL_FIELD:
342     ierr = DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
343                                         (void (**)(PetscInt, PetscInt, PetscInt,
344                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
345                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
346                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
347   case DM_BC_ESSENTIAL_BD_FIELD:
348     ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
349                                           (void (**)(PetscInt, PetscInt, PetscInt,
350                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
351                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
352                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
353   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
354   }
355   PetscFunctionReturn(0);
356 }
357 
358 static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
359 {
360   PetscReal      *points;
361   PetscInt       f, numPoints;
362   PetscErrorCode ierr;
363 
364   PetscFunctionBegin;
365   numPoints = 0;
366   for (f = 0; f < Nf; ++f) {
367     if (funcs[f]) {
368       PetscQuadrature fAllPoints;
369       PetscInt        fNumPoints;
370 
371       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
372       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
373       numPoints += fNumPoints;
374     }
375   }
376   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
377   numPoints = 0;
378   for (f = 0; f < Nf; ++f) {
379     PetscInt spDim;
380 
381     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
382     if (funcs[f]) {
383       PetscQuadrature fAllPoints;
384       PetscInt        qdim, fNumPoints, q;
385       const PetscReal *fPoints;
386 
387       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
388       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
389       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);
390       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
391       numPoints += fNumPoints;
392     }
393   }
394   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
395   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
396   PetscFunctionReturn(0);
397 }
398 
399 static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *ds)
400 {
401   DM              plex;
402   DMEnclosureType enc;
403   DMLabel         depthLabel;
404   PetscInt        dim, cdepth, ls = -1, i;
405   PetscErrorCode  ierr;
406 
407   PetscFunctionBegin;
408   if (lStart) *lStart = -1;
409   if (!label) PetscFunctionReturn(0);
410   ierr = DMGetEnclosureRelation(dm, odm, &enc);CHKERRQ(ierr);
411   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
412   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
413   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
414   cdepth = dim - height;
415   for (i = 0; i < numIds; ++i) {
416     IS              pointIS;
417     const PetscInt *points;
418     PetscInt        pdepth, point;
419 
420     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
421     if (!pointIS) continue; /* No points with that id on this process */
422     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
423     ierr = DMGetEnclosurePoint(dm, odm, enc, points[0], &point);CHKERRQ(ierr);
424     ierr = DMLabelGetValue(depthLabel, point, &pdepth);CHKERRQ(ierr);
425     if (pdepth == cdepth) {
426       ls = point;
427       if (ds) {ierr = DMGetCellDS(dm, ls, ds);CHKERRQ(ierr);}
428     }
429     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
430     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
431     if (ls >= 0) break;
432   }
433   if (lStart) *lStart = ls;
434   ierr = DMDestroy(&plex);CHKERRQ(ierr);
435   PetscFunctionReturn(0);
436 }
437 
438 /*
439   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
440 
441   There are several different scenarios:
442 
443   1) Volumetric mesh with volumetric auxiliary data
444 
445      Here minHeight=0 since we loop over cells.
446 
447   2) Boundary mesh with boundary auxiliary data
448 
449      Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation.
450 
451   3) Volumetric mesh with boundary auxiliary data
452 
453      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.
454 
455   The maxHeight is used to support enforcement of constraints in DMForest.
456 
457   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
458 
459   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.
460     - 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
461       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
462       dual spaces for the boundary from our input spaces.
463     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
464 
465   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
466 
467   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.
468 */
469 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
470                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
471                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
472                                                   InsertMode mode, Vec localX)
473 {
474   DM                 plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
475   DMEnclosureType    encIn, encAux;
476   PetscDS            ds = NULL, dsIn = NULL, dsAux = NULL;
477   Vec                localA = NULL, tv;
478   PetscSection       section;
479   PetscDualSpace    *sp, *cellsp;
480   PetscTabulation *T = NULL, *TAux = NULL;
481   PetscInt          *Nc;
482   PetscInt           dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfIn, NfAux = 0, f;
483   PetscBool         *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, transform;
484   DMField            coordField;
485   DMLabel            depthLabel;
486   PetscQuadrature    allPoints = NULL;
487   PetscErrorCode     ierr;
488 
489   PetscFunctionBegin;
490   if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);}
491   else        {dmIn = dm;}
492   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
493   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
494   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
495   ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr);
496   ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr);
497   ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
498   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
499   ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr);
500   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
501   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
502   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
503   /* Auxiliary information can only be used with interpolation of field functions */
504   if (dmAux) {
505     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
506     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
507       if (!localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
508       if (!minHeight) {
509         DMLabel spmap;
510 
511         /* If dmAux is a surface, then force the projection to take place over a surface */
512         ierr = DMPlexGetSubpointMap(plexAux, &spmap);CHKERRQ(ierr);
513         if (spmap) {
514           ierr = DMPlexGetVTKCellHeight(plexAux, &minHeight);CHKERRQ(ierr);
515           auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
516         }
517       }
518     }
519   }
520   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
521   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
522   ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr);
523   maxHeight = PetscMax(maxHeight, minHeight);
524   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
525   ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr);
526   if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);}
527   ierr = DMGetFirstLabelEntry_Private(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr);
528   if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);}
529   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
530   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
531   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
532   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
533   if (dmAux) {
534     ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr);
535     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
536   }
537   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
538   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
539   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
540   else               {cellsp = sp;}
541   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
542   /* Get cell dual spaces */
543   for (f = 0; f < Nf; ++f) {
544     PetscObject  obj;
545     PetscClassId id;
546 
547     ierr = PetscDSGetDiscretization(ds, f, &obj);CHKERRQ(ierr);
548     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
549     if (id == PETSCFE_CLASSID) {
550       PetscFE fe = (PetscFE) obj;
551 
552       hasFE   = PETSC_TRUE;
553       isFE[f] = PETSC_TRUE;
554       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
555     } else if (id == PETSCFV_CLASSID) {
556       PetscFV fv = (PetscFV) obj;
557 
558       hasFV   = PETSC_TRUE;
559       isFE[f] = PETSC_FALSE;
560       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
561     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
562   }
563   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
564   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
565     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
566     PetscFE          fem, subfem;
567     PetscBool        isfe;
568     const PetscReal *points;
569     PetscInt         numPoints;
570 
571     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
572     for (f = 0; f < Nf; ++f) {
573       if (!effectiveHeight) {sp[f] = cellsp[f];}
574       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
575     }
576     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
577     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
578     ierr = PetscMalloc2(NfIn, &T, NfAux, &TAux);CHKERRQ(ierr);
579     for (f = 0; f < NfIn; ++f) {
580       ierr = PetscDSIsFE_Internal(dsIn, f, &isfe);CHKERRQ(ierr);
581       if (!isfe) continue;
582       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
583       if (!effectiveHeight) {subfem = fem;}
584       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
585       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]);CHKERRQ(ierr);
586     }
587     for (f = 0; f < NfAux; ++f) {
588       ierr = PetscDSIsFE_Internal(dsAux, f, &isfe);CHKERRQ(ierr);
589       if (!isfe) continue;
590       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
591       if (!effectiveHeight || auxBd) {subfem = fem;}
592       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
593       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]);CHKERRQ(ierr);
594     }
595   }
596   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
597   for (h = minHeight; h <= maxHeight; h++) {
598     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
599     PetscDS      dsEff         = ds;
600     PetscScalar *values;
601     PetscBool   *fieldActive;
602     PetscInt     maxDegree;
603     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
604     IS           heightIS;
605 
606     /* Note we assume that dm and dmIn share the same topology */
607     ierr = DMPlexGetHeightStratum(plex, h, &pStart, &pEnd);CHKERRQ(ierr);
608     ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
609     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
610     if (!h) {
611       PetscInt cEndInterior;
612 
613       ierr = DMPlexGetInteriorCellStratum(dm, NULL, &cEndInterior);CHKERRQ(ierr);
614       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
615     }
616     if (pEnd <= pStart) {
617       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
618       continue;
619     }
620     /* Compute totDim, the number of dofs in the closure of a point at this height */
621     totDim = 0;
622     for (f = 0; f < Nf; ++f) {
623       if (!effectiveHeight) {
624         sp[f] = cellsp[f];
625       } else {
626         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
627         if (!sp[f]) continue;
628       }
629       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
630       totDim += spDim;
631     }
632     ierr = DMPlexVecGetClosure(plex, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr);
633     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim);
634     if (!totDim) {
635       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
636       continue;
637     }
638     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(ds, effectiveHeight, &dsEff);CHKERRQ(ierr);}
639     /* Loop over points at this height */
640     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
641     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
642     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
643     if (label) {
644       PetscInt i;
645 
646       for (i = 0; i < numIds; ++i) {
647         IS              pointIS, isectIS;
648         const PetscInt *points;
649         PetscInt        n;
650         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
651         PetscQuadrature quad = NULL;
652 
653         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
654         if (!pointIS) continue; /* No points with that id on this process */
655         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
656         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
657         if (!isectIS) continue;
658         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
659         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
660         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
661         if (maxDegree <= 1) {
662           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
663         }
664         if (!quad) {
665           if (!h && allPoints) {
666             quad = allPoints;
667             allPoints = NULL;
668           } else {
669             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
670           }
671         }
672         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
673         for (p = 0; p < n; ++p) {
674           const PetscInt  point = points[p];
675 
676           ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
677           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
678           ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values);
679           if (ierr) {
680             PetscErrorCode ierr2;
681             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
682             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
683             CHKERRQ(ierr);
684           }
685           if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
686           ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
687         }
688         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
689         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
690         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
691         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
692         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
693       }
694     } else {
695       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
696       PetscQuadrature quad = NULL;
697       IS              pointIS;
698 
699       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
700       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
701       if (maxDegree <= 1) {
702         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
703       }
704       if (!quad) {
705         if (!h && allPoints) {
706           quad = allPoints;
707           allPoints = NULL;
708         } else {
709           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
710         }
711       }
712       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
713       for (p = pStart; p < pEnd; ++p) {
714         ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
715         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
716         ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values);
717         if (ierr) {
718           PetscErrorCode ierr2;
719           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
720           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
721           CHKERRQ(ierr);
722         }
723         if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
724         ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
725       }
726       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
727       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
728       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
729       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
730     }
731     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
732     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
733     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
734   }
735   /* Cleanup */
736   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
737     PetscInt  effectiveHeight = auxBd ? minHeight : 0;
738     PetscFE   fem, subfem;
739     PetscBool isfe;
740 
741     for (f = 0; f < NfIn; ++f) {
742       ierr = PetscDSIsFE_Internal(dsIn, f, &isfe);CHKERRQ(ierr);
743       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
744       if (!effectiveHeight) {subfem = fem;}
745       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
746       ierr = PetscTabulationDestroy(&T[f]);CHKERRQ(ierr);
747     }
748     for (f = 0; f < NfAux; ++f) {
749       ierr = PetscDSIsFE_Internal(dsAux, f, &isfe);CHKERRQ(ierr);
750       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
751       if (!effectiveHeight || auxBd) {subfem = fem;}
752       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
753       ierr = PetscTabulationDestroy(&TAux[f]);CHKERRQ(ierr);
754     }
755     ierr = PetscFree2(T, TAux);CHKERRQ(ierr);
756   }
757   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
758   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
759   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
760   ierr = DMDestroy(&plex);CHKERRQ(ierr);
761   ierr = DMDestroy(&plexIn);CHKERRQ(ierr);
762   if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);}
763   PetscFunctionReturn(0);
764 }
765 
766 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
767 {
768   PetscErrorCode ierr;
769 
770   PetscFunctionBegin;
771   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
772   PetscFunctionReturn(0);
773 }
774 
775 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)
776 {
777   PetscErrorCode ierr;
778 
779   PetscFunctionBegin;
780   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
781   PetscFunctionReturn(0);
782 }
783 
784 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
785                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
786                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
787                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
788                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
789                                         InsertMode mode, Vec localX)
790 {
791   PetscErrorCode ierr;
792 
793   PetscFunctionBegin;
794   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
795   PetscFunctionReturn(0);
796 }
797 
798 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
799                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
800                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
801                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
802                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
803                                              InsertMode mode, Vec localX)
804 {
805   PetscErrorCode ierr;
806 
807   PetscFunctionBegin;
808   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
809   PetscFunctionReturn(0);
810 }
811