xref: /petsc/src/dm/impls/plex/plexproject.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c) !
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 Argument:
11 . dm   - the DM
12 
13   Output Argument:
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 Arguments:
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 derviatives 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   The maxHeight is used to support enforcement of constraints in DMForest.
522 
523   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
524 
525   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.
526     - 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
527       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
528       dual spaces for the boundary from our input spaces.
529     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
530 
531   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
532 
533   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.
534 */
535 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
536                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
537                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
538                                                   InsertMode mode, Vec localX)
539 {
540   DM                 plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
541   DMEnclosureType    encIn, encAux;
542   PetscDS            ds = NULL, dsIn = NULL, dsAux = NULL;
543   Vec                localA = NULL, tv;
544   IS                 fieldIS;
545   PetscSection       section;
546   PetscDualSpace    *sp, *cellsp;
547   PetscTabulation *T = NULL, *TAux = NULL;
548   PetscInt          *Nc;
549   PetscInt           dim, dimEmbed, depth, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
550   PetscBool         *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, isHybrid = PETSC_FALSE, transform;
551   DMField            coordField;
552   DMLabel            depthLabel;
553   PetscQuadrature    allPoints = NULL;
554   PetscErrorCode     ierr;
555 
556   PetscFunctionBegin;
557   if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);}
558   else        {dmIn = dm;}
559   ierr = DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, &localA);CHKERRQ(ierr);
560   if (localA) {ierr = VecGetDM(localA, &dmAux);CHKERRQ(ierr);} else {dmAux = NULL;}
561   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
562   ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr);
563   ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr);
564   ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
565   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
566   ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr);
567   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
568   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
569   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
570   /* Auxiliary information can only be used with interpolation of field functions */
571   if (dmAux) {
572     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
573     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
574       if (!localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
575       if (!minHeight) {
576         DMLabel spmap;
577 
578         /* If dmAux is a surface, then force the projection to take place over a surface */
579         ierr = DMPlexGetSubpointMap(plexAux, &spmap);CHKERRQ(ierr);
580         if (spmap) {
581           ierr = DMPlexGetVTKCellHeight(plexAux, &minHeight);CHKERRQ(ierr);
582           auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
583         }
584       }
585     }
586   }
587   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
588   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
589   ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr);
590   maxHeight = PetscMax(maxHeight, minHeight);
591   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
592   ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr);
593   if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);}
594   ierr = DMGetFirstLabelEntry_Private(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr);
595   if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);}
596   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
597   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
598   ierr = DMGetNumFields(dm, &NfTot);CHKERRQ(ierr);
599   ierr = DMFindRegionNum(dm, ds, &regionNum);CHKERRQ(ierr);
600   ierr = DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL);CHKERRQ(ierr);
601   ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr);
602   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
603   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
604   if (dmAux) {
605     ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr);
606     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
607   }
608   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
609   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
610   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
611   else               {cellsp = sp;}
612   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
613   /* Get cell dual spaces */
614   for (f = 0; f < Nf; ++f) {
615     PetscDiscType disctype;
616 
617     ierr = PetscDSGetDiscType_Internal(ds, f, &disctype);CHKERRQ(ierr);
618     if (disctype == PETSC_DISC_FE) {
619       PetscFE fe;
620 
621       isFE[f] = PETSC_TRUE;
622       hasFE   = PETSC_TRUE;
623       ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
624       ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
625     } else if (disctype == PETSC_DISC_FV) {
626       PetscFV fv;
627 
628       isFE[f] = PETSC_FALSE;
629       hasFV   = PETSC_TRUE;
630       ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fv);CHKERRQ(ierr);
631       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
632     } else {
633       isFE[f]   = PETSC_FALSE;
634       cellsp[f] = NULL;
635     }
636   }
637   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
638   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
639     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
640     PetscFE          fem, subfem;
641     PetscDiscType    disctype;
642     const PetscReal *points;
643     PetscInt         numPoints;
644 
645     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
646     for (f = 0; f < Nf; ++f) {
647       if (!effectiveHeight) {sp[f] = cellsp[f];}
648       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
649     }
650     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
651     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
652     ierr = PetscMalloc2(NfIn, &T, NfAux, &TAux);CHKERRQ(ierr);
653     for (f = 0; f < NfIn; ++f) {
654       ierr = PetscDSGetDiscType_Internal(dsIn, f, &disctype);CHKERRQ(ierr);
655       if (disctype != PETSC_DISC_FE) continue;
656       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
657       if (!effectiveHeight) {subfem = fem;}
658       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
659       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]);CHKERRQ(ierr);
660     }
661     for (f = 0; f < NfAux; ++f) {
662       ierr = PetscDSGetDiscType_Internal(dsAux, f, &disctype);CHKERRQ(ierr);
663       if (disctype != PETSC_DISC_FE) continue;
664       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
665       if (!effectiveHeight || auxBd) {subfem = fem;}
666       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
667       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]);CHKERRQ(ierr);
668     }
669   }
670   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
671   for (h = minHeight; h <= maxHeight; h++) {
672     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
673     PetscDS      dsEff         = ds;
674     PetscScalar *values;
675     PetscBool   *fieldActive;
676     PetscInt     maxDegree;
677     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
678     IS           heightIS;
679 
680     /* Note we assume that dm and dmIn share the same topology */
681     ierr = DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd);CHKERRQ(ierr);
682     ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
683     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
684     if (pEnd <= pStart) {
685       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
686       continue;
687     }
688     /* Compute totDim, the number of dofs in the closure of a point at this height */
689     totDim = 0;
690     for (f = 0; f < Nf; ++f) {
691       if (!effectiveHeight) {
692         sp[f] = cellsp[f];
693       } else {
694         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
695       }
696       if (!sp[f]) continue;
697       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
698       totDim += spDim;
699       if (isHybrid && (f < Nf-1)) totDim += spDim;
700     }
701     ierr = DMPlexVecGetClosure(plex, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr);
702     if (numValues != totDim) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point (%D) closure size %D != dual space dimension %D", lStart < 0 ? pStart : lStart, numValues, totDim);
703     if (!totDim) {
704       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
705       continue;
706     }
707     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(ds, effectiveHeight, &dsEff);CHKERRQ(ierr);}
708     /* Loop over points at this height */
709     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
710     ierr = DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive);CHKERRQ(ierr);
711     {
712       const PetscInt *fields;
713 
714       ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);
715       for (f = 0; f < NfTot; ++f) {fieldActive[f] = PETSC_FALSE;}
716       for (f = 0; f < Nf; ++f) {fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;}
717       ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);
718     }
719     if (label) {
720       PetscInt i;
721 
722       for (i = 0; i < numIds; ++i) {
723         IS              pointIS, isectIS;
724         const PetscInt *points;
725         PetscInt        n;
726         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
727         PetscQuadrature quad = NULL;
728 
729         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
730         if (!pointIS) continue; /* No points with that id on this process */
731         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
732         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
733         if (!isectIS) continue;
734         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
735         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
736         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
737         if (maxDegree <= 1) {
738           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
739         }
740         if (!quad) {
741           if (!h && allPoints) {
742             quad = allPoints;
743             allPoints = NULL;
744           } else {
745             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,isHybrid ? dim-effectiveHeight-1 : dim-effectiveHeight,funcs,&quad);CHKERRQ(ierr);
746           }
747         }
748         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
749         for (p = 0; p < n; ++p) {
750           const PetscInt  point = points[p];
751 
752           ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
753           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
754           ierr = DMPlexSetActivePoint(dm, point);CHKERRQ(ierr);
755           ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values);
756           if (ierr) {
757             PetscErrorCode ierr2;
758             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
759             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
760             CHKERRQ(ierr);
761           }
762           if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
763           ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode);CHKERRQ(ierr);
764         }
765         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
766         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
767         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
768         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
769         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
770       }
771     } else {
772       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
773       PetscQuadrature quad = NULL;
774       IS              pointIS;
775 
776       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
777       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
778       if (maxDegree <= 1) {
779         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
780       }
781       if (!quad) {
782         if (!h && allPoints) {
783           quad = allPoints;
784           allPoints = NULL;
785         } else {
786           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&quad);CHKERRQ(ierr);
787         }
788       }
789       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
790       for (p = pStart; p < pEnd; ++p) {
791         ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
792         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
793         ierr = DMPlexSetActivePoint(dm, p);CHKERRQ(ierr);
794         ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values);
795         if (ierr) {
796           PetscErrorCode ierr2;
797           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
798           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
799           CHKERRQ(ierr);
800         }
801         if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
802         ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode);CHKERRQ(ierr);
803       }
804       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
805       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
806       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
807       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
808     }
809     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
810     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
811     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
812   }
813   /* Cleanup */
814   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
815     PetscInt effectiveHeight = auxBd ? minHeight : 0;
816     PetscFE  fem, subfem;
817 
818     for (f = 0; f < NfIn; ++f) {
819       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
820       if (!effectiveHeight) {subfem = fem;}
821       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
822       ierr = PetscTabulationDestroy(&T[f]);CHKERRQ(ierr);
823     }
824     for (f = 0; f < NfAux; ++f) {
825       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
826       if (!effectiveHeight || auxBd) {subfem = fem;}
827       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
828       ierr = PetscTabulationDestroy(&TAux[f]);CHKERRQ(ierr);
829     }
830     ierr = PetscFree2(T, TAux);CHKERRQ(ierr);
831   }
832   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
833   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
834   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
835   ierr = DMDestroy(&plex);CHKERRQ(ierr);
836   ierr = DMDestroy(&plexIn);CHKERRQ(ierr);
837   if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);}
838   PetscFunctionReturn(0);
839 }
840 
841 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
842 {
843   PetscErrorCode ierr;
844 
845   PetscFunctionBegin;
846   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
847   PetscFunctionReturn(0);
848 }
849 
850 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)
851 {
852   PetscErrorCode ierr;
853 
854   PetscFunctionBegin;
855   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
856   PetscFunctionReturn(0);
857 }
858 
859 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
860                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
861                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
862                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
863                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
864                                         InsertMode mode, Vec localX)
865 {
866   PetscErrorCode ierr;
867 
868   PetscFunctionBegin;
869   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
870   PetscFunctionReturn(0);
871 }
872 
873 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
874                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
875                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
876                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
877                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
878                                              InsertMode mode, Vec localX)
879 {
880   PetscErrorCode ierr;
881 
882   PetscFunctionBegin;
883   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
884   PetscFunctionReturn(0);
885 }
886 
887 PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
888                                                void (**funcs)(PetscInt, PetscInt, PetscInt,
889                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
890                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
891                                                               PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
892                                                InsertMode mode, Vec localX)
893 {
894   PetscErrorCode ierr;
895 
896   PetscFunctionBegin;
897   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);
898   PetscFunctionReturn(0);
899 }
900