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