xref: /petsc/src/dm/impls/plex/plexproject.c (revision a7631e779ee37285eb510131acb81bacea428671)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #include <petsc/private/petscfeimpl.h>
4 
5 static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscReal time, PetscFECellGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
6                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
7                                                   PetscScalar values[])
8 {
9   PetscDS        prob;
10   PetscInt       Nf, *Nc, f, totDim, spDim, d, v;
11   PetscErrorCode ierr;
12 
13   PetscFunctionBeginHot;
14   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
15   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
16   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
17   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
18   /* Get values for closure */
19   for (f = 0, v = 0; f < Nf; ++f) {
20     void * const ctx = ctxs ? ctxs[f] : NULL;
21 
22     if (!sp[f]) continue;
23     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
24     for (d = 0; d < spDim; ++d, ++v) {
25       if (funcs[f]) {
26         if (isFE[f]) {ierr = PetscDualSpaceApply(sp[f], d, time, fegeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
27         else         {ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);}
28       } else {
29         values[v] = 0.0;
30       }
31     }
32   }
33   PetscFunctionReturn(0);
34 }
35 
36 static PetscErrorCode DMProjectPoint_Field_Private(DM dm, DM dmAux, PetscReal time, Vec localU, Vec localA, PetscFECellGeom *fegeom, PetscDualSpace sp[], PetscInt p,
37                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
38                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
39                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
40                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
41                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
42                                                    PetscScalar values[])
43 {
44   PetscDS            prob, probAux = NULL;
45   PetscSection       section, sectionAux = NULL;
46   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL, *bc;
47   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
48   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
49   const PetscScalar *constants;
50   PetscReal         *x;
51   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
52   PetscInt           dim, dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, c, v, tp = 0;
53   PetscErrorCode     ierr;
54 
55   PetscFunctionBeginHot;
56   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
57   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
58   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
59   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
60   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
61   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
62   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
63   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
64   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
65   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
66   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
67   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
68   if (dmAux) {
69     DMLabel  spmap;
70     PetscInt subp;
71 
72     ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
73     if (spmap) {
74       IS              subpointIS;
75       const PetscInt *subpoints;
76       PetscInt        numSubpoints;
77 
78       ierr = DMPlexCreateSubpointIS(dmAux, &subpointIS);CHKERRQ(ierr);
79       ierr = ISGetLocalSize(subpointIS, &numSubpoints);CHKERRQ(ierr);
80       ierr = ISGetIndices(subpointIS, &subpoints);CHKERRQ(ierr);
81       ierr = PetscFindInt(p, numSubpoints, subpoints, &subp);CHKERRQ(ierr);
82       if (subp < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d not found in submesh", p);
83       ierr = ISRestoreIndices(subpointIS, &subpoints);CHKERRQ(ierr);
84       ierr = ISDestroy(&subpointIS);CHKERRQ(ierr);
85     } else subp = p;
86     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
87     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
88     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
89     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
90     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
91     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
92     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
93     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
94     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
95     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
96     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
97   }
98   /* Get values for closure */
99   for (f = 0, v = 0; f < Nf; ++f) {
100     if (!sp[f]) continue;
101     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
102     for (d = 0; d < spDim; ++d, ++v) {
103       if (funcs[f]) {
104         PetscQuadrature  quad;
105         const PetscReal *points, *weights;
106         PetscInt         numPoints, q;
107 
108         ierr = PetscDualSpaceGetFunctional(sp[f], d, &quad);CHKERRQ(ierr);
109         ierr = PetscQuadratureGetData(quad, NULL, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
110         for (q = 0; q < numPoints; ++q, ++tp) {
111           CoordinatesRefToReal(dim, dim, fegeom->v0, fegeom->J, &points[q*dim], x);
112           EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, fegeom->invJ, coefficients, coefficients_t, u, u_x, u_t);
113           if (probAux) {EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, fegeom->invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);}
114           (*funcs[f])(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, x, numConstants, constants, bc);
115           for (c = 0; c < Nc[f]; ++c) values[v] += bc[c]*weights[c];
116         }
117       } else {
118         values[v] = 0.0;
119       }
120     }
121   }
122   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
123   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
124   PetscFunctionReturn(0);
125 }
126 
127 static PetscErrorCode DMProjectPoint_Private(DM dm, DM dmAux, PetscInt h, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
128                                              PetscDualSpace sp[], PetscInt p, PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
129                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
130 {
131   PetscFECellGeom fegeom;
132   PetscFVCellGeom fvgeom;
133   PetscInt        dim, dimEmbed;
134   PetscErrorCode  ierr;
135 
136   PetscFunctionBeginHot;
137   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
138   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
139   if (hasFE) {
140     ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, fegeom.v0, fegeom.J, fegeom.invJ, &fegeom.detJ);CHKERRQ(ierr);
141     fegeom.dim      = dim - h;
142     fegeom.dimEmbed = dimEmbed;
143   }
144   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
145   switch (type) {
146   case DM_BC_ESSENTIAL:
147   case DM_BC_NATURAL:
148     ierr = DMProjectPoint_Func_Private(dm, time, &fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
149   case DM_BC_ESSENTIAL_FIELD:
150   case DM_BC_NATURAL_FIELD:
151     ierr = DMProjectPoint_Field_Private(dm, dmAux, time, localU, localA, &fegeom, sp, p,
152                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
153                                         (void (**)(PetscInt, PetscInt, PetscInt,
154                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
155                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
156                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
157   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
158   }
159   PetscFunctionReturn(0);
160 }
161 
162 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
163                                                   DMLabel label, PetscInt numIds, const PetscInt ids[],
164                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
165                                                   InsertMode mode, Vec localX)
166 {
167   DM              dmAux = NULL;
168   PetscDS         prob, probAux = NULL;
169   Vec             localA = NULL;
170   PetscSection    section;
171   PetscDualSpace *sp, *cellsp;
172   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
173   PetscInt       *Nc;
174   PetscInt        dim, dimEmbed, minHeight, maxHeight, h, Nf, NfAux = 0, f;
175   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE;
176   PetscErrorCode  ierr;
177 
178   PetscFunctionBegin;
179   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
180   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
181   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
182   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
183   if (!minHeight && dmAux) {
184     DMLabel spmap;
185 
186     /* If dmAux is a surface, then force the projection to take place over a surface */
187     ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
188     if (spmap) {minHeight = 1; auxBd = PETSC_TRUE;}
189   }
190   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
191   maxHeight = PetscMax(maxHeight, minHeight);
192   if (maxHeight < 0 || maxHeight > dim) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %d not in [0, %d)\n", maxHeight, dim);}
193   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
194   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
195   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
196   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
197   if (dmAux) {
198     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
199     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
200   }
201   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
202   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
203   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
204   else               {cellsp = sp;}
205   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
206   /* Get cell dual spaces */
207   for (f = 0; f < Nf; ++f) {
208     PetscObject  obj;
209     PetscClassId id;
210 
211     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
212     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
213     if (id == PETSCFE_CLASSID) {
214       PetscFE fe = (PetscFE) obj;
215 
216       hasFE   = PETSC_TRUE;
217       isFE[f] = PETSC_TRUE;
218       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
219     } else if (id == PETSCFV_CLASSID) {
220       PetscFV fv = (PetscFV) obj;
221 
222       hasFV   = PETSC_TRUE;
223       isFE[f] = PETSC_FALSE;
224       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
225     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
226   }
227   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
228     PetscFE    fem;
229     PetscReal *points;
230     PetscInt   numPoints, spDim, d;
231 
232     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
233     numPoints = 0;
234     for (f = 0; f < Nf; ++f) {
235       ierr = PetscDualSpaceGetDimension(cellsp[f], &spDim);CHKERRQ(ierr);
236       for (d = 0; d < spDim; ++d) {
237         if (funcs[f]) {
238           PetscQuadrature quad;
239           PetscInt        Nq;
240 
241           ierr = PetscDualSpaceGetFunctional(cellsp[f], d, &quad);CHKERRQ(ierr);
242           ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
243           numPoints += Nq;
244         }
245       }
246     }
247     ierr = PetscMalloc1(numPoints*dim, &points);CHKERRQ(ierr);
248     numPoints = 0;
249     for (f = 0; f < Nf; ++f) {
250       ierr = PetscDualSpaceGetDimension(cellsp[f], &spDim);CHKERRQ(ierr);
251       for (d = 0; d < spDim; ++d) {
252         if (funcs[f]) {
253           PetscQuadrature  quad;
254           const PetscReal *qpoints;
255           PetscInt         Nq, q;
256 
257           ierr = PetscDualSpaceGetFunctional(cellsp[f], d, &quad);CHKERRQ(ierr);
258           ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, NULL);CHKERRQ(ierr);
259           for (q = 0; q < Nq*dim; ++q) points[numPoints*dim+q] = qpoints[q];
260           numPoints += Nq;
261         }
262       }
263     }
264     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
265     for (f = 0; f < Nf; ++f) {
266       if (!isFE[f]) continue;
267       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
268       ierr = PetscFEGetTabulation(fem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
269     }
270     for (f = 0; f < NfAux; ++f) {
271       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
272       ierr = PetscFEGetTabulation(fem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
273     }
274     ierr = PetscFree(points);CHKERRQ(ierr);
275   }
276   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
277   for (h = minHeight; h <= maxHeight; h++) {
278     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
279     PetscScalar *values;
280     PetscBool   *fieldActive;
281     PetscInt     pStart, pEnd, p, spDim, totDim, numValues;
282 
283     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
284     if (!h) {
285       PetscInt cEndInterior;
286 
287       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
288       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
289     }
290     if (pEnd <= pStart) continue;
291     /* Compute totDim, the number of dofs in the closure of a point at this height */
292     totDim = 0;
293     for (f = 0; f < Nf; ++f) {
294       if (!effectiveHeight) {
295         sp[f] = cellsp[f];
296       } else {
297         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
298         if (!sp[f]) continue;
299       }
300       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
301       totDim += spDim;
302     }
303     ierr = DMPlexVecGetClosure(dm, section, localX, pStart, &numValues, NULL);CHKERRQ(ierr);
304     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim);
305     if (!totDim) continue;
306     /* Loop over points at this height */
307     ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
308     ierr = DMGetWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
309     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
310     if (label) {
311       PetscInt i;
312 
313       for (i = 0; i < numIds; ++i) {
314         IS              pointIS;
315         const PetscInt *points;
316         PetscInt        n;
317 
318         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
319         if (!pointIS) continue; /* No points with that id on this process */
320         ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
321         ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
322         for (p = 0; p < n; ++p) {
323           const PetscInt  point = points[p];
324 
325           if ((point < pStart) || (point >= pEnd)) continue;
326           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
327           ierr = DMProjectPoint_Private(dm, dmAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
328           if (ierr) {
329             PetscErrorCode ierr2;
330             ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2);
331             ierr2 = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr2);
332             CHKERRQ(ierr);
333           }
334           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
335         }
336         ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
337         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
338       }
339     } else {
340       for (p = pStart; p < pEnd; ++p) {
341         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
342         ierr = DMProjectPoint_Private(dm, dmAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
343         if (ierr) {
344           PetscErrorCode ierr2;
345           ierr2 = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr2);
346           ierr2 = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr2);
347           CHKERRQ(ierr);
348         }
349         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, values, mode);CHKERRQ(ierr);
350       }
351     }
352     ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
353     ierr = DMRestoreWorkArray(dm, Nf, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
354   }
355   /* Cleanup */
356   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
357     PetscFE fem;
358 
359     for (f = 0; f < Nf; ++f) {
360       if (!isFE[f]) continue;
361       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
362       ierr = PetscFERestoreTabulation(fem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
363     }
364     for (f = 0; f < NfAux; ++f) {
365       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
366       ierr = PetscFERestoreTabulation(fem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
367     }
368     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
369   }
370   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
371   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
372   PetscFunctionReturn(0);
373 }
374 
375 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
376 {
377   PetscErrorCode ierr;
378 
379   PetscFunctionBegin;
380   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
381   PetscFunctionReturn(0);
382 }
383 
384 PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
385 {
386   PetscErrorCode ierr;
387 
388   PetscFunctionBegin;
389   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
390   PetscFunctionReturn(0);
391 }
392 
393 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
394                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
395                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
396                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
397                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
398                                         InsertMode mode, Vec localX)
399 {
400   PetscErrorCode ierr;
401 
402   PetscFunctionBegin;
403   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
404   PetscFunctionReturn(0);
405 }
406 
407 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], Vec localU,
408                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
409                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
410                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
411                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
412                                              InsertMode mode, Vec localX)
413 {
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin;
417   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
418   PetscFunctionReturn(0);
419 }
420