xref: /petsc/src/dm/impls/plex/plexproject.c (revision 702c4f2fb2b40d2bce4bfb3c8519157ae61fced8)
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, PetscDS prob, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
6                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
7                                                   PetscScalar values[])
8 {
9   PetscInt       coordDim, Nf, *Nc, f, totDim, spDim, d, v, tp;
10   PetscBool      isAffine, transform;
11   PetscErrorCode ierr;
12 
13   PetscFunctionBeginHot;
14   ierr = DMGetCoordinateDim(dm,&coordDim);CHKERRQ(ierr);
15   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
16   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
17   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
18   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
19   /* Get values for closure */
20   isAffine = fegeom->isAffine;
21   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
22     void * const ctx = ctxs ? ctxs[f] : NULL;
23 
24     if (!sp[f]) continue;
25     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
26     if (funcs[f]) {
27       if (isFE[f]) {
28         PetscQuadrature   allPoints;
29         PetscInt          q, dim, numPoints;
30         const PetscReal   *points;
31         PetscScalar       *pointEval;
32         PetscReal         *x;
33         DM                rdm;
34 
35         ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr);
36         ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
37         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
38         ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
39         ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
40         for (q = 0; q < numPoints; q++, tp++) {
41           const PetscReal *v0;
42 
43           if (isAffine) {
44             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
45             v0 = x;
46           } else {
47             v0 = &fegeom->v[tp*coordDim];
48           }
49           if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;}
50           ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
51         }
52         /* Transform point evaluations pointEval[q,c] */
53         ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr);
54         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
55         ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
56         ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
57         v += spDim;
58       } else {
59         for (d = 0; d < spDim; ++d, ++v) {
60           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
61         }
62       }
63     } else {
64       for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;}
65     }
66   }
67   PetscFunctionReturn(0);
68 }
69 
70 static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
71                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
72                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
73                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
74                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
75                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
76                                                    PetscScalar values[])
77 {
78   PetscSection       section, sectionAux = NULL;
79   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
80   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
81   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
82   const PetscScalar *constants;
83   PetscReal         *x;
84   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
85   PetscFEGeom        fegeom;
86   const PetscInt     dE = cgeom->dimEmbed;
87   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
88   PetscBool          isAffine, transform;
89   PetscErrorCode     ierr;
90 
91   PetscFunctionBeginHot;
92   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
93   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
94   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
95   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
96   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
97   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
98   ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
99   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
100   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
101   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
102   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
103   if (dmAux) {
104     PetscInt subp;
105 
106     ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, p, &subp);CHKERRQ(ierr);
107     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
108     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
109     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
110     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
111     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
112     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
113     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
114     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
115     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
116   }
117   /* Get values for closure */
118   isAffine = cgeom->isAffine;
119   if (isAffine) {
120     fegeom.v    = x;
121     fegeom.xi   = cgeom->xi;
122     fegeom.J    = cgeom->J;
123     fegeom.invJ = cgeom->invJ;
124     fegeom.detJ = cgeom->detJ;
125   }
126   for (f = 0, v = 0; f < Nf; ++f) {
127     PetscQuadrature   allPoints;
128     PetscInt          q, dim, numPoints;
129     const PetscReal   *points;
130     PetscScalar       *pointEval;
131     DM                dm;
132 
133     if (!sp[f]) continue;
134     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
135     if (!funcs[f]) {
136       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
137       continue;
138     }
139     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
140     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
141     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
142     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
143     for (q = 0; q < numPoints; ++q, ++tp) {
144       if (isAffine) {
145         CoordinatesRefToReal(dE, dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
146       } else {
147         fegeom.v    = &cgeom->v[tp*dE];
148         fegeom.J    = &cgeom->J[tp*dE*dE];
149         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
150         fegeom.detJ = &cgeom->detJ[tp];
151       }
152       ierr = PetscFEEvaluateFieldJets_Internal(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
153       if (probAux) {ierr = PetscFEEvaluateFieldJets_Internal(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
154       if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);}
155       (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, numConstants, constants, &pointEval[Nc[f]*q]);
156     }
157     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
158     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
159     v += spDim;
160   }
161   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
162   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
163   PetscFunctionReturn(0);
164 }
165 
166 static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
167                                                      PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
168                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
169                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
170                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
171                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
172                                                      PetscScalar values[])
173 {
174   PetscSection       section, sectionAux = NULL;
175   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
176   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
177   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
178   const PetscScalar *constants;
179   PetscReal         *x;
180   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
181   PetscFEGeom        fegeom, cgeom;
182   const PetscInt     dE = fgeom->dimEmbed;
183   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
184   PetscBool          isAffine;
185   PetscErrorCode     ierr;
186 
187   PetscFunctionBeginHot;
188   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
189   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
190   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
191   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
192   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
193   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
194   ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
195   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
196   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
197   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
198   if (dmAux) {
199     DMLabel  spmap;
200     PetscInt subp = p;
201 
202     /* If dm is a submesh, do not get subpoint */
203     ierr = DMPlexGetSubpointMap(dm, &spmap);CHKERRQ(ierr);
204     if (!spmap) {ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr);}
205     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
206     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
207     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
208     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
209     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
210     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
211     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
212     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
213     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
214   }
215   /* Get values for closure */
216   isAffine = fgeom->isAffine;
217   fegeom.n  = 0;
218   fegeom.J  = 0;
219   fegeom.v  = 0;
220   fegeom.xi = 0;
221   if (isAffine) {
222     fegeom.v    = x;
223     fegeom.xi   = fgeom->xi;
224     fegeom.J    = fgeom->J;
225     fegeom.invJ = fgeom->invJ;
226     fegeom.detJ = fgeom->detJ;
227     fegeom.n    = fgeom->n;
228 
229     cgeom.J     = fgeom->suppJ[0];
230     cgeom.invJ  = fgeom->suppInvJ[0];
231     cgeom.detJ  = fgeom->suppDetJ[0];
232   }
233   for (f = 0, v = 0; f < Nf; ++f) {
234     PetscQuadrature   allPoints;
235     PetscInt          q, dim, numPoints;
236     const PetscReal   *points;
237     PetscScalar       *pointEval;
238     DM                dm;
239 
240     if (!sp[f]) continue;
241     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
242     if (!funcs[f]) {
243       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
244       continue;
245     }
246     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
247     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
248     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
249     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
250     for (q = 0; q < numPoints; ++q, ++tp) {
251       if (isAffine) {
252         CoordinatesRefToReal(dE, dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
253       } else {
254         fegeom.v    = &fgeom->v[tp*dE];
255         fegeom.J    = &fgeom->J[tp*dE*dE];
256         fegeom.invJ = &fgeom->invJ[tp*dE*dE];
257         fegeom.detJ = &fgeom->detJ[tp];
258         fegeom.n    = &fgeom->n[tp*dE];
259 
260         cgeom.J     = &fgeom->suppJ[0][tp*dE*dE];
261         cgeom.invJ  = &fgeom->suppInvJ[0][tp*dE*dE];
262         cgeom.detJ  = &fgeom->suppDetJ[0][tp];
263       }
264       ierr = PetscFEEvaluateFieldJets_Internal(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
265       if (probAux) {ierr = PetscFEEvaluateFieldJets_Internal(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
266       (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, fegeom.n, numConstants, constants, &pointEval[Nc[f]*q]);
267     }
268     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
269     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
270     v += spDim;
271   }
272   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
273   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
274   PetscFunctionReturn(0);
275 }
276 
277 static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS prob, PetscFEGeom *fegeom, DM dmAux, PetscDS probAux, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
278                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
279                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
280                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
281 {
282   PetscFVCellGeom fvgeom;
283   PetscInt        dim, dimEmbed;
284   PetscErrorCode  ierr;
285 
286   PetscFunctionBeginHot;
287   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
288   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
289   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
290   switch (type) {
291   case DM_BC_ESSENTIAL:
292   case DM_BC_NATURAL:
293     ierr = DMProjectPoint_Func_Private(dm, prob, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
294   case DM_BC_ESSENTIAL_FIELD:
295   case DM_BC_NATURAL_FIELD:
296     ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
297                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
298                                         (void (**)(PetscInt, PetscInt, PetscInt,
299                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
300                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
301                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
302   case DM_BC_ESSENTIAL_BD_FIELD:
303     ierr = DMProjectPoint_BdField_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
304                                           basisTab, basisDerTab, basisTabAux, basisDerTabAux,
305                                           (void (**)(PetscInt, PetscInt, PetscInt,
306                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
307                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
308                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
309   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
310   }
311   PetscFunctionReturn(0);
312 }
313 
314 static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
315 {
316   PetscReal      *points;
317   PetscInt       f, numPoints;
318   PetscErrorCode ierr;
319 
320   PetscFunctionBegin;
321   numPoints = 0;
322   for (f = 0; f < Nf; ++f) {
323     if (funcs[f]) {
324       PetscQuadrature fAllPoints;
325       PetscInt        fNumPoints;
326 
327       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
328       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
329       numPoints += fNumPoints;
330     }
331   }
332   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
333   numPoints = 0;
334   for (f = 0; f < Nf; ++f) {
335     PetscInt spDim;
336 
337     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
338     if (funcs[f]) {
339       PetscQuadrature fAllPoints;
340       PetscInt        qdim, fNumPoints, q;
341       const PetscReal *fPoints;
342 
343       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
344       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
345       if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %D for dual basis does not match input dimension %D", qdim, dim);
346       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
347       numPoints += fNumPoints;
348     }
349   }
350   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
351   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
352   PetscFunctionReturn(0);
353 }
354 
355 static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *prob)
356 {
357   DMLabel        depthLabel;
358   PetscInt       dim, cdepth, ls = -1, i;
359   PetscErrorCode ierr;
360 
361   PetscFunctionBegin;
362   if (lStart) *lStart = -1;
363   if (!label) PetscFunctionReturn(0);
364   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
365   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
366   cdepth = dim - height;
367   for (i = 0; i < numIds; ++i) {
368     IS              pointIS;
369     const PetscInt *points;
370     PetscInt        pdepth;
371 
372     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
373     if (!pointIS) continue; /* No points with that id on this process */
374     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
375     ierr = DMLabelGetValue(depthLabel, points[0], &pdepth);CHKERRQ(ierr);
376     if (pdepth == cdepth) {
377       ls = points[0];
378       if (prob) {ierr = DMGetCellDS(dm, ls, prob);CHKERRQ(ierr);}
379     }
380     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
381     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
382     if (ls >= 0) break;
383   }
384   if (lStart) *lStart = ls;
385   PetscFunctionReturn(0);
386 }
387 
388 /*
389   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
390 
391   There are several different scenarios:
392 
393   1) Volumetric mesh with volumetric auxiliary data
394 
395      Here minHeight=0 since we loop over cells.
396 
397   2) Boundary mesh with boundary auxiliary data
398 
399      Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation.
400 
401   3) Volumetric mesh with boundary auxiliary data
402 
403      Here minHeight=1 and auxbd=PETSC_TRUE since we loop over faces and use data only supported on those faces. This is common when imposing Dirichlet boundary conditions.
404 
405   The maxHeight is used to support enforcement of constraints in DMForest.
406 
407   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
408 
409   If we are using an input field (DM_BC_ESSENTIAL_FIELD or DM_BC_NATURAL_FIELD), we need to evaluate it at all the quadrature points of the dual basis functionals.
410     - We use effectiveHeight to mean the height above our incoming DS. For example, if the DS is for a submesh then the effective height is zero, whereas if the DS
411       is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When th effective height is nonzero, we need to extract
412       dual spaces for the boundary from our input spaces.
413     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
414 
415   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
416 
417   If we have a label, we iterate over those points. This will probably break the maxHeight functionality since we do not check the height of those points.
418 */
419 static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
420                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
421                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
422                                                   InsertMode mode, Vec localX)
423 {
424   DM              dmAux = NULL, tdm;
425   PetscDS         prob = NULL, probAux = NULL;
426   Vec             localA = NULL, tv;
427   PetscSection    section;
428   PetscDualSpace *sp, *cellsp;
429   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
430   PetscInt       *Nc;
431   PetscInt        dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f;
432   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, transform;
433   DMField         coordField;
434   DMLabel         depthLabel;
435   PetscQuadrature allPoints = NULL;
436   PetscErrorCode  ierr;
437 
438   PetscFunctionBegin;
439   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
440   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
441   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
442   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
443   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
444   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
445   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
446   /* Auxiliary information can only be used with interpolation of field functions */
447   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
448     if (dmAux && !localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
449     if (!minHeight && dmAux) {
450       DMLabel spmap;
451 
452       /* If dmAux is a surface, then force the projection to take place over a surface */
453       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
454       if (spmap) {
455         ierr = DMPlexGetVTKCellHeight(dmAux, &minHeight);CHKERRQ(ierr);
456         auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
457       }
458     }
459   }
460   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
461   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
462   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
463   maxHeight = PetscMax(maxHeight, minHeight);
464   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
465   ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, 0, NULL, &prob);CHKERRQ(ierr);
466   if (!prob) {ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);}
467   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
468   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
469   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
470   if (dmAux) {
471     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
472     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
473   }
474   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
475   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
476   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
477   else               {cellsp = sp;}
478   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
479   /* Get cell dual spaces */
480   for (f = 0; f < Nf; ++f) {
481     PetscObject  obj;
482     PetscClassId id;
483 
484     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
485     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
486     if (id == PETSCFE_CLASSID) {
487       PetscFE fe = (PetscFE) obj;
488 
489       hasFE   = PETSC_TRUE;
490       isFE[f] = PETSC_TRUE;
491       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
492     } else if (id == PETSCFV_CLASSID) {
493       PetscFV fv = (PetscFV) obj;
494 
495       hasFV   = PETSC_TRUE;
496       isFE[f] = PETSC_FALSE;
497       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
498     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
499   }
500   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
501   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
502     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
503     PetscFE          fem, subfem;
504     const PetscReal *points;
505     PetscInt         numPoints;
506 
507     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
508     for (f = 0; f < Nf; ++f) {
509       if (!effectiveHeight) {sp[f] = cellsp[f];}
510       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
511     }
512     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
513     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
514     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
515     for (f = 0; f < Nf; ++f) {
516       if (!isFE[f]) continue;
517       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
518       if (!effectiveHeight) {subfem = fem;}
519       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
520       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
521     }
522     for (f = 0; f < NfAux; ++f) {
523       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
524       if (!effectiveHeight || auxBd) {subfem = fem;}
525       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
526       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
527     }
528   }
529   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
530   for (h = minHeight; h <= maxHeight; h++) {
531     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
532     PetscDS      probEff         = prob;
533     PetscScalar *values;
534     PetscBool   *fieldActive;
535     PetscInt     maxDegree;
536     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
537     IS           heightIS;
538 
539     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
540     ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
541     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
542     if (!h) {
543       PetscInt cEndInterior;
544 
545       ierr = DMPlexGetInteriorCellStratum(dm, NULL, &cEndInterior);CHKERRQ(ierr);
546       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
547     }
548     if (pEnd <= pStart) {
549       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
550       continue;
551     }
552     /* Compute totDim, the number of dofs in the closure of a point at this height */
553     totDim = 0;
554     for (f = 0; f < Nf; ++f) {
555       if (!effectiveHeight) {
556         sp[f] = cellsp[f];
557       } else {
558         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
559         if (!sp[f]) continue;
560       }
561       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
562       totDim += spDim;
563     }
564     ierr = DMPlexVecGetClosure(dm, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr);
565     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim);
566     if (!totDim) {
567       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
568       continue;
569     }
570     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);}
571     /* Loop over points at this height */
572     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
573     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
574     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
575     if (label) {
576       PetscInt i;
577 
578       for (i = 0; i < numIds; ++i) {
579         IS              pointIS, isectIS;
580         const PetscInt *points;
581         PetscInt        n;
582         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
583         PetscQuadrature quad = NULL;
584 
585         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
586         if (!pointIS) continue; /* No points with that id on this process */
587         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
588         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
589         if (!isectIS) continue;
590         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
591         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
592         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
593         if (maxDegree <= 1) {
594           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
595         }
596         if (!quad) {
597           if (!h && allPoints) {
598             quad = allPoints;
599             allPoints = NULL;
600           } else {
601             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
602           }
603         }
604         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
605         for (p = 0; p < n; ++p) {
606           const PetscInt  point = points[p];
607 
608           ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
609           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
610           ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
611           if (ierr) {
612             PetscErrorCode ierr2;
613             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
614             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
615             CHKERRQ(ierr);
616           }
617           if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
618           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
619         }
620         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
621         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
622         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
623         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
624         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
625       }
626     } else {
627       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
628       PetscQuadrature quad = NULL;
629       IS              pointIS;
630 
631       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
632       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
633       if (maxDegree <= 1) {
634         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
635       }
636       if (!quad) {
637         if (!h && allPoints) {
638           quad = allPoints;
639           allPoints = NULL;
640         } else {
641           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
642         }
643       }
644       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
645       for (p = pStart; p < pEnd; ++p) {
646         ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
647         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
648         ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
649         if (ierr) {
650           PetscErrorCode ierr2;
651           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
652           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
653           CHKERRQ(ierr);
654         }
655         if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
656         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
657       }
658       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
659       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
660       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
661       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
662     }
663     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
664     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
665     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
666   }
667   /* Cleanup */
668   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
669     PetscInt effectiveHeight = auxBd ? minHeight : 0;
670     PetscFE  fem, subfem;
671 
672     for (f = 0; f < Nf; ++f) {
673       if (!isFE[f]) continue;
674       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
675       if (!effectiveHeight) {subfem = fem;}
676       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
677       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
678     }
679     for (f = 0; f < NfAux; ++f) {
680       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
681       if (!effectiveHeight || auxBd) {subfem = fem;}
682       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
683       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
684     }
685     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
686   }
687   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
688   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
689   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
690   PetscFunctionReturn(0);
691 }
692 
693 PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
694 {
695   PetscErrorCode ierr;
696 
697   PetscFunctionBegin;
698   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
699   PetscFunctionReturn(0);
700 }
701 
702 PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
703 {
704   PetscErrorCode ierr;
705 
706   PetscFunctionBegin;
707   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
708   PetscFunctionReturn(0);
709 }
710 
711 PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
712                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
713                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
714                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
715                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
716                                         InsertMode mode, Vec localX)
717 {
718   PetscErrorCode ierr;
719 
720   PetscFunctionBegin;
721   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
722   PetscFunctionReturn(0);
723 }
724 
725 PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
726                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
727                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
728                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
729                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
730                                              InsertMode mode, Vec localX)
731 {
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
736   PetscFunctionReturn(0);
737 }
738