Actual source code: plexproject.c

  1: #include <petsc/private/dmpleximpl.h>

  3: #include <petsc/private/petscfeimpl.h>

  5: /*@
  6:   DMPlexGetActivePoint - Get the point on which projection is currently working

  8:   Not Collective

 10:   Input Parameter:
 11: . dm - the `DM`

 13:   Output Parameter:
 14: . point - The mesh point involved in the current projection

 16:   Level: developer

 18: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`
 19: @*/
 20: PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point)
 21: {
 22:   PetscFunctionBeginHot;
 23:   *point = ((DM_Plex *)dm->data)->activePoint;
 24:   PetscFunctionReturn(PETSC_SUCCESS);
 25: }

 27: /*@
 28:   DMPlexSetActivePoint - Set the point on which projection is currently working

 30:   Not Collective

 32:   Input Parameters:
 33: + dm    - the `DM`
 34: - point - The mesh point involved in the current projection

 36:   Level: developer

 38: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetActivePoint()`
 39: @*/
 40: PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point)
 41: {
 42:   PetscFunctionBeginHot;
 43:   ((DM_Plex *)dm->data)->activePoint = point;
 44:   PetscFunctionReturn(PETSC_SUCCESS);
 45: }

 47: /*
 48:   DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point

 50:   Input Parameters:
 51: + dm     - The output `DM`
 52: . ds     - The output `DS`
 53: . dmIn   - The input `DM`
 54: . dsIn   - The input `DS`
 55: . time   - The time for this evaluation
 56: . fegeom - The FE geometry for this point
 57: . fvgeom - The FV geometry for this point
 58: . isFE   - Flag indicating whether each output field has an FE discretization
 59: . sp     - The output `PetscDualSpace` for each field
 60: . funcs  - The evaluation function for each field
 61: - ctxs   - The user context for each field

 63:   Output Parameter:
 64: . values - The value for each dual basis vector in the output dual space

 66:   Level: developer

 68: .seealso:[](ch_unstructured), `DM`, `DMPLEX`, `PetscDS`, `PetscFEGeom`, `PetscFVCellGeom`, `PetscDualSpace`
 69: */
 70: static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, PetscScalar values[])
 71: {
 72:   PetscInt  coordDim, Nf, *Nc, f, spDim, d, v, tp;
 73:   PetscBool isAffine, isCohesive, transform;

 75:   PetscFunctionBeginHot;
 76:   PetscCall(DMGetCoordinateDim(dmIn, &coordDim));
 77:   PetscCall(DMHasBasisTransform(dmIn, &transform));
 78:   PetscCall(PetscDSGetNumFields(ds, &Nf));
 79:   PetscCall(PetscDSGetComponents(ds, &Nc));
 80:   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
 81:   /* Get values for closure */
 82:   isAffine = fegeom->isAffine;
 83:   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
 84:     void *const ctx = ctxs ? ctxs[f] : NULL;
 85:     PetscBool   cohesive;

 87:     if (!sp[f]) continue;
 88:     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
 89:     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
 90:     if (funcs[f]) {
 91:       if (isFE[f]) {
 92:         PetscQuadrature  allPoints;
 93:         PetscInt         q, dim, numPoints;
 94:         const PetscReal *points;
 95:         PetscScalar     *pointEval;
 96:         PetscReal       *x;
 97:         DM               rdm;

 99:         PetscCall(PetscDualSpaceGetDM(sp[f], &rdm));
100:         PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
101:         PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
102:         PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
103:         PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x));
104:         PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f]));
105:         for (q = 0; q < numPoints; q++, tp++) {
106:           const PetscReal *v0;

108:           if (isAffine) {
109:             const PetscReal *refpoint    = &points[q * dim];
110:             PetscReal        injpoint[3] = {0., 0., 0.};

112:             if (dim != fegeom->dim) {
113:               if (isCohesive) {
114:                 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
115:                 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
116:                 refpoint = injpoint;
117:               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim);
118:             }
119:             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x);
120:             v0 = x;
121:           } else {
122:             v0 = &fegeom->v[tp * coordDim];
123:           }
124:           if (transform) {
125:             PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx));
126:             v0 = x;
127:           }
128:           PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx));
129:         }
130:         /* Transform point evaluations pointEval[q,c] */
131:         PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
132:         PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
133:         PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x));
134:         PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
135:         v += spDim;
136:         if (isCohesive && !cohesive) {
137:           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
138:         }
139:       } else {
140:         for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
141:       }
142:     } else {
143:       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
144:       if (isCohesive && !cohesive) {
145:         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
146:       }
147:     }
148:   }
149:   PetscFunctionReturn(PETSC_SUCCESS);
150: }

152: /*
153:   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point

155:   Input Parameters:
156: + dm             - The output DM
157: . ds             - The output DS
158: . dmIn           - The input DM
159: . dsIn           - The input DS
160: . dmAux          - The auxiliary DM, which is always for the input space
161: . dsAux          - The auxiliary DS, which is always for the input space
162: . time           - The time for this evaluation
163: . localU         - The local solution
164: . localA         - The local auziliary fields
165: . cgeom          - The FE geometry for this point
166: . sp             - The output PetscDualSpace for each field
167: . p              - The point in the output DM
168: . T              - Input basis and derivatives for each field tabulated on the quadrature points
169: . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
170: . funcs          - The evaluation function for each field
171: - ctxs           - The user context for each field

173:   Output Parameter:
174: . values         - The value for each dual basis vector in the output dual space

176:   Level: developer

178:   Note:
179:   Not supported for FV

181: .seealso: `DMProjectPoint_Field_Private()`
182: */
183: static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[])
184: {
185:   PetscSection       section, sectionAux = NULL;
186:   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
187:   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
188:   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
189:   const PetscScalar *constants;
190:   PetscReal         *x;
191:   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
192:   PetscFEGeom        fegeom;
193:   const PetscInt     dE = cgeom->dimEmbed, *cone, *ornt;
194:   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
195:   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
196:   DMPolytopeType     qct;

198:   PetscFunctionBeginHot;
199:   PetscCall(PetscDSGetNumFields(ds, &Nf));
200:   PetscCall(PetscDSGetComponents(ds, &Nc));
201:   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
202:   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
203:   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
204:   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
205:   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
206:   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
207:   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
208:   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
209:   PetscCall(DMHasBasisTransform(dmIn, &transform));
210:   PetscCall(DMGetLocalSection(dmIn, &section));
211:   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
212:   // Get cohesive cell hanging off face
213:   if (isCohesiveIn) {
214:     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
215:     if ((qct != DM_POLYTOPE_POINT_PRISM_TENSOR) && (qct != DM_POLYTOPE_SEG_PRISM_TENSOR) && (qct != DM_POLYTOPE_TRI_PRISM_TENSOR) && (qct != DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
216:       DMPolytopeType  ct;
217:       const PetscInt *support;
218:       PetscInt        Ns, s;

220:       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
221:       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
222:       for (s = 0; s < Ns; ++s) {
223:         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
224:         if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
225:           inp = support[s];
226:           break;
227:         }
228:       }
229:       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
230:       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
231:       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
232:       face[0] = 0;
233:       face[1] = 0;
234:     }
235:   }
236:   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
237:   if (dmAux) {
238:     PetscInt subp;

240:     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
241:     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
242:     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
243:     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
244:     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
245:     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
246:     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
247:   }
248:   /* Get values for closure */
249:   isAffine        = cgeom->isAffine;
250:   fegeom.dim      = cgeom->dim;
251:   fegeom.dimEmbed = cgeom->dimEmbed;
252:   if (isAffine) {
253:     fegeom.v    = x;
254:     fegeom.xi   = cgeom->xi;
255:     fegeom.J    = cgeom->J;
256:     fegeom.invJ = cgeom->invJ;
257:     fegeom.detJ = cgeom->detJ;
258:   }
259:   for (f = 0, v = 0; f < Nf; ++f) {
260:     PetscQuadrature  allPoints;
261:     PetscInt         q, dim, numPoints;
262:     const PetscReal *points;
263:     PetscScalar     *pointEval;
264:     PetscBool        cohesive;
265:     DM               dm;

267:     if (!sp[f]) continue;
268:     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
269:     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
270:     if (!funcs[f]) {
271:       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
272:       if (isCohesive && !cohesive) {
273:         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
274:       }
275:       continue;
276:     }
277:     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
278:     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
279:     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
280:     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
281:     for (q = 0; q < numPoints; ++q, ++tp) {
282:       PetscInt qpt[2];

284:       if (isCohesiveIn) {
285:         PetscCall(PetscDSPermuteQuadPoint(dsIn, ornt[0], f, q, &qpt[0]));
286:         PetscCall(PetscDSPermuteQuadPoint(dsIn, DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0), f, q, &qpt[1]));
287:       }
288:       if (isAffine) {
289:         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x);
290:       } else {
291:         fegeom.v    = &cgeom->v[tp * dE];
292:         fegeom.J    = &cgeom->J[tp * dE * dE];
293:         fegeom.invJ = &cgeom->invJ[tp * dE * dE];
294:         fegeom.detJ = &cgeom->detJ[tp];
295:       }
296:       if (coefficients) {
297:         if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
298:         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
299:       }
300:       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
301:       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
302:       (*funcs[f])(dE, NfIn, 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]);
303:     }
304:     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
305:     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
306:     v += spDim;
307:     /* TODO: For now, set both sides equal, but this should use info from other support cell */
308:     if (isCohesive && !cohesive) {
309:       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
310:     }
311:   }
312:   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
313:   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
314:   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
315:   PetscFunctionReturn(PETSC_SUCCESS);
316: }

318: static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[])
319: {
320:   PetscSection       section, sectionAux = NULL;
321:   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
322:   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
323:   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
324:   const PetscScalar *constants;
325:   PetscReal         *x;
326:   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
327:   PetscFEGeom        fegeom, cgeom;
328:   const PetscInt     dE = fgeom->dimEmbed;
329:   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
330:   PetscBool          isAffine;

332:   PetscFunctionBeginHot;
333:   PetscCheck(dm == dmIn, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
334:   PetscCall(PetscDSGetNumFields(ds, &Nf));
335:   PetscCall(PetscDSGetComponents(ds, &Nc));
336:   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
337:   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
338:   PetscCall(PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x));
339:   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
340:   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
341:   PetscCall(DMGetLocalSection(dm, &section));
342:   PetscCall(DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients));
343:   if (dmAux) {
344:     PetscInt subp;

346:     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
347:     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
348:     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
349:     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
350:     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
351:     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
352:     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
353:   }
354:   /* Get values for closure */
355:   isAffine       = fgeom->isAffine;
356:   fegeom.n       = NULL;
357:   fegeom.J       = NULL;
358:   fegeom.v       = NULL;
359:   fegeom.xi      = NULL;
360:   cgeom.dim      = fgeom->dim;
361:   cgeom.dimEmbed = fgeom->dimEmbed;
362:   if (isAffine) {
363:     fegeom.v    = x;
364:     fegeom.xi   = fgeom->xi;
365:     fegeom.J    = fgeom->J;
366:     fegeom.invJ = fgeom->invJ;
367:     fegeom.detJ = fgeom->detJ;
368:     fegeom.n    = fgeom->n;

370:     cgeom.J    = fgeom->suppJ[0];
371:     cgeom.invJ = fgeom->suppInvJ[0];
372:     cgeom.detJ = fgeom->suppDetJ[0];
373:   }
374:   for (f = 0, v = 0; f < Nf; ++f) {
375:     PetscQuadrature  allPoints;
376:     PetscInt         q, dim, numPoints;
377:     const PetscReal *points;
378:     PetscScalar     *pointEval;
379:     DM               dm;

381:     if (!sp[f]) continue;
382:     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
383:     if (!funcs[f]) {
384:       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
385:       continue;
386:     }
387:     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
388:     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
389:     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
390:     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
391:     for (q = 0; q < numPoints; ++q, ++tp) {
392:       if (isAffine) {
393:         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x);
394:       } else {
395:         fegeom.v    = &fgeom->v[tp * dE];
396:         fegeom.J    = &fgeom->J[tp * dE * dE];
397:         fegeom.invJ = &fgeom->invJ[tp * dE * dE];
398:         fegeom.detJ = &fgeom->detJ[tp];
399:         fegeom.n    = &fgeom->n[tp * dE];

401:         cgeom.J    = &fgeom->suppJ[0][tp * dE * dE];
402:         cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE];
403:         cgeom.detJ = &fgeom->suppDetJ[0][tp];
404:       }
405:       /* TODO We should use cgeom here, instead of fegeom, however the geometry coming in through fgeom does not have the support cell geometry */
406:       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
407:       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
408:       (*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]);
409:     }
410:     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
411:     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
412:     v += spDim;
413:   }
414:   PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients));
415:   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
416:   PetscFunctionReturn(PETSC_SUCCESS);
417: }

419: static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscFEGeom *fegeom, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[], PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
420: {
421:   PetscFVCellGeom fvgeom;
422:   PetscInt        dim, dimEmbed;

424:   PetscFunctionBeginHot;
425:   PetscCall(DMGetDimension(dm, &dim));
426:   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
427:   if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
428:   switch (type) {
429:   case DM_BC_ESSENTIAL:
430:   case DM_BC_NATURAL:
431:     PetscCall(DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode(**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))funcs, ctxs, values));
432:     break;
433:   case DM_BC_ESSENTIAL_FIELD:
434:   case DM_BC_NATURAL_FIELD:
435:     PetscCall(DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
436:     break;
437:   case DM_BC_ESSENTIAL_BD_FIELD:
438:     PetscCall(DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
439:     break;
440:   default:
441:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type);
442:   }
443:   PetscFunctionReturn(PETSC_SUCCESS);
444: }

446: static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
447: {
448:   PetscReal *points;
449:   PetscInt   f, numPoints;

451:   PetscFunctionBegin;
452:   if (!dim) {
453:     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
454:     PetscFunctionReturn(PETSC_SUCCESS);
455:   }
456:   numPoints = 0;
457:   for (f = 0; f < Nf; ++f) {
458:     if (funcs[f]) {
459:       PetscQuadrature fAllPoints;
460:       PetscInt        fNumPoints;

462:       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
463:       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
464:       numPoints += fNumPoints;
465:     }
466:   }
467:   PetscCall(PetscMalloc1(dim * numPoints, &points));
468:   numPoints = 0;
469:   for (f = 0; f < Nf; ++f) {
470:     if (funcs[f]) {
471:       PetscQuadrature  fAllPoints;
472:       PetscInt         qdim, fNumPoints, q;
473:       const PetscReal *fPoints;

475:       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
476:       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
477:       PetscCheck(qdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %" PetscInt_FMT " for dual basis does not match input dimension %" PetscInt_FMT, qdim, dim);
478:       for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q];
479:       numPoints += fNumPoints;
480:     }
481:   }
482:   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
483:   PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL));
484:   PetscFunctionReturn(PETSC_SUCCESS);
485: }

487: /*@C
488:   DMGetFirstLabeledPoint - Find first labeled `point` in `odm` such that the corresponding point in `dm` has the specified `height`. Return `point` and the corresponding `ds`.

490:   Input Parameters:
491: + dm     - the `DM`
492: . odm    - the enclosing `DM`
493: . label  - label for `DM` domain, or `NULL` for whole domain
494: . numIds - the number of `ids`
495: . ids    - An array of the label ids in sequence for the domain
496: - height - Height of target cells in `DMPLEX` topology

498:   Output Parameters:
499: + point - the first labeled point
500: - ds    - the ds corresponding to the first labeled point

502:   Level: developer

504: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS`
505: @*/
506: PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
507: {
508:   DM              plex;
509:   DMEnclosureType enc;
510:   PetscInt        ls = -1;

512:   PetscFunctionBegin;
513:   if (point) *point = -1;
514:   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
515:   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
516:   PetscCall(DMConvert(dm, DMPLEX, &plex));
517:   for (PetscInt i = 0; i < numIds; ++i) {
518:     IS       labelIS;
519:     PetscInt num_points, pStart, pEnd;
520:     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
521:     if (!labelIS) continue; /* No points with that id on this process */
522:     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
523:     PetscCall(ISGetSize(labelIS, &num_points));
524:     if (num_points) {
525:       const PetscInt *points;
526:       PetscCall(ISGetIndices(labelIS, &points));
527:       for (PetscInt i = 0; i < num_points; i++) {
528:         PetscInt point;
529:         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
530:         if (pStart <= point && point < pEnd) {
531:           ls = point;
532:           if (ds) PetscCall(DMGetCellDS(dm, ls, ds, NULL));
533:           if (ls >= 0) break;
534:         }
535:       }
536:       PetscCall(ISRestoreIndices(labelIS, &points));
537:     }
538:     PetscCall(ISDestroy(&labelIS));
539:     if (ls >= 0) break;
540:   }
541:   if (point) *point = ls;
542:   PetscCall(DMDestroy(&plex));
543:   PetscFunctionReturn(PETSC_SUCCESS);
544: }

546: /*
547:   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM

549:   There are several different scenarios:

551:   1) Volumetric mesh with volumetric auxiliary data

553:      Here minHeight=0 since we loop over cells.

555:   2) Boundary mesh with boundary auxiliary data

557:      Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation.

559:   3) Volumetric mesh with boundary auxiliary data

561:      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.

563:   4) Volumetric input mesh with boundary output mesh

565:      Here we must get a subspace for the input DS

567:   The maxHeight is used to support enforcement of constraints in DMForest.

569:   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.

571:   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.
572:     - 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
573:       is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When the effective height is nonzero, we need to extract
574:       dual spaces for the boundary from our input spaces.
575:     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.

577:   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration

579:   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.
580: */
581: static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, InsertMode mode, Vec localX)
582: {
583:   DM               plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
584:   DMEnclosureType  encIn, encAux;
585:   PetscDS          ds = NULL, dsIn = NULL, dsAux = NULL;
586:   Vec              localA = NULL, tv;
587:   IS               fieldIS;
588:   PetscSection     section;
589:   PetscDualSpace  *sp, *cellsp, *spIn, *cellspIn;
590:   PetscTabulation *T = NULL, *TAux = NULL;
591:   PetscInt        *Nc;
592:   PetscInt         dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
593:   PetscBool       *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform;
594:   DMField          coordField;
595:   DMLabel          depthLabel;
596:   PetscQuadrature  allPoints = NULL;

598:   PetscFunctionBegin;
599:   if (localU) PetscCall(VecGetDM(localU, &dmIn));
600:   else dmIn = dm;
601:   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
602:   if (localA) PetscCall(VecGetDM(localA, &dmAux));
603:   else dmAux = NULL;
604:   PetscCall(DMConvert(dm, DMPLEX, &plex));
605:   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
606:   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
607:   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
608:   PetscCall(DMGetDimension(dm, &dim));
609:   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
610:   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
611:   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
612:   PetscCall(DMHasBasisTransform(dm, &transform));
613:   /* Auxiliary information can only be used with interpolation of field functions */
614:   if (dmAux) {
615:     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
616:     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) PetscCheck(localA, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector");
617:   }
618:   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL));
619:   PetscCall(DMGetCoordinateField(dm, &coordField));
620:   /**** No collective calls below this point ****/
621:   /* Determine height for iteration of all meshes */
622:   {
623:     DMPolytopeType ct, ctIn, ctAux;
624:     PetscInt       minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux;
625:     PetscInt       dim = -1, dimIn = -1, dimAux = -1;

627:     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
628:     if (pEnd > pStart) {
629:       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
630:       p = lStart < 0 ? pStart : lStart;
631:       PetscCall(DMPlexGetCellType(plex, p, &ct));
632:       dim = DMPolytopeTypeGetDim(ct);
633:       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
634:       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
635:       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
636:       dimIn = DMPolytopeTypeGetDim(ctIn);
637:       if (dmAux) {
638:         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
639:         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
640:         if (pStartAux < pEndAux) {
641:           PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
642:           dimAux = DMPolytopeTypeGetDim(ctAux);
643:         }
644:       } else dimAux = dim;
645:     } else {
646:       PetscCall(DMDestroy(&plex));
647:       PetscCall(DMDestroy(&plexIn));
648:       if (dmAux) PetscCall(DMDestroy(&plexAux));
649:       PetscFunctionReturn(PETSC_SUCCESS);
650:     }
651:     if (dim < 0) {
652:       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;

654:       /* Fall back to determination based on being a submesh */
655:       PetscCall(DMPlexGetSubpointMap(plex, &spmap));
656:       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
657:       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
658:       dim    = spmap ? 1 : 0;
659:       dimIn  = spmapIn ? 1 : 0;
660:       dimAux = spmapAux ? 1 : 0;
661:     }
662:     {
663:       PetscInt dimProj   = PetscMin(PetscMin(dim, dimIn), (dimAux < 0 ? PETSC_MAX_INT : dimAux));
664:       PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;

666:       PetscCheck(PetscAbsInt(dimProj - dim) <= 1 && PetscAbsInt(dimProj - dimIn) <= 1 && PetscAbsInt(dimProj - dimAuxEff) <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not currently support differences of more than 1 in dimension");
667:       if (dimProj < dim) minHeight = 1;
668:       htInc    = dim - dimProj;
669:       htIncIn  = dimIn - dimProj;
670:       htIncAux = dimAuxEff - dimProj;
671:     }
672:   }
673:   PetscCall(DMPlexGetDepth(plex, &depth));
674:   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
675:   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
676:   maxHeight = PetscMax(maxHeight, minHeight);
677:   PetscCheck(maxHeight >= 0 && maxHeight <= dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", maxHeight, dim);
678:   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, 0, NULL, &ds));
679:   if (!ds) PetscCall(DMGetDS(dm, &ds));
680:   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn));
681:   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
682:   PetscCall(PetscDSGetNumFields(ds, &Nf));
683:   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
684:   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
685:   if (isCohesiveIn) --htIncIn; // Should be rearranged
686:   PetscCall(DMGetNumFields(dm, &NfTot));
687:   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
688:   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL));
689:   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
690:   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
691:   PetscCall(DMGetLocalSection(dm, &section));
692:   if (dmAux) {
693:     PetscCall(DMGetDS(dmAux, &dsAux));
694:     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
695:   }
696:   PetscCall(PetscDSGetComponents(ds, &Nc));
697:   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
698:   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
699:   else {
700:     cellsp   = sp;
701:     cellspIn = spIn;
702:   }
703:   /* Get cell dual spaces */
704:   for (f = 0; f < Nf; ++f) {
705:     PetscDiscType disctype;

707:     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
708:     if (disctype == PETSC_DISC_FE) {
709:       PetscFE fe;

711:       isFE[f] = PETSC_TRUE;
712:       hasFE   = PETSC_TRUE;
713:       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
714:       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
715:     } else if (disctype == PETSC_DISC_FV) {
716:       PetscFV fv;

718:       isFE[f] = PETSC_FALSE;
719:       hasFV   = PETSC_TRUE;
720:       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv));
721:       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
722:     } else {
723:       isFE[f]   = PETSC_FALSE;
724:       cellsp[f] = NULL;
725:     }
726:   }
727:   for (f = 0; f < NfIn; ++f) {
728:     PetscDiscType disctype;

730:     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
731:     if (disctype == PETSC_DISC_FE) {
732:       PetscFE fe;

734:       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
735:       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
736:     } else if (disctype == PETSC_DISC_FV) {
737:       PetscFV fv;

739:       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv));
740:       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
741:     } else {
742:       cellspIn[f] = NULL;
743:     }
744:   }
745:   for (f = 0; f < Nf; ++f) {
746:     if (!htInc) {
747:       sp[f] = cellsp[f];
748:     } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
749:   }
750:   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
751:     PetscFE          fem, subfem;
752:     PetscDiscType    disctype;
753:     const PetscReal *points;
754:     PetscInt         numPoints;

756:     PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
757:     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints));
758:     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
759:     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
760:     for (f = 0; f < NfIn; ++f) {
761:       if (!htIncIn) {
762:         spIn[f] = cellspIn[f];
763:       } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));

765:       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
766:       if (disctype != PETSC_DISC_FE) continue;
767:       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem));
768:       if (!htIncIn) {
769:         subfem = fem;
770:       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
771:       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]));
772:     }
773:     for (f = 0; f < NfAux; ++f) {
774:       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
775:       if (disctype != PETSC_DISC_FE) continue;
776:       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem));
777:       if (!htIncAux) {
778:         subfem = fem;
779:       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
780:       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]));
781:     }
782:   }
783:   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
784:   for (h = minHeight; h <= maxHeight; h++) {
785:     PetscInt     hEff     = h - minHeight + htInc;
786:     PetscInt     hEffIn   = h - minHeight + htIncIn;
787:     PetscInt     hEffAux  = h - minHeight + htIncAux;
788:     PetscDS      dsEff    = ds;
789:     PetscDS      dsEffIn  = dsIn;
790:     PetscDS      dsEffAux = dsAux;
791:     PetscScalar *values;
792:     PetscBool   *fieldActive;
793:     PetscInt     maxDegree;
794:     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
795:     IS           heightIS;

797:     if (h > minHeight) {
798:       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
799:     }
800:     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
801:     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
802:     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
803:     if (pEnd <= pStart) {
804:       PetscCall(ISDestroy(&heightIS));
805:       continue;
806:     }
807:     /* Compute totDim, the number of dofs in the closure of a point at this height */
808:     totDim = 0;
809:     for (f = 0; f < Nf; ++f) {
810:       PetscBool cohesive;

812:       if (!sp[f]) continue;
813:       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
814:       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
815:       totDim += spDim;
816:       if (isCohesive && !cohesive) totDim += spDim;
817:     }
818:     p = lStart < 0 ? pStart : lStart;
819:     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
820:     PetscCheck(numValues == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The output section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT " in [%" PetscInt_FMT ", %" PetscInt_FMT "]", p, numValues, totDim, h, minHeight, maxHeight);
821:     if (!totDim) {
822:       PetscCall(ISDestroy(&heightIS));
823:       continue;
824:     }
825:     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
826:     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
827:     if (localU) {
828:       PetscInt totDimIn, pIn, numValuesIn;

830:       totDimIn = 0;
831:       for (f = 0; f < NfIn; ++f) {
832:         PetscBool cohesive;

834:         if (!spIn[f]) continue;
835:         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
836:         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
837:         totDimIn += spDim;
838:         if (isCohesiveIn && !cohesive) totDimIn += spDim;
839:       }
840:       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
841:       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
842:       // TODO We could check that pIn is a cohesive cell for this check
843:       PetscCheck(isCohesiveIn || (numValuesIn == totDimIn), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The input section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT, pIn, numValuesIn, totDimIn, htIncIn);
844:       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
845:     }
846:     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
847:     /* Loop over points at this height */
848:     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
849:     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
850:     {
851:       const PetscInt *fields;

853:       PetscCall(ISGetIndices(fieldIS, &fields));
854:       for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE;
855:       for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
856:       PetscCall(ISRestoreIndices(fieldIS, &fields));
857:     }
858:     if (label) {
859:       PetscInt i;

861:       for (i = 0; i < numIds; ++i) {
862:         IS              pointIS, isectIS;
863:         const PetscInt *points;
864:         PetscInt        n;
865:         PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
866:         PetscQuadrature quad = NULL;

868:         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
869:         if (!pointIS) continue; /* No points with that id on this process */
870:         PetscCall(ISIntersect(pointIS, heightIS, &isectIS));
871:         PetscCall(ISDestroy(&pointIS));
872:         if (!isectIS) continue;
873:         PetscCall(ISGetLocalSize(isectIS, &n));
874:         PetscCall(ISGetIndices(isectIS, &points));
875:         PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree));
876:         if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad));
877:         if (!quad) {
878:           if (!h && allPoints) {
879:             quad      = allPoints;
880:             allPoints = NULL;
881:           } else {
882:             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad));
883:           }
884:         }
885:         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
886:         for (p = 0; p < n; ++p) {
887:           const PetscInt point = points[p];

889:           PetscCall(PetscArrayzero(values, numValues));
890:           PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom));
891:           PetscCall(DMPlexSetActivePoint(dm, point));
892:           PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values));
893:           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
894:           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
895:         }
896:         PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom));
897:         PetscCall(PetscFEGeomDestroy(&fegeom));
898:         PetscCall(PetscQuadratureDestroy(&quad));
899:         PetscCall(ISRestoreIndices(isectIS, &points));
900:         PetscCall(ISDestroy(&isectIS));
901:       }
902:     } else {
903:       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
904:       PetscQuadrature quad = NULL;
905:       IS              pointIS;

907:       PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS));
908:       PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
909:       if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
910:       if (!quad) {
911:         if (!h && allPoints) {
912:           quad      = allPoints;
913:           allPoints = NULL;
914:         } else {
915:           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad));
916:         }
917:       }
918:       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
919:       for (p = pStart; p < pEnd; ++p) {
920:         PetscCall(PetscArrayzero(values, numValues));
921:         PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom));
922:         PetscCall(DMPlexSetActivePoint(dm, p));
923:         PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values));
924:         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
925:         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
926:       }
927:       PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom));
928:       PetscCall(PetscFEGeomDestroy(&fegeom));
929:       PetscCall(PetscQuadratureDestroy(&quad));
930:       PetscCall(ISDestroy(&pointIS));
931:     }
932:     PetscCall(ISDestroy(&heightIS));
933:     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
934:     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
935:   }
936:   /* Cleanup */
937:   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
938:     for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f]));
939:     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
940:     PetscCall(PetscFree2(T, TAux));
941:   }
942:   PetscCall(PetscQuadratureDestroy(&allPoints));
943:   PetscCall(PetscFree3(isFE, sp, spIn));
944:   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
945:   PetscCall(DMDestroy(&plex));
946:   PetscCall(DMDestroy(&plexIn));
947:   if (dmAux) PetscCall(DMDestroy(&plexAux));
948:   PetscFunctionReturn(PETSC_SUCCESS);
949: }

951: PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
952: {
953:   PetscFunctionBegin;
954:   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
955:   PetscFunctionReturn(PETSC_SUCCESS);
956: }

958: 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)
959: {
960:   PetscFunctionBegin;
961:   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
962:   PetscFunctionReturn(PETSC_SUCCESS);
963: }

965: PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
966: {
967:   PetscFunctionBegin;
968:   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
969:   PetscFunctionReturn(PETSC_SUCCESS);
970: }

972: PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
973: {
974:   PetscFunctionBegin;
975:   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
976:   PetscFunctionReturn(PETSC_SUCCESS);
977: }

979: PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
980: {
981:   PetscFunctionBegin;
982:   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX));
983:   PetscFunctionReturn(PETSC_SUCCESS);
984: }