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