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