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