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