xref: /petsc/src/dm/impls/plex/plexproject.c (revision a3b724e88d046d3fedf46cc78f622355d8be981f)
147923291SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
247923291SMatthew G. Knepley 
38c6c5593SMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
48c6c5593SMatthew G. Knepley 
51b32699bSMatthew G. Knepley /*@
61b32699bSMatthew G. Knepley   DMPlexGetActivePoint - Get the point on which projection is currently working
71b32699bSMatthew G. Knepley 
820f4b53cSBarry Smith   Not Collective
91b32699bSMatthew G. Knepley 
104165533cSJose E. Roman   Input Parameter:
11a1cb98faSBarry Smith . dm - the `DM`
121b32699bSMatthew G. Knepley 
134165533cSJose E. Roman   Output Parameter:
141b32699bSMatthew G. Knepley . point - The mesh point involved in the current projection
151b32699bSMatthew G. Knepley 
161b32699bSMatthew G. Knepley   Level: developer
171b32699bSMatthew G. Knepley 
181cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`
191b32699bSMatthew G. Knepley @*/
20d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point)
21d71ae5a4SJacob Faibussowitsch {
221b32699bSMatthew G. Knepley   PetscFunctionBeginHot;
231b32699bSMatthew G. Knepley   *point = ((DM_Plex *)dm->data)->activePoint;
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251b32699bSMatthew G. Knepley }
261b32699bSMatthew G. Knepley 
271b32699bSMatthew G. Knepley /*@
281b32699bSMatthew G. Knepley   DMPlexSetActivePoint - Set the point on which projection is currently working
291b32699bSMatthew G. Knepley 
3020f4b53cSBarry Smith   Not Collective
311b32699bSMatthew G. Knepley 
324165533cSJose E. Roman   Input Parameters:
33a1cb98faSBarry Smith + dm    - the `DM`
341b32699bSMatthew G. Knepley - point - The mesh point involved in the current projection
351b32699bSMatthew G. Knepley 
361b32699bSMatthew G. Knepley   Level: developer
371b32699bSMatthew G. Knepley 
381cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetActivePoint()`
391b32699bSMatthew G. Knepley @*/
40d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point)
41d71ae5a4SJacob Faibussowitsch {
421b32699bSMatthew G. Knepley   PetscFunctionBeginHot;
431b32699bSMatthew G. Knepley   ((DM_Plex *)dm->data)->activePoint = point;
443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
451b32699bSMatthew G. Knepley }
461b32699bSMatthew G. Knepley 
47a6e0b375SMatthew G. Knepley /*
48a6e0b375SMatthew G. Knepley   DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point
49a6e0b375SMatthew G. Knepley 
50a6e0b375SMatthew G. Knepley   Input Parameters:
51a1cb98faSBarry Smith + dm     - The output `DM`
52a1cb98faSBarry Smith . ds     - The output `DS`
53a1cb98faSBarry Smith . dmIn   - The input `DM`
54a1cb98faSBarry Smith . dsIn   - The input `DS`
55a6e0b375SMatthew G. Knepley . time   - The time for this evaluation
56a6e0b375SMatthew G. Knepley . fegeom - The FE geometry for this point
57a6e0b375SMatthew G. Knepley . fvgeom - The FV geometry for this point
58a6e0b375SMatthew G. Knepley . isFE   - Flag indicating whether each output field has an FE discretization
59a1cb98faSBarry Smith . sp     - The output `PetscDualSpace` for each field
60a6e0b375SMatthew G. Knepley . funcs  - The evaluation function for each field
61a6e0b375SMatthew G. Knepley - ctxs   - The user context for each field
62a6e0b375SMatthew G. Knepley 
63a6e0b375SMatthew G. Knepley   Output Parameter:
64a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space
65a6e0b375SMatthew G. Knepley 
66a6e0b375SMatthew G. Knepley   Level: developer
67a6e0b375SMatthew G. Knepley 
681cc06b55SBarry Smith .seealso:[](ch_unstructured), `DM`, `DMPLEX`, `PetscDS`, `PetscFEGeom`, `PetscFVCellGeom`, `PetscDualSpace`
69a6e0b375SMatthew G. Knepley */
70d71ae5a4SJacob Faibussowitsch 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[])
71d71ae5a4SJacob Faibussowitsch {
72a6e0b375SMatthew G. Knepley   PetscInt  coordDim, Nf, *Nc, f, spDim, d, v, tp;
735fedec97SMatthew G. Knepley   PetscBool isAffine, isCohesive, transform;
748c6c5593SMatthew G. Knepley 
758c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
769566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmIn, &coordDim));
779566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
809566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
818c6c5593SMatthew G. Knepley   /* Get values for closure */
82c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
83c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
848c6c5593SMatthew G. Knepley     void *const ctx = ctxs ? ctxs[f] : NULL;
855fedec97SMatthew G. Knepley     PetscBool   cohesive;
868c6c5593SMatthew G. Knepley 
878c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
899566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
908c6c5593SMatthew G. Knepley     if (funcs[f]) {
91c330f8ffSToby Isaac       if (isFE[f]) {
92c330f8ffSToby Isaac         PetscQuadrature  allPoints;
93c330f8ffSToby Isaac         PetscInt         q, dim, numPoints;
94c330f8ffSToby Isaac         const PetscReal *points;
95c330f8ffSToby Isaac         PetscScalar     *pointEval;
96c330f8ffSToby Isaac         PetscReal       *x;
97ca3d3a14SMatthew G. Knepley         DM               rdm;
98c330f8ffSToby Isaac 
999566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDM(sp[f], &rdm));
1009566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
1019566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
1029566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
1039566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x));
104737e23dcSMatthew G. Knepley         PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f]));
105c330f8ffSToby Isaac         for (q = 0; q < numPoints; q++, tp++) {
106c330f8ffSToby Isaac           const PetscReal *v0;
107c330f8ffSToby Isaac 
108c330f8ffSToby Isaac           if (isAffine) {
109665f567fSMatthew G. Knepley             const PetscReal *refpoint    = &points[q * dim];
110665f567fSMatthew G. Knepley             PetscReal        injpoint[3] = {0., 0., 0.};
111665f567fSMatthew G. Knepley 
112665f567fSMatthew G. Knepley             if (dim != fegeom->dim) {
1135fedec97SMatthew G. Knepley               if (isCohesive) {
114665f567fSMatthew G. Knepley                 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
115665f567fSMatthew G. Knepley                 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
116665f567fSMatthew G. Knepley                 refpoint = injpoint;
11763a3b9bcSJacob Faibussowitsch               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim);
118665f567fSMatthew G. Knepley             }
119665f567fSMatthew G. Knepley             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x);
120c330f8ffSToby Isaac             v0 = x;
1218c6c5593SMatthew G. Knepley           } else {
122c330f8ffSToby Isaac             v0 = &fegeom->v[tp * coordDim];
1238c6c5593SMatthew G. Knepley           }
1249371c9d4SSatish Balay           if (transform) {
1259371c9d4SSatish Balay             PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx));
1269371c9d4SSatish Balay             v0 = x;
1279371c9d4SSatish Balay           }
1289566063dSJacob Faibussowitsch           PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx));
129c330f8ffSToby Isaac         }
1304bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
1319566063dSJacob Faibussowitsch         PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
1329566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
1339566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x));
1349566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
135c330f8ffSToby Isaac         v += spDim;
1365fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) {
13727f02ce8SMatthew G. Knepley           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
13827f02ce8SMatthew G. Knepley         }
139c330f8ffSToby Isaac       } else {
14048a46eb9SPierre Jolivet         for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
141c330f8ffSToby Isaac       }
142c330f8ffSToby Isaac     } else {
14327f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
1445fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
14527f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
14627f02ce8SMatthew G. Knepley       }
1478c6c5593SMatthew G. Knepley     }
1489c3cf19fSMatthew G. Knepley   }
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1508c6c5593SMatthew G. Knepley }
1518c6c5593SMatthew G. Knepley 
152a6e0b375SMatthew G. Knepley /*
153a6e0b375SMatthew G. Knepley   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
154a6e0b375SMatthew G. Knepley 
155a6e0b375SMatthew G. Knepley   Input Parameters:
156a6e0b375SMatthew G. Knepley + dm             - The output DM
157a6e0b375SMatthew G. Knepley . ds             - The output DS
158a6e0b375SMatthew G. Knepley . dmIn           - The input DM
159a6e0b375SMatthew G. Knepley . dsIn           - The input DS
160a6e0b375SMatthew G. Knepley . dmAux          - The auxiliary DM, which is always for the input space
161a6e0b375SMatthew G. Knepley . dsAux          - The auxiliary DS, which is always for the input space
162a6e0b375SMatthew G. Knepley . time           - The time for this evaluation
163a6e0b375SMatthew G. Knepley . localU         - The local solution
164a6e0b375SMatthew G. Knepley . localA         - The local auziliary fields
165a6e0b375SMatthew G. Knepley . cgeom          - The FE geometry for this point
166a6e0b375SMatthew G. Knepley . sp             - The output PetscDualSpace for each field
167a6e0b375SMatthew G. Knepley . p              - The point in the output DM
168a5b23f4aSJose E. Roman . T              - Input basis and derivatives for each field tabulated on the quadrature points
169ef0bb6c7SMatthew G. Knepley . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
170a6e0b375SMatthew G. Knepley . funcs          - The evaluation function for each field
171a6e0b375SMatthew G. Knepley - ctxs           - The user context for each field
172a6e0b375SMatthew G. Knepley 
173a6e0b375SMatthew G. Knepley   Output Parameter:
174a6e0b375SMatthew G. Knepley . values         - The value for each dual basis vector in the output dual space
175a6e0b375SMatthew G. Knepley 
176a6e0b375SMatthew G. Knepley   Level: developer
177a6e0b375SMatthew G. Knepley 
178a1cb98faSBarry Smith   Note:
179a1cb98faSBarry Smith   Not supported for FV
180a1cb98faSBarry Smith 
181db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()`
182a6e0b375SMatthew G. Knepley */
183d71ae5a4SJacob Faibussowitsch 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[])
184d71ae5a4SJacob Faibussowitsch {
1858c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
1864bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
1878c6c5593SMatthew G. Knepley   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
1888c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
189191494d9SMatthew G. Knepley   const PetscScalar *constants;
1908c6c5593SMatthew G. Knepley   PetscReal         *x;
1918e6d238bSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
1924bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
1938e6d238bSMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed, *cone, *ornt;
194ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
1958e6d238bSMatthew G. Knepley   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
1968e6d238bSMatthew G. Knepley   DMPolytopeType     qct;
1978c6c5593SMatthew G. Knepley 
1988c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
2019566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
2029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
2038e6d238bSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
2049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
2059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
2069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
2079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
2089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
2099566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
2109566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmIn, &section));
2119566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
2128e6d238bSMatthew G. Knepley   // Get cohesive cell hanging off face
2138e6d238bSMatthew G. Knepley   if (isCohesiveIn) {
2148e6d238bSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
2158e6d238bSMatthew G. Knepley     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)) {
2168e6d238bSMatthew G. Knepley       DMPolytopeType  ct;
2178e6d238bSMatthew G. Knepley       const PetscInt *support;
2188e6d238bSMatthew G. Knepley       PetscInt        Ns, s;
2198e6d238bSMatthew G. Knepley 
2208e6d238bSMatthew G. Knepley       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
2218e6d238bSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
2228e6d238bSMatthew G. Knepley       for (s = 0; s < Ns; ++s) {
2238e6d238bSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
2248e6d238bSMatthew G. Knepley         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)) {
2258e6d238bSMatthew G. Knepley           inp = support[s];
2268e6d238bSMatthew G. Knepley           break;
2278e6d238bSMatthew G. Knepley         }
2288e6d238bSMatthew G. Knepley       }
2298e6d238bSMatthew G. Knepley       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
2308e6d238bSMatthew G. Knepley       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
2318e6d238bSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
2328e6d238bSMatthew G. Knepley       face[0] = 0;
2338e6d238bSMatthew G. Knepley       face[1] = 0;
2348e6d238bSMatthew G. Knepley     }
2358e6d238bSMatthew G. Knepley   }
236eb8f539aSJed Brown   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
2378c6c5593SMatthew G. Knepley   if (dmAux) {
23844171101SMatthew G. Knepley     PetscInt subp;
2391ba34526SMatthew G. Knepley 
2409566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
2419566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2429566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
2439566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2449566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2459566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
2469566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
2478c6c5593SMatthew G. Knepley   }
2488c6c5593SMatthew G. Knepley   /* Get values for closure */
2494bee2e38SMatthew G. Knepley   isAffine        = cgeom->isAffine;
25027f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
25127f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
2524bee2e38SMatthew G. Knepley   if (isAffine) {
2534bee2e38SMatthew G. Knepley     fegeom.v    = x;
2544bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
2554bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
2564bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
2574bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
2584bee2e38SMatthew G. Knepley   }
2598c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
260c330f8ffSToby Isaac     PetscQuadrature  allPoints;
261c330f8ffSToby Isaac     PetscInt         q, dim, numPoints;
262c330f8ffSToby Isaac     const PetscReal *points;
263c330f8ffSToby Isaac     PetscScalar     *pointEval;
2645fedec97SMatthew G. Knepley     PetscBool        cohesive;
265c330f8ffSToby Isaac     DM               dm;
266c330f8ffSToby Isaac 
2678c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2689566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
2699566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
270c330f8ffSToby Isaac     if (!funcs[f]) {
271be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
2725fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
27327f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
27427f02ce8SMatthew G. Knepley       }
275c330f8ffSToby Isaac       continue;
276c330f8ffSToby Isaac     }
2779566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
2789566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
2799566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
2809566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
2810c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
2828e6d238bSMatthew G. Knepley       PetscInt qpt[2];
2838e6d238bSMatthew G. Knepley 
2848e6d238bSMatthew G. Knepley       if (isCohesiveIn) {
2858e6d238bSMatthew G. Knepley         PetscCall(PetscDSPermuteQuadPoint(dsIn, ornt[0], f, q, &qpt[0]));
2868e6d238bSMatthew G. Knepley         PetscCall(PetscDSPermuteQuadPoint(dsIn, DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0), f, q, &qpt[1]));
2878e6d238bSMatthew G. Knepley       }
288c330f8ffSToby Isaac       if (isAffine) {
289ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x);
2901c531cf8SMatthew G. Knepley       } else {
2914bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp * dE];
2924bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp * dE * dE];
2934bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp * dE * dE];
2944bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
2958c6c5593SMatthew G. Knepley       }
2968e6d238bSMatthew G. Knepley       if (coefficients) {
2978e6d238bSMatthew G. Knepley         if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
2988e6d238bSMatthew G. Knepley         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
2998e6d238bSMatthew G. Knepley       }
3009566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
3019566063dSJacob Faibussowitsch       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
302a6e0b375SMatthew G. Knepley       (*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]);
3031c531cf8SMatthew G. Knepley     }
3049566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
3059566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
306c330f8ffSToby Isaac     v += spDim;
30727f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
3085fedec97SMatthew G. Knepley     if (isCohesive && !cohesive) {
30927f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
31027f02ce8SMatthew G. Knepley     }
3118c6c5593SMatthew G. Knepley   }
312eb8f539aSJed Brown   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
3139566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
3148e6d238bSMatthew G. Knepley   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3168c6c5593SMatthew G. Knepley }
3178c6c5593SMatthew G. Knepley 
31879f2dbaeSMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, 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[])
319d71ae5a4SJacob Faibussowitsch {
320b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
3214bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
322b18799e0SMatthew G. Knepley   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
323b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
324b18799e0SMatthew G. Knepley   const PetscScalar *constants;
325b18799e0SMatthew G. Knepley   PetscReal         *x;
32679f2dbaeSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
3279f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
32879f2dbaeSMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed, *cone, *ornt;
32979f2dbaeSMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
33079f2dbaeSMatthew G. Knepley   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
33179f2dbaeSMatthew G. Knepley   DMPolytopeType     qct;
332b18799e0SMatthew G. Knepley 
333b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
3349566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3359566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
33679f2dbaeSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
33779f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
33879f2dbaeSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
33979f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
34079f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
34179f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
34279f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
34379f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
34479f2dbaeSMatthew G. Knepley   PetscCall(DMHasBasisTransform(dmIn, &transform));
34579f2dbaeSMatthew G. Knepley   PetscCall(DMGetLocalSection(dmIn, &section));
34679f2dbaeSMatthew G. Knepley   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
34779f2dbaeSMatthew G. Knepley   // Get cohesive cell hanging off face
34879f2dbaeSMatthew G. Knepley   if (isCohesiveIn) {
34979f2dbaeSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
35079f2dbaeSMatthew G. Knepley     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)) {
35179f2dbaeSMatthew G. Knepley       DMPolytopeType  ct;
35279f2dbaeSMatthew G. Knepley       const PetscInt *support;
35379f2dbaeSMatthew G. Knepley       PetscInt        Ns, s;
35479f2dbaeSMatthew G. Knepley 
35579f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
35679f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
35779f2dbaeSMatthew G. Knepley       for (s = 0; s < Ns; ++s) {
35879f2dbaeSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
35979f2dbaeSMatthew G. Knepley         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)) {
36079f2dbaeSMatthew G. Knepley           inp = support[s];
36179f2dbaeSMatthew G. Knepley           break;
36279f2dbaeSMatthew G. Knepley         }
36379f2dbaeSMatthew G. Knepley       }
36479f2dbaeSMatthew G. Knepley       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
36579f2dbaeSMatthew G. Knepley       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
36679f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
36779f2dbaeSMatthew G. Knepley       face[0] = 0;
36879f2dbaeSMatthew G. Knepley       face[1] = 0;
36979f2dbaeSMatthew G. Knepley     }
37079f2dbaeSMatthew G. Knepley   }
37179f2dbaeSMatthew G. Knepley   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
372b18799e0SMatthew G. Knepley   if (dmAux) {
373a6e0b375SMatthew G. Knepley     PetscInt subp;
374b18799e0SMatthew G. Knepley 
3759566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
3769566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3779566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
3789566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3799566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3809566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
3819566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
382b18799e0SMatthew G. Knepley   }
383b18799e0SMatthew G. Knepley   /* Get values for closure */
3844bee2e38SMatthew G. Knepley   isAffine       = fgeom->isAffine;
385ea78f98cSLisandro Dalcin   fegeom.n       = NULL;
386ea78f98cSLisandro Dalcin   fegeom.J       = NULL;
387ea78f98cSLisandro Dalcin   fegeom.v       = NULL;
388ea78f98cSLisandro Dalcin   fegeom.xi      = NULL;
389a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
390a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
3914bee2e38SMatthew G. Knepley   if (isAffine) {
3924bee2e38SMatthew G. Knepley     fegeom.v    = x;
3934bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
3944bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
3954bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
3964bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
3974bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
3989f209ee4SMatthew G. Knepley 
3999f209ee4SMatthew G. Knepley     cgeom.J    = fgeom->suppJ[0];
4009f209ee4SMatthew G. Knepley     cgeom.invJ = fgeom->suppInvJ[0];
4019f209ee4SMatthew G. Knepley     cgeom.detJ = fgeom->suppDetJ[0];
4024bee2e38SMatthew G. Knepley   }
403b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
404b18799e0SMatthew G. Knepley     PetscQuadrature  allPoints;
405b18799e0SMatthew G. Knepley     PetscInt         q, dim, numPoints;
406b18799e0SMatthew G. Knepley     const PetscReal *points;
407b18799e0SMatthew G. Knepley     PetscScalar     *pointEval;
40879f2dbaeSMatthew G. Knepley     PetscBool        cohesive;
409b18799e0SMatthew G. Knepley     DM               dm;
410b18799e0SMatthew G. Knepley 
411b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
41279f2dbaeSMatthew G. Knepley     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
4139566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
414b18799e0SMatthew G. Knepley     if (!funcs[f]) {
415b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
41679f2dbaeSMatthew G. Knepley       if (isCohesive && !cohesive) {
41779f2dbaeSMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
41879f2dbaeSMatthew G. Knepley       }
419b18799e0SMatthew G. Knepley       continue;
420b18799e0SMatthew G. Knepley     }
4219566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
4229566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
4239566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
4249566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
425b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
42679f2dbaeSMatthew G. Knepley       PetscInt qpt[2];
42779f2dbaeSMatthew G. Knepley 
42879f2dbaeSMatthew G. Knepley       if (isCohesiveIn) {
42979f2dbaeSMatthew G. Knepley         // These points are not integration quadratures, but dual space quadratures
430d8b4a066SPierre Jolivet         // If they had multiple points we should match them from both sides, similar to hybrid residual eval
43179f2dbaeSMatthew G. Knepley         qpt[0] = qpt[1] = q;
43279f2dbaeSMatthew G. Knepley       }
433b18799e0SMatthew G. Knepley       if (isAffine) {
434ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x);
435b18799e0SMatthew G. Knepley       } else {
4363fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp * dE];
4379f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp * dE * dE];
4389f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp * dE * dE];
4394bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
4404bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp * dE];
4419f209ee4SMatthew G. Knepley 
4429f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][tp * dE * dE];
4439f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE];
4449f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][tp];
445b18799e0SMatthew G. Knepley       }
446a6e0b375SMatthew G. Knepley       /* TODO We should use cgeom here, instead of fegeom, however the geometry coming in through fgeom does not have the support cell geometry */
44779f2dbaeSMatthew G. Knepley       if (coefficients) {
44879f2dbaeSMatthew G. Knepley         if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
44979f2dbaeSMatthew G. Knepley         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
45079f2dbaeSMatthew G. Knepley       }
4519566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
45279f2dbaeSMatthew G. Knepley       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
45379f2dbaeSMatthew G. Knepley       (*funcs[f])(dE, NfIn, 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]);
454b18799e0SMatthew G. Knepley     }
4559566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
4569566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
457b18799e0SMatthew G. Knepley     v += spDim;
45879f2dbaeSMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
45979f2dbaeSMatthew G. Knepley     if (isCohesive && !cohesive) {
46079f2dbaeSMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
461b18799e0SMatthew G. Knepley     }
46279f2dbaeSMatthew G. Knepley   }
46379f2dbaeSMatthew G. Knepley   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
4649566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
46579f2dbaeSMatthew G. Knepley   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
4663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
467b18799e0SMatthew G. Knepley }
468b18799e0SMatthew G. Knepley 
469d71ae5a4SJacob Faibussowitsch 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[])
470d71ae5a4SJacob Faibussowitsch {
4718c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
4728c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
4738c6c5593SMatthew G. Knepley 
4748c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
4759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4769566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
4779566063dSJacob Faibussowitsch   if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
4788c6c5593SMatthew G. Knepley   switch (type) {
4798c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
480d71ae5a4SJacob Faibussowitsch   case DM_BC_NATURAL:
481d71ae5a4SJacob Faibussowitsch     PetscCall(DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode(**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))funcs, ctxs, values));
482d71ae5a4SJacob Faibussowitsch     break;
4838c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
4848c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
4859371c9d4SSatish Balay     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));
486d0609cedSBarry Smith     break;
487b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
48879f2dbaeSMatthew G. Knepley     PetscCall(DMProjectPoint_BdField_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[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
489d0609cedSBarry Smith     break;
490d71ae5a4SJacob Faibussowitsch   default:
491d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type);
4928c6c5593SMatthew G. Knepley   }
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4948c6c5593SMatthew G. Knepley }
4958c6c5593SMatthew G. Knepley 
496d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
497d71ae5a4SJacob Faibussowitsch {
498df5c1128SToby Isaac   PetscReal *points;
499df5c1128SToby Isaac   PetscInt   f, numPoints;
500df5c1128SToby Isaac 
501df5c1128SToby Isaac   PetscFunctionBegin;
50219552267SMatthew G. Knepley   if (!dim) {
50319552267SMatthew G. Knepley     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5043ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
50519552267SMatthew G. Knepley   }
506df5c1128SToby Isaac   numPoints = 0;
507df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
508df5c1128SToby Isaac     if (funcs[f]) {
509df5c1128SToby Isaac       PetscQuadrature fAllPoints;
510df5c1128SToby Isaac       PetscInt        fNumPoints;
511df5c1128SToby Isaac 
5129566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5139566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
514df5c1128SToby Isaac       numPoints += fNumPoints;
515df5c1128SToby Isaac     }
516df5c1128SToby Isaac   }
5179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * numPoints, &points));
518df5c1128SToby Isaac   numPoints = 0;
519df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
520df5c1128SToby Isaac     if (funcs[f]) {
521df5c1128SToby Isaac       PetscQuadrature  fAllPoints;
52254ee0cceSMatthew G. Knepley       PetscInt         qdim, fNumPoints, q;
523df5c1128SToby Isaac       const PetscReal *fPoints;
524df5c1128SToby Isaac 
5259566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5269566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
52763a3b9bcSJacob Faibussowitsch       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);
528df5c1128SToby Isaac       for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q];
529df5c1128SToby Isaac       numPoints += fNumPoints;
530df5c1128SToby Isaac     }
531df5c1128SToby Isaac   }
5329566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5339566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL));
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
535df5c1128SToby Isaac }
536df5c1128SToby Isaac 
5375f15299fSJeremy L Thompson /*@C
53820f4b53cSBarry Smith   DMGetFirstLabeledPoint - Find first labeled `point` in `odm` such that the corresponding point in `dm` has the specified `height`. Return `point` and the corresponding `ds`.
5395f15299fSJeremy L Thompson 
5405f15299fSJeremy L Thompson   Input Parameters:
5412fe279fdSBarry Smith + dm     - the `DM`
5422fe279fdSBarry Smith . odm    - the enclosing `DM`
5432fe279fdSBarry Smith . label  - label for `DM` domain, or `NULL` for whole domain
5442fe279fdSBarry Smith . numIds - the number of `ids`
5452fe279fdSBarry Smith . ids    - An array of the label ids in sequence for the domain
5462fe279fdSBarry Smith - height - Height of target cells in `DMPLEX` topology
5475f15299fSJeremy L Thompson 
5485f15299fSJeremy L Thompson   Output Parameters:
5492fe279fdSBarry Smith + point - the first labeled point
550*a3b724e8SBarry Smith - ds    - the `PetscDS` corresponding to the first labeled point
5515f15299fSJeremy L Thompson 
5525f15299fSJeremy L Thompson   Level: developer
553a1cb98faSBarry Smith 
5541cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS`
555817da375SSatish Balay @*/
556d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
557d71ae5a4SJacob Faibussowitsch {
558a6e0b375SMatthew G. Knepley   DM              plex;
559a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
560e9f0ba4eSJed Brown   PetscInt        ls = -1;
561e5e52638SMatthew G. Knepley 
562e5e52638SMatthew G. Knepley   PetscFunctionBegin;
5635f15299fSJeremy L Thompson   if (point) *point = -1;
5643ba16761SJacob Faibussowitsch   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
5659566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
5669566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
567e9f0ba4eSJed Brown   for (PetscInt i = 0; i < numIds; ++i) {
568e9f0ba4eSJed Brown     IS       labelIS;
569e9f0ba4eSJed Brown     PetscInt num_points, pStart, pEnd;
5709566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
571e9f0ba4eSJed Brown     if (!labelIS) continue; /* No points with that id on this process */
5729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
5739566063dSJacob Faibussowitsch     PetscCall(ISGetSize(labelIS, &num_points));
574e9f0ba4eSJed Brown     if (num_points) {
575e5e52638SMatthew G. Knepley       const PetscInt *points;
5769566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(labelIS, &points));
577e9f0ba4eSJed Brown       for (PetscInt i = 0; i < num_points; i++) {
578e9f0ba4eSJed Brown         PetscInt point;
5799566063dSJacob Faibussowitsch         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
580e9f0ba4eSJed Brown         if (pStart <= point && point < pEnd) {
581a6e0b375SMatthew G. Knepley           ls = point;
58279f2dbaeSMatthew G. Knepley           if (ds) {
58379f2dbaeSMatthew G. Knepley             // If this is a face of a cohesive cell, then prefer that DS
58479f2dbaeSMatthew G. Knepley             if (height == 1) {
58579f2dbaeSMatthew G. Knepley               const PetscInt *supp;
58679f2dbaeSMatthew G. Knepley               PetscInt        suppSize;
58779f2dbaeSMatthew G. Knepley               DMPolytopeType  ct;
58879f2dbaeSMatthew G. Knepley 
58979f2dbaeSMatthew G. Knepley               PetscCall(DMPlexGetSupport(dm, ls, &supp));
59079f2dbaeSMatthew G. Knepley               PetscCall(DMPlexGetSupportSize(dm, ls, &suppSize));
59179f2dbaeSMatthew G. Knepley               for (PetscInt s = 0; s < suppSize; ++s) {
59279f2dbaeSMatthew G. Knepley                 PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
59379f2dbaeSMatthew G. Knepley                 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)) {
59479f2dbaeSMatthew G. Knepley                   ls = supp[s];
59579f2dbaeSMatthew G. Knepley                   break;
59679f2dbaeSMatthew G. Knepley                 }
59779f2dbaeSMatthew G. Knepley               }
59879f2dbaeSMatthew G. Knepley             }
59979f2dbaeSMatthew G. Knepley             PetscCall(DMGetCellDS(dm, ls, ds, NULL));
60079f2dbaeSMatthew G. Knepley           }
6018e6d238bSMatthew G. Knepley           if (ls >= 0) break;
602e5e52638SMatthew G. Knepley         }
603e9f0ba4eSJed Brown       }
6049566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(labelIS, &points));
605e9f0ba4eSJed Brown     }
6069566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&labelIS));
607e5e52638SMatthew G. Knepley     if (ls >= 0) break;
608e5e52638SMatthew G. Knepley   }
6095f15299fSJeremy L Thompson   if (point) *point = ls;
6109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
612e5e52638SMatthew G. Knepley }
613e5e52638SMatthew G. Knepley 
6140de51b56SMatthew G. Knepley /*
6150de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
6160de51b56SMatthew G. Knepley 
6170de51b56SMatthew G. Knepley   There are several different scenarios:
6180de51b56SMatthew G. Knepley 
6190de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
6200de51b56SMatthew G. Knepley 
6210de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
6220de51b56SMatthew G. Knepley 
6230de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
6240de51b56SMatthew G. Knepley 
6250de51b56SMatthew G. Knepley      Here minHeight=1 since we loop over faces. This normally happens since we hang cells off of our boundary meshes to facilitate computation.
6260de51b56SMatthew G. Knepley 
6270de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
6280de51b56SMatthew G. Knepley 
6290de51b56SMatthew G. Knepley      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.
6300de51b56SMatthew G. Knepley 
631636b9322SMatthew G. Knepley   4) Volumetric input mesh with boundary output mesh
632636b9322SMatthew G. Knepley 
633636b9322SMatthew G. Knepley      Here we must get a subspace for the input DS
634636b9322SMatthew G. Knepley 
6350de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
6360de51b56SMatthew G. Knepley 
6370de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
6380de51b56SMatthew G. Knepley 
6390de51b56SMatthew G. Knepley   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.
6400de51b56SMatthew G. Knepley     - 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
6415f790a90SMatthew G. Knepley       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
6420de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
6430de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
6440de51b56SMatthew G. Knepley 
6450de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
6460de51b56SMatthew G. Knepley 
6470de51b56SMatthew G. Knepley   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.
6480de51b56SMatthew G. Knepley */
649d71ae5a4SJacob Faibussowitsch 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)
650d71ae5a4SJacob Faibussowitsch {
651a6e0b375SMatthew G. Knepley   DM               plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
652a6e0b375SMatthew G. Knepley   DMEnclosureType  encIn, encAux;
653a6e0b375SMatthew G. Knepley   PetscDS          ds = NULL, dsIn = NULL, dsAux = NULL;
654ca3d3a14SMatthew G. Knepley   Vec              localA = NULL, tv;
655aa0eca99SMatthew G. Knepley   IS               fieldIS;
65647923291SMatthew G. Knepley   PetscSection     section;
657636b9322SMatthew G. Knepley   PetscDualSpace  *sp, *cellsp, *spIn, *cellspIn;
658ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
6598c6c5593SMatthew G. Knepley   PetscInt        *Nc;
66079f2dbaeSMatthew G. Knepley   PetscInt         dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, minHeightIn, minHeightAux, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
6618e6d238bSMatthew G. Knepley   PetscBool       *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform;
662c330f8ffSToby Isaac   DMField          coordField;
663c330f8ffSToby Isaac   DMLabel          depthLabel;
664c330f8ffSToby Isaac   PetscQuadrature  allPoints = NULL;
66547923291SMatthew G. Knepley 
66647923291SMatthew G. Knepley   PetscFunctionBegin;
6679566063dSJacob Faibussowitsch   if (localU) PetscCall(VecGetDM(localU, &dmIn));
668ad540459SPierre Jolivet   else dmIn = dm;
6699566063dSJacob Faibussowitsch   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
6709371c9d4SSatish Balay   if (localA) PetscCall(VecGetDM(localA, &dmAux));
671ad540459SPierre Jolivet   else dmAux = NULL;
6729566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
6739566063dSJacob Faibussowitsch   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
6749566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
6759566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
6769566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
6779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
6789566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
6799566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
6809566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
6810de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
682a6e0b375SMatthew G. Knepley   if (dmAux) {
6839566063dSJacob Faibussowitsch     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
684ad540459SPierre Jolivet     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");
685636b9322SMatthew G. Knepley   }
6865ee220baSJed Brown   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plexIn, PETSC_TRUE, localU, time, NULL, NULL, NULL));
687e04ae0b4SMatthew G. Knepley   PetscCall(DMGetCoordinateField(dm, &coordField));
688e44f6aebSMatthew G. Knepley   PetscCheck(coordField, PETSC_COMM_SELF, PETSC_ERR_USER, "DM must have a coordinate field");
689e04ae0b4SMatthew G. Knepley   /**** No collective calls below this point ****/
690636b9322SMatthew G. Knepley   /* Determine height for iteration of all meshes */
691636b9322SMatthew G. Knepley   {
692636b9322SMatthew G. Knepley     DMPolytopeType ct, ctIn, ctAux;
69379f2dbaeSMatthew G. Knepley     PetscInt       lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux;
69419552267SMatthew G. Knepley     PetscInt       dim = -1, dimIn = -1, dimAux = -1;
69588aca1feSMatthew G. Knepley 
6969566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
697636b9322SMatthew G. Knepley     if (pEnd > pStart) {
6989566063dSJacob Faibussowitsch       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
699636b9322SMatthew G. Knepley       p = lStart < 0 ? pStart : lStart;
7009566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plex, p, &ct));
701636b9322SMatthew G. Knepley       dim = DMPolytopeTypeGetDim(ct);
7029566063dSJacob Faibussowitsch       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
7039566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
7049566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
705636b9322SMatthew G. Knepley       dimIn = DMPolytopeTypeGetDim(ctIn);
706636b9322SMatthew G. Knepley       if (dmAux) {
7079566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
70819552267SMatthew G. Knepley         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
70919552267SMatthew G. Knepley         if (pStartAux < pEndAux) {
7109566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
711636b9322SMatthew G. Knepley           dimAux = DMPolytopeTypeGetDim(ctAux);
71219552267SMatthew G. Knepley         }
713636b9322SMatthew G. Knepley       } else dimAux = dim;
714e04ae0b4SMatthew G. Knepley     } else {
715e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plex));
716e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plexIn));
717e04ae0b4SMatthew G. Knepley       if (dmAux) PetscCall(DMDestroy(&plexAux));
7183ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
719e39fcb4eSStefano Zampini     }
720636b9322SMatthew G. Knepley     if (dim < 0) {
721636b9322SMatthew G. Knepley       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
722636b9322SMatthew G. Knepley 
723636b9322SMatthew G. Knepley       /* Fall back to determination based on being a submesh */
7249566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plex, &spmap));
7259566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
7269566063dSJacob Faibussowitsch       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
727636b9322SMatthew G. Knepley       dim    = spmap ? 1 : 0;
728636b9322SMatthew G. Knepley       dimIn  = spmapIn ? 1 : 0;
729636b9322SMatthew G. Knepley       dimAux = spmapAux ? 1 : 0;
73088aca1feSMatthew G. Knepley     }
731636b9322SMatthew G. Knepley     {
73219552267SMatthew G. Knepley       PetscInt dimProj   = PetscMin(PetscMin(dim, dimIn), (dimAux < 0 ? PETSC_MAX_INT : dimAux));
73319552267SMatthew G. Knepley       PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;
734636b9322SMatthew G. Knepley 
73519552267SMatthew G. Knepley       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");
736636b9322SMatthew G. Knepley       if (dimProj < dim) minHeight = 1;
737636b9322SMatthew G. Knepley       htInc    = dim - dimProj;
738636b9322SMatthew G. Knepley       htIncIn  = dimIn - dimProj;
73919552267SMatthew G. Knepley       htIncAux = dimAuxEff - dimProj;
7400de51b56SMatthew G. Knepley     }
741a6e0b375SMatthew G. Knepley   }
7429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
7439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
7449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
7452af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
74663a3b9bcSJacob Faibussowitsch   PetscCheck(maxHeight >= 0 && maxHeight <= dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", maxHeight, dim);
74779f2dbaeSMatthew G. Knepley   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, NULL, &ds));
7489566063dSJacob Faibussowitsch   if (!ds) PetscCall(DMGetDS(dm, &ds));
74979f2dbaeSMatthew G. Knepley   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, minHeight, NULL, &dsIn));
7509566063dSJacob Faibussowitsch   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
7519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
7529566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
7538e6d238bSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
7548e6d238bSMatthew G. Knepley   if (isCohesiveIn) --htIncIn; // Should be rearranged
7559566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &NfTot));
7569566063dSJacob Faibussowitsch   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
75707218a29SMatthew G. Knepley   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL));
7589566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
7599566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
7609566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
7610c898c18SMatthew G. Knepley   if (dmAux) {
7629566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmAux, &dsAux));
7639566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
7640c898c18SMatthew G. Knepley   }
7659566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
7669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
7679566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
7689371c9d4SSatish Balay   else {
7699371c9d4SSatish Balay     cellsp   = sp;
7709371c9d4SSatish Balay     cellspIn = spIn;
7719371c9d4SSatish Balay   }
7728c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
77347923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
774665f567fSMatthew G. Knepley     PetscDiscType disctype;
77547923291SMatthew G. Knepley 
7769566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
777665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
778665f567fSMatthew G. Knepley       PetscFE fe;
77947923291SMatthew G. Knepley 
78047923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
781665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
7829566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
7839566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
784665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
785665f567fSMatthew G. Knepley       PetscFV fv;
78647923291SMatthew G. Knepley 
78747923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
788665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
7899566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv));
7909566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
791665f567fSMatthew G. Knepley     } else {
792665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
793665f567fSMatthew G. Knepley       cellsp[f] = NULL;
794665f567fSMatthew G. Knepley     }
79547923291SMatthew G. Knepley   }
796636b9322SMatthew G. Knepley   for (f = 0; f < NfIn; ++f) {
797636b9322SMatthew G. Knepley     PetscDiscType disctype;
798636b9322SMatthew G. Knepley 
7999566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
800636b9322SMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
801636b9322SMatthew G. Knepley       PetscFE fe;
802636b9322SMatthew G. Knepley 
8039566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
8049566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
805636b9322SMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
806636b9322SMatthew G. Knepley       PetscFV fv;
807636b9322SMatthew G. Knepley 
8089566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv));
8099566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
810636b9322SMatthew G. Knepley     } else {
811636b9322SMatthew G. Knepley       cellspIn[f] = NULL;
812636b9322SMatthew G. Knepley     }
813636b9322SMatthew G. Knepley   }
814636b9322SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
8159371c9d4SSatish Balay     if (!htInc) {
8169371c9d4SSatish Balay       sp[f] = cellsp[f];
8179371c9d4SSatish Balay     } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
818636b9322SMatthew G. Knepley   }
819ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
82074fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
821665f567fSMatthew G. Knepley     PetscDiscType    disctype;
8224a3e9fdbSToby Isaac     const PetscReal *points;
8234a3e9fdbSToby Isaac     PetscInt         numPoints;
8240c898c18SMatthew G. Knepley 
82508401ef6SPierre Jolivet     PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
8269566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints));
8279566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
8289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
829a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
8309371c9d4SSatish Balay       if (!htIncIn) {
8319371c9d4SSatish Balay         spIn[f] = cellspIn[f];
8329371c9d4SSatish Balay       } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
833636b9322SMatthew G. Knepley 
8349566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
835665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
8369566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem));
8379371c9d4SSatish Balay       if (!htIncIn) {
8389371c9d4SSatish Balay         subfem = fem;
8399371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
8409566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]));
8410c898c18SMatthew G. Knepley     }
8420c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
8439566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
844665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
8459566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem));
8469371c9d4SSatish Balay       if (!htIncAux) {
8479371c9d4SSatish Balay         subfem = fem;
8489371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
8499566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]));
8500c898c18SMatthew G. Knepley     }
8510c898c18SMatthew G. Knepley   }
85247923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
8532af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
854636b9322SMatthew G. Knepley     PetscInt     hEff     = h - minHeight + htInc;
855636b9322SMatthew G. Knepley     PetscInt     hEffIn   = h - minHeight + htIncIn;
856636b9322SMatthew G. Knepley     PetscInt     hEffAux  = h - minHeight + htIncAux;
857a6e0b375SMatthew G. Knepley     PetscDS      dsEff    = ds;
858636b9322SMatthew G. Knepley     PetscDS      dsEffIn  = dsIn;
859636b9322SMatthew G. Knepley     PetscDS      dsEffAux = dsAux;
8608c6c5593SMatthew G. Knepley     PetscScalar *values;
861b7260050SToby Isaac     PetscBool   *fieldActive;
862b7260050SToby Isaac     PetscInt     maxDegree;
863e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
864c330f8ffSToby Isaac     IS           heightIS;
8658c6c5593SMatthew G. Knepley 
866636b9322SMatthew G. Knepley     if (h > minHeight) {
8679566063dSJacob Faibussowitsch       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
868636b9322SMatthew G. Knepley     }
8699566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
8709566063dSJacob Faibussowitsch     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
8719566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
872c330f8ffSToby Isaac     if (pEnd <= pStart) {
8739566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
874c330f8ffSToby Isaac       continue;
875c330f8ffSToby Isaac     }
87647923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
87747923291SMatthew G. Knepley     totDim = 0;
87847923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
8795fedec97SMatthew G. Knepley       PetscBool cohesive;
8805fedec97SMatthew G. Knepley 
881665f567fSMatthew G. Knepley       if (!sp[f]) continue;
8829566063dSJacob Faibussowitsch       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
8839566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
8849c3cf19fSMatthew G. Knepley       totDim += spDim;
8855fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) totDim += spDim;
88647923291SMatthew G. Knepley     }
887636b9322SMatthew G. Knepley     p = lStart < 0 ? pStart : lStart;
8889566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
88963a3b9bcSJacob Faibussowitsch     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);
890c330f8ffSToby Isaac     if (!totDim) {
8919566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
892c330f8ffSToby Isaac       continue;
893c330f8ffSToby Isaac     }
8949566063dSJacob Faibussowitsch     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
895636b9322SMatthew G. Knepley     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
896636b9322SMatthew G. Knepley     if (localU) {
897636b9322SMatthew G. Knepley       PetscInt totDimIn, pIn, numValuesIn;
898636b9322SMatthew G. Knepley 
899636b9322SMatthew G. Knepley       totDimIn = 0;
900636b9322SMatthew G. Knepley       for (f = 0; f < NfIn; ++f) {
9015fedec97SMatthew G. Knepley         PetscBool cohesive;
9025fedec97SMatthew G. Knepley 
903636b9322SMatthew G. Knepley         if (!spIn[f]) continue;
9049566063dSJacob Faibussowitsch         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
9059566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
906636b9322SMatthew G. Knepley         totDimIn += spDim;
9078e6d238bSMatthew G. Knepley         if (isCohesiveIn && !cohesive) totDimIn += spDim;
908636b9322SMatthew G. Knepley       }
9099566063dSJacob Faibussowitsch       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
9109566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
9118e6d238bSMatthew G. Knepley       // TODO We could check that pIn is a cohesive cell for this check
9128e6d238bSMatthew G. Knepley       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);
9139566063dSJacob Faibussowitsch       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
914636b9322SMatthew G. Knepley     }
9159566063dSJacob Faibussowitsch     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
91647923291SMatthew G. Knepley     /* Loop over points at this height */
9179566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
9189566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
919aa0eca99SMatthew G. Knepley     {
920aa0eca99SMatthew G. Knepley       const PetscInt *fields;
921aa0eca99SMatthew G. Knepley 
9229566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
923ad540459SPierre Jolivet       for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE;
924ad540459SPierre Jolivet       for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
9259566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
926aa0eca99SMatthew G. Knepley     }
9278c6c5593SMatthew G. Knepley     if (label) {
9288c6c5593SMatthew G. Knepley       PetscInt i;
92947923291SMatthew G. Knepley 
93047923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
931c330f8ffSToby Isaac         IS              pointIS, isectIS;
93247923291SMatthew G. Knepley         const PetscInt *points;
9338c6c5593SMatthew G. Knepley         PetscInt        n;
934c330f8ffSToby Isaac         PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
935c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
93647923291SMatthew G. Knepley 
9379566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
93847923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
9399566063dSJacob Faibussowitsch         PetscCall(ISIntersect(pointIS, heightIS, &isectIS));
9409566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&pointIS));
941c330f8ffSToby Isaac         if (!isectIS) continue;
9429566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(isectIS, &n));
9439566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(isectIS, &points));
9449566063dSJacob Faibussowitsch         PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree));
94548a46eb9SPierre Jolivet         if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad));
946c330f8ffSToby Isaac         if (!quad) {
947c330f8ffSToby Isaac           if (!h && allPoints) {
948c330f8ffSToby Isaac             quad      = allPoints;
949c330f8ffSToby Isaac             allPoints = NULL;
950c330f8ffSToby Isaac           } else {
9519566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad));
952c330f8ffSToby Isaac           }
953c330f8ffSToby Isaac         }
95479f2dbaeSMatthew G. Knepley         PetscBool computeFaceGeom = htInc && h == minHeight ? PETSC_TRUE : PETSC_FALSE;
95579f2dbaeSMatthew G. Knepley 
95679f2dbaeSMatthew G. Knepley         if (n) {
95779f2dbaeSMatthew G. Knepley           PetscInt depth, dep;
95879f2dbaeSMatthew G. Knepley 
95979f2dbaeSMatthew G. Knepley           PetscCall(DMPlexGetDepth(dm, &depth));
96079f2dbaeSMatthew G. Knepley           PetscCall(DMPlexGetPointDepth(dm, points[0], &dep));
96179f2dbaeSMatthew G. Knepley           if (dep < depth && h == minHeight) computeFaceGeom = PETSC_TRUE;
96279f2dbaeSMatthew G. Knepley         }
96379f2dbaeSMatthew G. Knepley         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, computeFaceGeom, &fegeom));
96447923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
96547923291SMatthew G. Knepley           const PetscInt point = points[p];
96647923291SMatthew G. Knepley 
9679566063dSJacob Faibussowitsch           PetscCall(PetscArrayzero(values, numValues));
9689566063dSJacob Faibussowitsch           PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom));
9699566063dSJacob Faibussowitsch           PetscCall(DMPlexSetActivePoint(dm, point));
970d0609cedSBarry Smith           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));
9719566063dSJacob Faibussowitsch           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
9729566063dSJacob Faibussowitsch           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
97347923291SMatthew G. Knepley         }
9749566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom));
9759566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&fegeom));
9769566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureDestroy(&quad));
9779566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(isectIS, &points));
9789566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isectIS));
97947923291SMatthew G. Knepley       }
9808c6c5593SMatthew G. Knepley     } else {
981c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
982c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
983c330f8ffSToby Isaac       IS              pointIS;
984c330f8ffSToby Isaac 
9859566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS));
9869566063dSJacob Faibussowitsch       PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
98748a46eb9SPierre Jolivet       if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
988c330f8ffSToby Isaac       if (!quad) {
989c330f8ffSToby Isaac         if (!h && allPoints) {
990c330f8ffSToby Isaac           quad      = allPoints;
991c330f8ffSToby Isaac           allPoints = NULL;
992c330f8ffSToby Isaac         } else {
9939566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad));
994c330f8ffSToby Isaac         }
995c330f8ffSToby Isaac       }
9969566063dSJacob Faibussowitsch       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
9978c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
9989566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(values, numValues));
9999566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom));
10009566063dSJacob Faibussowitsch         PetscCall(DMPlexSetActivePoint(dm, p));
1001d0609cedSBarry Smith         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));
10029566063dSJacob Faibussowitsch         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
10039566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
10048c6c5593SMatthew G. Knepley       }
10059566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom));
10069566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomDestroy(&fegeom));
10079566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
10089566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
10098c6c5593SMatthew G. Knepley     }
10109566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&heightIS));
10119566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
10129566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
101347923291SMatthew G. Knepley   }
10148c6c5593SMatthew G. Knepley   /* Cleanup */
1015ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
10169566063dSJacob Faibussowitsch     for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f]));
10179566063dSJacob Faibussowitsch     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
10189566063dSJacob Faibussowitsch     PetscCall(PetscFree2(T, TAux));
10190c898c18SMatthew G. Knepley   }
10209566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&allPoints));
10219566063dSJacob Faibussowitsch   PetscCall(PetscFree3(isFE, sp, spIn));
10229566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
10239566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
10249566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plexIn));
10259566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMDestroy(&plexAux));
10263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
102747923291SMatthew G. Knepley }
10288c6c5593SMatthew G. Knepley 
1029d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
1030d71ae5a4SJacob Faibussowitsch {
10318c6c5593SMatthew G. Knepley   PetscFunctionBegin;
10329566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
10333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10348c6c5593SMatthew G. Knepley }
10358c6c5593SMatthew G. Knepley 
1036d71ae5a4SJacob Faibussowitsch 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)
1037d71ae5a4SJacob Faibussowitsch {
10388c6c5593SMatthew G. Knepley   PetscFunctionBegin;
10399566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
10403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104147923291SMatthew G. Knepley }
104247923291SMatthew G. Knepley 
1043d71ae5a4SJacob Faibussowitsch 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)
1044d71ae5a4SJacob Faibussowitsch {
104547923291SMatthew G. Knepley   PetscFunctionBegin;
10469566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
10473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104847923291SMatthew G. Knepley }
104947923291SMatthew G. Knepley 
1050d71ae5a4SJacob Faibussowitsch 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)
1051d71ae5a4SJacob Faibussowitsch {
10528c6c5593SMatthew G. Knepley   PetscFunctionBegin;
10539566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
10543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105547923291SMatthew G. Knepley }
1056ece3a9fcSMatthew G. Knepley 
1057d71ae5a4SJacob Faibussowitsch 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)
1058d71ae5a4SJacob Faibussowitsch {
1059ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
10609566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX));
10613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1062ece3a9fcSMatthew G. Knepley }
1063