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