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