xref: /petsc/src/dm/impls/plex/plexproject.c (revision 3fed57382a69ddef8b301aa4ce9b2f05bf867c00)
1 #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2 
3 #include <petsc/private/petscfeimpl.h>
4 
5 /*@
6   DMPlexGetActivePoint - Get the point on which projection is currently working
7 
8   Not Collective
9 
10   Input Parameter:
11 . dm - the `DM`
12 
13   Output Parameter:
14 . point - The mesh point involved in the current projection
15 
16   Level: developer
17 
18 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`
19 @*/
20 PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point)
21 {
22   PetscFunctionBeginHot;
23   *point = ((DM_Plex *)dm->data)->activePoint;
24   PetscFunctionReturn(PETSC_SUCCESS);
25 }
26 
27 /*@
28   DMPlexSetActivePoint - Set the point on which projection is currently working
29 
30   Not Collective
31 
32   Input Parameters:
33 + dm    - the `DM`
34 - point - The mesh point involved in the current projection
35 
36   Level: developer
37 
38 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetActivePoint()`
39 @*/
40 PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point)
41 {
42   PetscFunctionBeginHot;
43   ((DM_Plex *)dm->data)->activePoint = point;
44   PetscFunctionReturn(PETSC_SUCCESS);
45 }
46 
47 /*
48   DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point
49 
50   Input Parameters:
51 + dm     - The output `DM`
52 . ds     - The output `DS`
53 . dmIn   - The input `DM`
54 . dsIn   - The input `DS`
55 . time   - The time for this evaluation
56 . fegeom - The FE geometry for this point
57 . fvgeom - The FV geometry for this point
58 . isFE   - Flag indicating whether each output field has an FE discretization
59 . sp     - The output `PetscDualSpace` for each field
60 . funcs  - The evaluation function for each field
61 - ctxs   - The user context for each field
62 
63   Output Parameter:
64 . values - The value for each dual basis vector in the output dual space
65 
66   Level: developer
67 
68 .seealso:[](ch_unstructured), `DM`, `DMPLEX`, `PetscDS`, `PetscFEGeom`, `PetscFVCellGeom`, `PetscDualSpace`
69 */
70 static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, PetscScalar values[])
71 {
72   PetscInt  debug = ((DM_Plex *)dm->data)->printProject;
73   PetscInt  coordDim, Nf, *Nc, f, spDim, d, v, tp;
74   PetscBool isAffine, isCohesive, transform;
75 
76   PetscFunctionBeginHot;
77   PetscCall(DMGetCoordinateDim(dmIn, &coordDim));
78   PetscCall(DMHasBasisTransform(dmIn, &transform));
79   PetscCall(PetscDSGetNumFields(ds, &Nf));
80   PetscCall(PetscDSGetComponents(ds, &Nc));
81   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
82   /* Get values for closure */
83   isAffine = fegeom->isAffine;
84   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
85     void *const ctx = ctxs ? ctxs[f] : NULL;
86     PetscBool   cohesive;
87 
88     if (!sp[f]) continue;
89     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
90     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
91     if (funcs[f]) {
92       if (isFE[f]) {
93         PetscQuadrature  allPoints;
94         PetscInt         q, dim, numPoints;
95         const PetscReal *points;
96         PetscScalar     *pointEval;
97         PetscReal       *x;
98         DM               rdm;
99 
100         PetscCall(PetscDualSpaceGetDM(sp[f], &rdm));
101         PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
102         PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
103         PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
104         PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x));
105         PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f]));
106         for (q = 0; q < numPoints; q++, tp++) {
107           const PetscReal *v0;
108 
109           if (isAffine) {
110             const PetscReal *refpoint    = &points[q * dim];
111             PetscReal        injpoint[3] = {0., 0., 0.};
112 
113             if (dim != fegeom->dim) {
114               if (isCohesive) {
115                 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
116                 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
117                 refpoint = injpoint;
118               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim);
119             }
120             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x);
121             v0 = x;
122           } else {
123             v0 = &fegeom->v[tp * coordDim];
124           }
125           if (transform) {
126             PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx));
127             v0 = x;
128           }
129           if (debug > 3) {
130             PetscInt ap;
131             PetscCall(DMPlexGetActivePoint(dm, &ap));
132             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Project point %" PetscInt_FMT ", analytic: ref (", ap));
133             for (PetscInt d = 0; d < dim; ++d) {
134               if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
135               PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)points[q * dim + d]));
136             }
137             PetscCall(PetscPrintf(PETSC_COMM_SELF, ") real ("));
138             for (PetscInt d = 0; d < dim; ++d) {
139               if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
140               PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)v0[d]));
141             }
142             PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
143           }
144           PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx));
145         }
146         /* Transform point evaluations pointEval[q,c] */
147         PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
148         PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
149         PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x));
150         PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
151         v += spDim;
152         if (isCohesive && !cohesive) {
153           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
154         }
155       } else {
156         for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
157       }
158     } else {
159       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
160       if (isCohesive && !cohesive) {
161         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
162       }
163     }
164   }
165   PetscFunctionReturn(PETSC_SUCCESS);
166 }
167 
168 /*
169   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
170 
171   Input Parameters:
172 + dm             - The output DM
173 . ds             - The output DS
174 . dmIn           - The input DM
175 . dsIn           - The input DS
176 . dmAux          - The auxiliary DM, which is always for the input space
177 . dsAux          - The auxiliary DS, which is always for the input space
178 . time           - The time for this evaluation
179 . localU         - The local solution
180 . localA         - The local auziliary fields
181 . cgeom          - The FE geometry for this point
182 . sp             - The output PetscDualSpace for each field
183 . p              - The point in the output DM
184 . T              - Input basis and derivatives for each field tabulated on the quadrature points
185 . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
186 . funcs          - The evaluation function for each field
187 - ctxs           - The user context for each field
188 
189   Output Parameter:
190 . values         - The value for each dual basis vector in the output dual space
191 
192   Level: developer
193 
194   Note:
195   Not supported for FV
196 
197 .seealso: `DMProjectPoint_Field_Private()`
198 */
199 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, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[])
200 {
201   PetscSection       section, sectionAux = NULL;
202   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
203   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
204   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
205   const PetscScalar *constants;
206   PetscReal         *x;
207   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
208   PetscFEGeom        fegeom;
209   const PetscInt     dE = cgeom->dimEmbed, *cone, *ornt;
210   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
211   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
212   DMPolytopeType     qct;
213 
214   PetscFunctionBeginHot;
215   PetscCall(PetscDSGetNumFields(ds, &Nf));
216   PetscCall(PetscDSGetComponents(ds, &Nc));
217   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
218   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
219   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
220   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
221   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
222   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
223   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
224   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
225   PetscCall(DMHasBasisTransform(dmIn, &transform));
226   PetscCall(DMGetLocalSection(dmIn, &section));
227   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
228   // Get cohesive cell hanging off face
229   if (isCohesiveIn) {
230     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
231     if ((qct != DM_POLYTOPE_POINT_PRISM_TENSOR) && (qct != DM_POLYTOPE_SEG_PRISM_TENSOR) && (qct != DM_POLYTOPE_TRI_PRISM_TENSOR) && (qct != DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
232       DMPolytopeType  ct;
233       const PetscInt *support;
234       PetscInt        Ns, s;
235 
236       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
237       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
238       for (s = 0; s < Ns; ++s) {
239         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
240         if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
241           inp = support[s];
242           break;
243         }
244       }
245       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
246       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
247       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
248       face[0] = 0;
249       face[1] = 0;
250     }
251   }
252   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
253   if (dmAux) {
254     PetscInt subp;
255 
256     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
257     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
258     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
259     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
260     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
261     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
262     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
263   }
264   /* Get values for closure */
265   isAffine        = cgeom->isAffine;
266   fegeom.dim      = cgeom->dim;
267   fegeom.dimEmbed = cgeom->dimEmbed;
268   if (isAffine) {
269     fegeom.v    = x;
270     fegeom.xi   = cgeom->xi;
271     fegeom.J    = cgeom->J;
272     fegeom.invJ = cgeom->invJ;
273     fegeom.detJ = cgeom->detJ;
274   }
275   for (f = 0, v = 0; f < Nf; ++f) {
276     PetscQuadrature  allPoints;
277     PetscInt         q, dim, numPoints;
278     const PetscReal *points;
279     PetscScalar     *pointEval;
280     PetscBool        cohesive;
281     DM               dm;
282 
283     if (!sp[f]) continue;
284     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
285     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
286     if (!funcs[f]) {
287       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
288       if (isCohesive && !cohesive) {
289         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
290       }
291       continue;
292     }
293     const PetscInt ***perms;
294     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
295     PetscCall(PetscDualSpaceGetSymmetries(sp[f], &perms, NULL));
296     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
297     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
298     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
299     for (q = 0; q < numPoints; ++q, ++tp) {
300       PetscInt qpt[2];
301 
302       if (isCohesiveIn) {
303         qpt[0] = perms ? perms[0][ornt[0]][q] : q;
304         qpt[1] = perms ? perms[0][DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0)][q] : q;
305       }
306       if (isAffine) {
307         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x);
308       } else {
309         fegeom.v    = &cgeom->v[tp * dE];
310         fegeom.J    = &cgeom->J[tp * dE * dE];
311         fegeom.invJ = &cgeom->invJ[tp * dE * dE];
312         fegeom.detJ = &cgeom->detJ[tp];
313       }
314       if (coefficients) {
315         if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
316         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
317       }
318       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
319       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
320       (*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]);
321     }
322     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
323     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
324     v += spDim;
325     /* TODO: For now, set both sides equal, but this should use info from other support cell */
326     if (isCohesive && !cohesive) {
327       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
328     }
329   }
330   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
331   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
332   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
333   PetscFunctionReturn(PETSC_SUCCESS);
334 }
335 
336 static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[])
337 {
338   PetscSection       section, sectionAux = NULL;
339   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
340   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
341   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
342   const PetscScalar *constants;
343   PetscReal         *x;
344   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
345   PetscFEGeom        fegeom, cgeom;
346   const PetscInt     dE = fgeom->dimEmbed, *cone, *ornt;
347   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
348   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
349   DMPolytopeType     qct;
350 
351   PetscFunctionBeginHot;
352   PetscCall(PetscDSGetNumFields(ds, &Nf));
353   PetscCall(PetscDSGetComponents(ds, &Nc));
354   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
355   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
356   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
357   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
358   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
359   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
360   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
361   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
362   PetscCall(DMHasBasisTransform(dmIn, &transform));
363   PetscCall(DMGetLocalSection(dmIn, &section));
364   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
365   // Get cohesive cell hanging off face
366   if (isCohesiveIn) {
367     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
368     if ((qct != DM_POLYTOPE_POINT_PRISM_TENSOR) && (qct != DM_POLYTOPE_SEG_PRISM_TENSOR) && (qct != DM_POLYTOPE_TRI_PRISM_TENSOR) && (qct != DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
369       DMPolytopeType  ct;
370       const PetscInt *support;
371       PetscInt        Ns, s;
372 
373       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
374       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
375       for (s = 0; s < Ns; ++s) {
376         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
377         if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
378           inp = support[s];
379           break;
380         }
381       }
382       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
383       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
384       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
385       face[0] = 0;
386       face[1] = 0;
387     }
388   }
389   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
390   if (dmAux) {
391     PetscInt subp;
392 
393     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
394     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
395     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
396     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
397     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
398     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
399     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
400   }
401   /* Get values for closure */
402   isAffine       = fgeom->isAffine;
403   fegeom.n       = NULL;
404   fegeom.J       = NULL;
405   fegeom.v       = NULL;
406   fegeom.xi      = NULL;
407   cgeom.dim      = fgeom->dim;
408   cgeom.dimEmbed = fgeom->dimEmbed;
409   if (isAffine) {
410     fegeom.v    = x;
411     fegeom.xi   = fgeom->xi;
412     fegeom.J    = fgeom->J;
413     fegeom.invJ = fgeom->invJ;
414     fegeom.detJ = fgeom->detJ;
415     fegeom.n    = fgeom->n;
416 
417     cgeom.J    = fgeom->suppJ[0];
418     cgeom.invJ = fgeom->suppInvJ[0];
419     cgeom.detJ = fgeom->suppDetJ[0];
420   }
421   for (f = 0, v = 0; f < Nf; ++f) {
422     PetscQuadrature  allPoints;
423     PetscInt         q, dim, numPoints;
424     const PetscReal *points;
425     PetscScalar     *pointEval;
426     PetscBool        cohesive;
427     DM               dm;
428 
429     if (!sp[f]) continue;
430     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
431     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
432     if (!funcs[f]) {
433       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
434       if (isCohesive && !cohesive) {
435         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
436       }
437       continue;
438     }
439     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
440     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
441     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
442     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
443     for (q = 0; q < numPoints; ++q, ++tp) {
444       PetscInt qpt[2];
445 
446       if (isCohesiveIn) {
447         // These points are not integration quadratures, but dual space quadratures
448         // If they had multiple points we should match them from both sides, similar to hybrid residual eval
449         qpt[0] = qpt[1] = q;
450       }
451       if (isAffine) {
452         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x);
453       } else {
454         fegeom.v    = &fgeom->v[tp * dE];
455         fegeom.J    = &fgeom->J[tp * dE * dE];
456         fegeom.invJ = &fgeom->invJ[tp * dE * dE];
457         fegeom.detJ = &fgeom->detJ[tp];
458         fegeom.n    = &fgeom->n[tp * dE];
459 
460         cgeom.J    = &fgeom->suppJ[0][tp * dE * dE];
461         cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE];
462         cgeom.detJ = &fgeom->suppDetJ[0][tp];
463       }
464       /* TODO We should use cgeom here, instead of fegeom, however the geometry coming in through fgeom does not have the support cell geometry */
465       if (coefficients) {
466         if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
467         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
468       }
469       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
470       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
471       (*funcs[f])(dE, NfIn, 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]);
472     }
473     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
474     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
475     v += spDim;
476     /* TODO: For now, set both sides equal, but this should use info from other support cell */
477     if (isCohesive && !cohesive) {
478       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
479     }
480   }
481   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
482   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
483   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
484   PetscFunctionReturn(PETSC_SUCCESS);
485 }
486 
487 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[], PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
488 {
489   PetscFVCellGeom fvgeom;
490   PetscInt        dim, dimEmbed;
491 
492   PetscFunctionBeginHot;
493   PetscCall(DMGetDimension(dm, &dim));
494   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
495   if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
496   switch (type) {
497   case DM_BC_ESSENTIAL:
498   case DM_BC_NATURAL:
499     PetscCall(DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))funcs, ctxs, values));
500     break;
501   case DM_BC_ESSENTIAL_FIELD:
502   case DM_BC_NATURAL_FIELD:
503     PetscCall(DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
504     break;
505   case DM_BC_ESSENTIAL_BD_FIELD:
506     PetscCall(DMProjectPoint_BdField_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
507     break;
508   default:
509     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type);
510   }
511   PetscFunctionReturn(PETSC_SUCCESS);
512 }
513 
514 static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
515 {
516   PetscReal *points;
517   PetscInt   f, numPoints;
518 
519   PetscFunctionBegin;
520   if (!dim) {
521     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
522     PetscFunctionReturn(PETSC_SUCCESS);
523   }
524   numPoints = 0;
525   for (f = 0; f < Nf; ++f) {
526     if (funcs[f]) {
527       PetscQuadrature fAllPoints;
528       PetscInt        fNumPoints;
529 
530       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
531       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
532       numPoints += fNumPoints;
533     }
534   }
535   PetscCall(PetscMalloc1(dim * numPoints, &points));
536   numPoints = 0;
537   for (f = 0; f < Nf; ++f) {
538     if (funcs[f]) {
539       PetscQuadrature  fAllPoints;
540       PetscInt         qdim, fNumPoints, q;
541       const PetscReal *fPoints;
542 
543       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
544       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
545       PetscCheck(qdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %" PetscInt_FMT " for dual basis does not match input dimension %" PetscInt_FMT, qdim, dim);
546       for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q];
547       numPoints += fNumPoints;
548     }
549   }
550   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
551   PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL));
552   PetscFunctionReturn(PETSC_SUCCESS);
553 }
554 
555 /*@C
556   DMGetFirstLabeledPoint - Find first labeled `point` in `odm` such that the corresponding point in `dm` has the specified `height`. Return `point` and the corresponding `ds`.
557 
558   Input Parameters:
559 + dm     - the `DM`
560 . odm    - the enclosing `DM`
561 . label  - label for `DM` domain, or `NULL` for whole domain
562 . numIds - the number of `ids`
563 . ids    - An array of the label ids in sequence for the domain
564 - height - Height of target cells in `DMPLEX` topology
565 
566   Output Parameters:
567 + point - the first labeled point
568 - ds    - the `PetscDS` corresponding to the first labeled point
569 
570   Level: developer
571 
572 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS`
573 @*/
574 PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
575 {
576   DM              plex;
577   DMEnclosureType enc;
578   PetscInt        ls = -1;
579 
580   PetscFunctionBegin;
581   if (point) *point = -1;
582   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
583   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
584   PetscCall(DMConvert(dm, DMPLEX, &plex));
585   for (PetscInt i = 0; i < numIds; ++i) {
586     IS       labelIS;
587     PetscInt num_points, pStart, pEnd;
588     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
589     if (!labelIS) continue; /* No points with that id on this process */
590     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
591     PetscCall(ISGetSize(labelIS, &num_points));
592     if (num_points) {
593       const PetscInt *points;
594       PetscCall(ISGetIndices(labelIS, &points));
595       for (PetscInt i = 0; i < num_points; i++) {
596         PetscInt point;
597         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
598         if (pStart <= point && point < pEnd) {
599           ls = point;
600           if (ds) {
601             // If this is a face of a cohesive cell, then prefer that DS
602             if (height == 1) {
603               const PetscInt *supp;
604               PetscInt        suppSize;
605               DMPolytopeType  ct;
606 
607               PetscCall(DMPlexGetSupport(dm, ls, &supp));
608               PetscCall(DMPlexGetSupportSize(dm, ls, &suppSize));
609               for (PetscInt s = 0; s < suppSize; ++s) {
610                 PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
611                 if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
612                   ls = supp[s];
613                   break;
614                 }
615               }
616             }
617             PetscCall(DMGetCellDS(dm, ls, ds, NULL));
618           }
619           if (ls >= 0) break;
620         }
621       }
622       PetscCall(ISRestoreIndices(labelIS, &points));
623     }
624     PetscCall(ISDestroy(&labelIS));
625     if (ls >= 0) break;
626   }
627   if (point) *point = ls;
628   PetscCall(DMDestroy(&plex));
629   PetscFunctionReturn(PETSC_SUCCESS);
630 }
631 
632 /*
633   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
634 
635   There are several different scenarios:
636 
637   1) Volumetric mesh with volumetric auxiliary data
638 
639      Here minHeight=0 since we loop over cells.
640 
641   2) Boundary mesh with boundary auxiliary data
642 
643      Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation.
644 
645   3) Volumetric mesh with boundary auxiliary data
646 
647      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.
648 
649   4) Volumetric input mesh with boundary output mesh
650 
651      Here we must get a subspace for the input DS
652 
653   The maxHeight is used to support enforcement of constraints in DMForest.
654 
655   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
656 
657   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.
658     - 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
659       is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When the effective height is nonzero, we need to extract
660       dual spaces for the boundary from our input spaces.
661     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
662 
663   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
664 
665   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.
666 */
667 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, InsertMode mode, Vec localX)
668 {
669   DM               plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
670   DMEnclosureType  encIn, encAux;
671   PetscDS          ds = NULL, dsIn = NULL, dsAux = NULL;
672   Vec              localA = NULL, tv;
673   IS               fieldIS;
674   PetscSection     section;
675   PetscDualSpace  *sp, *cellsp, *spIn, *cellspIn;
676   PetscTabulation *T = NULL, *TAux = NULL;
677   PetscInt        *Nc;
678   PetscInt         dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, minHeightIn, minHeightAux, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
679   PetscBool       *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform;
680   DMField          coordField;
681   DMLabel          depthLabel;
682   PetscQuadrature  allPoints = NULL;
683 
684   PetscFunctionBegin;
685   if (localU) PetscCall(VecGetDM(localU, &dmIn));
686   else dmIn = dm;
687   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
688   if (localA) PetscCall(VecGetDM(localA, &dmAux));
689   else dmAux = NULL;
690   PetscCall(DMConvert(dm, DMPLEX, &plex));
691   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
692   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
693   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
694   PetscCall(DMGetDimension(dm, &dim));
695   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
696   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
697   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
698   PetscCall(DMHasBasisTransform(dm, &transform));
699   /* Auxiliary information can only be used with interpolation of field functions */
700   if (dmAux) {
701     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
702     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) PetscCheck(localA, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector");
703   }
704   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plexIn, PETSC_TRUE, localU, time, NULL, NULL, NULL));
705   PetscCall(DMGetCoordinateField(dm, &coordField));
706   PetscCheck(coordField, PETSC_COMM_SELF, PETSC_ERR_USER, "DM must have a coordinate field");
707   /**** No collective calls below this point ****/
708   /* Determine height for iteration of all meshes */
709   {
710     DMPolytopeType ct, ctIn, ctAux;
711     PetscInt       lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux;
712     PetscInt       dim = -1, dimIn = -1, dimAux = -1;
713 
714     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
715     if (pEnd > pStart) {
716       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
717       p = lStart < 0 ? pStart : lStart;
718       PetscCall(DMPlexGetCellType(plex, p, &ct));
719       dim = DMPolytopeTypeGetDim(ct);
720       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
721       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
722       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
723       dimIn = DMPolytopeTypeGetDim(ctIn);
724       if (dmAux) {
725         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
726         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
727         if (pStartAux < pEndAux) {
728           PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
729           dimAux = DMPolytopeTypeGetDim(ctAux);
730         }
731       } else dimAux = dim;
732     } else {
733       PetscCall(DMDestroy(&plex));
734       PetscCall(DMDestroy(&plexIn));
735       if (dmAux) PetscCall(DMDestroy(&plexAux));
736       PetscFunctionReturn(PETSC_SUCCESS);
737     }
738     if (dim < 0) {
739       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
740 
741       /* Fall back to determination based on being a submesh */
742       PetscCall(DMPlexGetSubpointMap(plex, &spmap));
743       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
744       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
745       dim    = spmap ? 1 : 0;
746       dimIn  = spmapIn ? 1 : 0;
747       dimAux = spmapAux ? 1 : 0;
748     }
749     {
750       PetscInt dimProj   = PetscMin(PetscMin(dim, dimIn), dimAux < 0 ? PETSC_INT_MAX : dimAux);
751       PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;
752 
753       PetscCheck(PetscAbsInt(dimProj - dim) <= 1 && PetscAbsInt(dimProj - dimIn) <= 1 && PetscAbsInt(dimProj - dimAuxEff) <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not currently support differences of more than 1 in dimension");
754       if (dimProj < dim) minHeight = 1;
755       htInc    = dim - dimProj;
756       htIncIn  = dimIn - dimProj;
757       htIncAux = dimAuxEff - dimProj;
758     }
759   }
760   PetscCall(DMPlexGetDepth(plex, &depth));
761   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
762   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
763   maxHeight = PetscMax(maxHeight, minHeight);
764   PetscCheck(maxHeight >= 0 && maxHeight <= dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", maxHeight, dim);
765   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, NULL, &ds));
766   if (!ds) PetscCall(DMGetDS(dm, &ds));
767   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, minHeight, NULL, &dsIn));
768   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
769   PetscCall(PetscDSGetNumFields(ds, &Nf));
770   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
771   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
772   if (isCohesiveIn) --htIncIn; // Should be rearranged
773   PetscCall(DMGetNumFields(dm, &NfTot));
774   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
775   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL));
776   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
777   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
778   PetscCall(DMGetLocalSection(dm, &section));
779   if (dmAux) {
780     PetscCall(DMGetDS(dmAux, &dsAux));
781     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
782   }
783   PetscCall(PetscDSGetComponents(ds, &Nc));
784   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
785   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
786   else {
787     cellsp   = sp;
788     cellspIn = spIn;
789   }
790   /* Get cell dual spaces */
791   for (f = 0; f < Nf; ++f) {
792     PetscDiscType disctype;
793 
794     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
795     if (disctype == PETSC_DISC_FE) {
796       PetscFE fe;
797 
798       isFE[f] = PETSC_TRUE;
799       hasFE   = PETSC_TRUE;
800       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
801       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
802     } else if (disctype == PETSC_DISC_FV) {
803       PetscFV fv;
804 
805       isFE[f] = PETSC_FALSE;
806       hasFV   = PETSC_TRUE;
807       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv));
808       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
809     } else {
810       isFE[f]   = PETSC_FALSE;
811       cellsp[f] = NULL;
812     }
813   }
814   for (f = 0; f < NfIn; ++f) {
815     PetscDiscType disctype;
816 
817     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
818     if (disctype == PETSC_DISC_FE) {
819       PetscFE fe;
820 
821       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
822       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
823     } else if (disctype == PETSC_DISC_FV) {
824       PetscFV fv;
825 
826       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv));
827       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
828     } else {
829       cellspIn[f] = NULL;
830     }
831   }
832   for (f = 0; f < Nf; ++f) {
833     if (!htInc) {
834       sp[f] = cellsp[f];
835     } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
836   }
837   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
838     PetscFE          fem, subfem;
839     PetscDiscType    disctype;
840     const PetscReal *points;
841     PetscInt         numPoints, k;
842 
843     PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
844     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints));
845     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
846     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
847     for (f = 0; f < NfIn; ++f) {
848       if (!htIncIn) {
849         spIn[f] = cellspIn[f];
850       } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
851 
852       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
853       if (disctype != PETSC_DISC_FE) continue;
854       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem));
855       if (!htIncIn) {
856         subfem = fem;
857       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
858       PetscCall(PetscDSGetJetDegree(dsIn, f, &k));
859       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &T[f]));
860     }
861     for (f = 0; f < NfAux; ++f) {
862       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
863       if (disctype != PETSC_DISC_FE) continue;
864       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem));
865       if (!htIncAux) {
866         subfem = fem;
867       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
868       PetscCall(PetscDSGetJetDegree(dsAux, f, &k));
869       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &TAux[f]));
870     }
871   }
872   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
873   for (h = minHeight; h <= maxHeight; h++) {
874     PetscInt     hEff     = h - minHeight + htInc;
875     PetscInt     hEffIn   = h - minHeight + htIncIn;
876     PetscInt     hEffAux  = h - minHeight + htIncAux;
877     PetscDS      dsEff    = ds;
878     PetscDS      dsEffIn  = dsIn;
879     PetscDS      dsEffAux = dsAux;
880     PetscScalar *values;
881     PetscBool   *fieldActive;
882     PetscInt     maxDegree;
883     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
884     IS           heightIS;
885 
886     if (h > minHeight) {
887       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
888     }
889     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
890     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
891     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
892     if (pEnd <= pStart) {
893       PetscCall(ISDestroy(&heightIS));
894       continue;
895     }
896     /* Compute totDim, the number of dofs in the closure of a point at this height */
897     totDim = 0;
898     for (f = 0; f < Nf; ++f) {
899       PetscBool cohesive;
900 
901       if (!sp[f]) continue;
902       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
903       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
904       totDim += spDim;
905       if (isCohesive && !cohesive) totDim += spDim;
906     }
907     p = lStart < 0 ? pStart : lStart;
908     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
909     PetscCheck(numValues == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The output section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT " in [%" PetscInt_FMT ", %" PetscInt_FMT "]", p, numValues, totDim, h, minHeight, maxHeight);
910     if (!totDim) {
911       PetscCall(ISDestroy(&heightIS));
912       continue;
913     }
914     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
915     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
916     if (localU) {
917       PetscInt totDimIn, pIn, numValuesIn;
918 
919       totDimIn = 0;
920       for (f = 0; f < NfIn; ++f) {
921         PetscBool cohesive;
922 
923         if (!spIn[f]) continue;
924         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
925         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
926         totDimIn += spDim;
927         if (isCohesiveIn && !cohesive) totDimIn += spDim;
928       }
929       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
930       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
931       // TODO We could check that pIn is a cohesive cell for this check
932       PetscCheck(isCohesiveIn || (numValuesIn == totDimIn), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The input section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT, pIn, numValuesIn, totDimIn, htIncIn);
933       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
934     }
935     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
936     /* Loop over points at this height */
937     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
938     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
939     {
940       const PetscInt *fields;
941 
942       PetscCall(ISGetIndices(fieldIS, &fields));
943       for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE;
944       for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
945       PetscCall(ISRestoreIndices(fieldIS, &fields));
946     }
947     if (label) {
948       PetscInt i;
949 
950       for (i = 0; i < numIds; ++i) {
951         IS              pointIS, isectIS;
952         const PetscInt *points;
953         PetscInt        n;
954         PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
955         PetscQuadrature quad = NULL;
956 
957         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
958         if (!pointIS) continue; /* No points with that id on this process */
959         PetscCall(ISIntersect(pointIS, heightIS, &isectIS));
960         PetscCall(ISDestroy(&pointIS));
961         if (!isectIS) continue;
962         PetscCall(ISGetLocalSize(isectIS, &n));
963         PetscCall(ISGetIndices(isectIS, &points));
964         PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree));
965         if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad));
966         if (!quad) {
967           if (!h && allPoints) {
968             quad      = allPoints;
969             allPoints = NULL;
970           } else {
971             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad));
972           }
973         }
974         PetscFEGeomMode geommode = htInc && h == minHeight ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC;
975 
976         if (n) {
977           PetscInt depth, dep;
978 
979           PetscCall(DMPlexGetDepth(dm, &depth));
980           PetscCall(DMPlexGetPointDepth(dm, points[0], &dep));
981           if (dep < depth && h == minHeight) geommode = PETSC_FEGEOM_BOUNDARY;
982         }
983         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, geommode, &fegeom));
984         for (p = 0; p < n; ++p) {
985           const PetscInt point = points[p];
986 
987           PetscCall(PetscArrayzero(values, numValues));
988           PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom));
989           PetscCall(DMPlexSetActivePoint(dm, point));
990           PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values));
991           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
992           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
993         }
994         PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom));
995         PetscCall(PetscFEGeomDestroy(&fegeom));
996         PetscCall(PetscQuadratureDestroy(&quad));
997         PetscCall(ISRestoreIndices(isectIS, &points));
998         PetscCall(ISDestroy(&isectIS));
999       }
1000     } else {
1001       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
1002       PetscQuadrature quad = NULL;
1003       IS              pointIS;
1004 
1005       PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS));
1006       PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
1007       if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
1008       if (!quad) {
1009         if (!h && allPoints) {
1010           quad      = allPoints;
1011           allPoints = NULL;
1012         } else {
1013           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad));
1014         }
1015       }
1016       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC, &fegeom));
1017       for (p = pStart; p < pEnd; ++p) {
1018         PetscCall(PetscArrayzero(values, numValues));
1019         PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom));
1020         PetscCall(DMPlexSetActivePoint(dm, p));
1021         PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values));
1022         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
1023         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
1024       }
1025       PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom));
1026       PetscCall(PetscFEGeomDestroy(&fegeom));
1027       PetscCall(PetscQuadratureDestroy(&quad));
1028       PetscCall(ISDestroy(&pointIS));
1029     }
1030     PetscCall(ISDestroy(&heightIS));
1031     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
1032     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
1033   }
1034   /* Cleanup */
1035   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
1036     for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f]));
1037     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
1038     PetscCall(PetscFree2(T, TAux));
1039   }
1040   PetscCall(PetscQuadratureDestroy(&allPoints));
1041   PetscCall(PetscFree3(isFE, sp, spIn));
1042   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
1043   PetscCall(DMDestroy(&plex));
1044   PetscCall(DMDestroy(&plexIn));
1045   if (dmAux) PetscCall(DMDestroy(&plexAux));
1046   PetscFunctionReturn(PETSC_SUCCESS);
1047 }
1048 
1049 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
1050 {
1051   PetscFunctionBegin;
1052   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
1053   PetscFunctionReturn(PETSC_SUCCESS);
1054 }
1055 
1056 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)
1057 {
1058   PetscFunctionBegin;
1059   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
1060   PetscFunctionReturn(PETSC_SUCCESS);
1061 }
1062 
1063 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
1064 {
1065   PetscFunctionBegin;
1066   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
1067   PetscFunctionReturn(PETSC_SUCCESS);
1068 }
1069 
1070 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
1071 {
1072   PetscFunctionBegin;
1073   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
1074   PetscFunctionReturn(PETSC_SUCCESS);
1075 }
1076 
1077 PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
1078 {
1079   PetscFunctionBegin;
1080   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX));
1081   PetscFunctionReturn(PETSC_SUCCESS);
1082 }
1083