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