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