xref: /petsc/src/dm/impls/plex/plexproject.c (revision 989fa6391339cf14e48187f6100e73b78ee6c140)
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 {
72a77a5016SMatthew G. Knepley   PetscInt  debug = ((DM_Plex *)dm->data)->printProject;
73a6e0b375SMatthew G. Knepley   PetscInt  coordDim, Nf, *Nc, f, spDim, d, v, tp;
745fedec97SMatthew G. Knepley   PetscBool isAffine, isCohesive, transform;
758c6c5593SMatthew G. Knepley 
768c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
779566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmIn, &coordDim));
789566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
819566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
828c6c5593SMatthew G. Knepley   /* Get values for closure */
83c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
84c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
858c6c5593SMatthew G. Knepley     void *const ctx = ctxs ? ctxs[f] : NULL;
865fedec97SMatthew G. Knepley     PetscBool   cohesive;
878c6c5593SMatthew G. Knepley 
888c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
899566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
909566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
918c6c5593SMatthew G. Knepley     if (funcs[f]) {
92c330f8ffSToby Isaac       if (isFE[f]) {
93c330f8ffSToby Isaac         PetscQuadrature  allPoints;
94c330f8ffSToby Isaac         PetscInt         q, dim, numPoints;
95c330f8ffSToby Isaac         const PetscReal *points;
96c330f8ffSToby Isaac         PetscScalar     *pointEval;
97c330f8ffSToby Isaac         PetscReal       *x;
98ca3d3a14SMatthew G. Knepley         DM               rdm;
99c330f8ffSToby Isaac 
1009566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDM(sp[f], &rdm));
1019566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
1029566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
1039566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
1049566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x));
105737e23dcSMatthew G. Knepley         PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f]));
106c330f8ffSToby Isaac         for (q = 0; q < numPoints; q++, tp++) {
107c330f8ffSToby Isaac           const PetscReal *v0;
108c330f8ffSToby Isaac 
109c330f8ffSToby Isaac           if (isAffine) {
110665f567fSMatthew G. Knepley             const PetscReal *refpoint    = &points[q * dim];
111665f567fSMatthew G. Knepley             PetscReal        injpoint[3] = {0., 0., 0.};
112665f567fSMatthew G. Knepley 
113665f567fSMatthew G. Knepley             if (dim != fegeom->dim) {
1145fedec97SMatthew G. Knepley               if (isCohesive) {
115665f567fSMatthew G. Knepley                 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
116665f567fSMatthew G. Knepley                 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
117665f567fSMatthew G. Knepley                 refpoint = injpoint;
11863a3b9bcSJacob Faibussowitsch               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim);
119665f567fSMatthew G. Knepley             }
120665f567fSMatthew G. Knepley             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x);
121c330f8ffSToby Isaac             v0 = x;
1228c6c5593SMatthew G. Knepley           } else {
123c330f8ffSToby Isaac             v0 = &fegeom->v[tp * coordDim];
1248c6c5593SMatthew G. Knepley           }
1259371c9d4SSatish Balay           if (transform) {
1269371c9d4SSatish Balay             PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx));
1279371c9d4SSatish Balay             v0 = x;
1289371c9d4SSatish Balay           }
129a77a5016SMatthew G. Knepley           if (debug > 3) {
130a77a5016SMatthew G. Knepley             PetscInt ap;
131a77a5016SMatthew G. Knepley             PetscCall(DMPlexGetActivePoint(dm, &ap));
132a77a5016SMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, "Project point %" PetscInt_FMT ", analytic: ref (", ap));
133a77a5016SMatthew G. Knepley             for (PetscInt d = 0; d < dim; ++d) {
134a77a5016SMatthew G. Knepley               if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
135a77a5016SMatthew G. Knepley               PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)points[q * dim + d]));
136a77a5016SMatthew G. Knepley             }
137a77a5016SMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, ") real ("));
138a77a5016SMatthew G. Knepley             for (PetscInt d = 0; d < dim; ++d) {
139a77a5016SMatthew G. Knepley               if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
140a77a5016SMatthew G. Knepley               PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)v0[d]));
141a77a5016SMatthew G. Knepley             }
142a77a5016SMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
143a77a5016SMatthew G. Knepley           }
1449566063dSJacob Faibussowitsch           PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx));
145c330f8ffSToby Isaac         }
1464bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
1479566063dSJacob Faibussowitsch         PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
1489566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
1499566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x));
1509566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
151c330f8ffSToby Isaac         v += spDim;
1525fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) {
15327f02ce8SMatthew G. Knepley           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
15427f02ce8SMatthew G. Knepley         }
155c330f8ffSToby Isaac       } else {
15648a46eb9SPierre Jolivet         for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
157c330f8ffSToby Isaac       }
158c330f8ffSToby Isaac     } else {
15927f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
1605fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
16127f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
16227f02ce8SMatthew G. Knepley       }
1638c6c5593SMatthew G. Knepley     }
1649c3cf19fSMatthew G. Knepley   }
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1668c6c5593SMatthew G. Knepley }
1678c6c5593SMatthew G. Knepley 
168a6e0b375SMatthew G. Knepley /*
169a6e0b375SMatthew G. Knepley   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
170a6e0b375SMatthew G. Knepley 
171a6e0b375SMatthew G. Knepley   Input Parameters:
172a6e0b375SMatthew G. Knepley + dm             - The output DM
173a6e0b375SMatthew G. Knepley . ds             - The output DS
174a6e0b375SMatthew G. Knepley . dmIn           - The input DM
175a6e0b375SMatthew G. Knepley . dsIn           - The input DS
176a6e0b375SMatthew G. Knepley . dmAux          - The auxiliary DM, which is always for the input space
177a6e0b375SMatthew G. Knepley . dsAux          - The auxiliary DS, which is always for the input space
178a6e0b375SMatthew G. Knepley . time           - The time for this evaluation
179a6e0b375SMatthew G. Knepley . localU         - The local solution
180a6e0b375SMatthew G. Knepley . localA         - The local auziliary fields
181a6e0b375SMatthew G. Knepley . cgeom          - The FE geometry for this point
182a6e0b375SMatthew G. Knepley . sp             - The output PetscDualSpace for each field
183a6e0b375SMatthew G. Knepley . p              - The point in the output DM
184a5b23f4aSJose E. Roman . T              - Input basis and derivatives for each field tabulated on the quadrature points
185ef0bb6c7SMatthew G. Knepley . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
186a6e0b375SMatthew G. Knepley . funcs          - The evaluation function for each field
187a6e0b375SMatthew G. Knepley - ctxs           - The user context for each field
188a6e0b375SMatthew G. Knepley 
189a6e0b375SMatthew G. Knepley   Output Parameter:
190a6e0b375SMatthew G. Knepley . values         - The value for each dual basis vector in the output dual space
191a6e0b375SMatthew G. Knepley 
192a6e0b375SMatthew G. Knepley   Level: developer
193a6e0b375SMatthew G. Knepley 
194a1cb98faSBarry Smith   Note:
195a1cb98faSBarry Smith   Not supported for FV
196a1cb98faSBarry Smith 
197db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()`
198a6e0b375SMatthew G. Knepley */
199d71ae5a4SJacob 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[])
200d71ae5a4SJacob Faibussowitsch {
2018c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
2024bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
2038c6c5593SMatthew G. Knepley   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
2048c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
205191494d9SMatthew G. Knepley   const PetscScalar *constants;
2068c6c5593SMatthew G. Knepley   PetscReal         *x;
2078e6d238bSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
208*989fa639SBrad Aagaard   PetscFEGeom        fegeom, fgeomN[2];
2098e6d238bSMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed, *cone, *ornt;
210ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
2118e6d238bSMatthew G. Knepley   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
2128e6d238bSMatthew G. Knepley   DMPolytopeType     qct;
2138c6c5593SMatthew G. Knepley 
2148c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
2159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
2179566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
2189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
2198e6d238bSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
2209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
2219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
2229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
2239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
2259566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
2269566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmIn, &section));
2279566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
2288e6d238bSMatthew G. Knepley   // Get cohesive cell hanging off face
2298e6d238bSMatthew G. Knepley   if (isCohesiveIn) {
2308e6d238bSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
2318e6d238bSMatthew 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)) {
2328e6d238bSMatthew G. Knepley       DMPolytopeType  ct;
2338e6d238bSMatthew G. Knepley       const PetscInt *support;
2348e6d238bSMatthew G. Knepley       PetscInt        Ns, s;
2358e6d238bSMatthew G. Knepley 
2368e6d238bSMatthew G. Knepley       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
2378e6d238bSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
2388e6d238bSMatthew G. Knepley       for (s = 0; s < Ns; ++s) {
2398e6d238bSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
2408e6d238bSMatthew 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)) {
2418e6d238bSMatthew G. Knepley           inp = support[s];
2428e6d238bSMatthew G. Knepley           break;
2438e6d238bSMatthew G. Knepley         }
2448e6d238bSMatthew G. Knepley       }
2458e6d238bSMatthew G. Knepley       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
2468e6d238bSMatthew G. Knepley       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
2478e6d238bSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
2488e6d238bSMatthew G. Knepley       face[0] = 0;
2498e6d238bSMatthew G. Knepley       face[1] = 0;
2508e6d238bSMatthew G. Knepley     }
2518e6d238bSMatthew G. Knepley   }
252eb8f539aSJed Brown   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
2538c6c5593SMatthew G. Knepley   if (dmAux) {
25444171101SMatthew G. Knepley     PetscInt subp;
2551ba34526SMatthew G. Knepley 
2569566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
2579566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2589566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
2599566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2609566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2619566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
2629566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
2638c6c5593SMatthew G. Knepley   }
2648c6c5593SMatthew G. Knepley   /* Get values for closure */
2654bee2e38SMatthew G. Knepley   isAffine        = cgeom->isAffine;
26627f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
26727f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
268*989fa639SBrad Aagaard   if (isCohesiveIn) {
269*989fa639SBrad Aagaard     fgeomN[0].dim      = cgeom->dim;
270*989fa639SBrad Aagaard     fgeomN[0].dimEmbed = cgeom->dimEmbed;
271*989fa639SBrad Aagaard     fgeomN[1].dim      = cgeom->dim;
272*989fa639SBrad Aagaard     fgeomN[1].dimEmbed = cgeom->dimEmbed;
273*989fa639SBrad Aagaard   }
2744bee2e38SMatthew G. Knepley   if (isAffine) {
2754bee2e38SMatthew G. Knepley     fegeom.v    = x;
2764bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
2774bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
2784bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
2794bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
280*989fa639SBrad Aagaard     if (isCohesiveIn) {
281*989fa639SBrad Aagaard       fgeomN[0].J    = cgeom->suppJ[0];
282*989fa639SBrad Aagaard       fgeomN[0].invJ = cgeom->suppInvJ[0];
283*989fa639SBrad Aagaard       fgeomN[0].detJ = cgeom->suppDetJ[0];
284*989fa639SBrad Aagaard       fgeomN[1].J    = cgeom->suppJ[1];
285*989fa639SBrad Aagaard       fgeomN[1].invJ = cgeom->suppInvJ[1];
286*989fa639SBrad Aagaard       fgeomN[1].detJ = cgeom->suppDetJ[1];
287*989fa639SBrad Aagaard     }
2884bee2e38SMatthew G. Knepley   }
2898c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
290c330f8ffSToby Isaac     PetscQuadrature  allPoints;
291c330f8ffSToby Isaac     PetscInt         q, dim, numPoints;
292c330f8ffSToby Isaac     const PetscReal *points;
293c330f8ffSToby Isaac     PetscScalar     *pointEval;
2945fedec97SMatthew G. Knepley     PetscBool        cohesive;
295c330f8ffSToby Isaac     DM               dm;
296c330f8ffSToby Isaac 
2978c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
2999566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
300c330f8ffSToby Isaac     if (!funcs[f]) {
301be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
3025fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
30327f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
30427f02ce8SMatthew G. Knepley       }
305c330f8ffSToby Isaac       continue;
306c330f8ffSToby Isaac     }
3076f0e67eaSMatthew G. Knepley     const PetscInt ***perms;
3089566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
3096f0e67eaSMatthew G. Knepley     PetscCall(PetscDualSpaceGetSymmetries(sp[f], &perms, NULL));
3109566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
3119566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
3129566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
3130c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
3148e6d238bSMatthew G. Knepley       PetscInt qpt[2];
3158e6d238bSMatthew G. Knepley 
3168e6d238bSMatthew G. Knepley       if (isCohesiveIn) {
3176f0e67eaSMatthew G. Knepley         qpt[0] = perms ? perms[0][ornt[0]][q] : q;
3186f0e67eaSMatthew G. Knepley         qpt[1] = perms ? perms[0][DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0)][q] : q;
3198e6d238bSMatthew G. Knepley       }
320c330f8ffSToby Isaac       if (isAffine) {
321ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x);
3221c531cf8SMatthew G. Knepley       } else {
3234bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp * dE];
3244bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp * dE * dE];
3254bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp * dE * dE];
3264bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
327*989fa639SBrad Aagaard         if (isCohesiveIn) {
328*989fa639SBrad Aagaard           fgeomN[0].J    = &cgeom->suppJ[0][tp * dE * dE];
329*989fa639SBrad Aagaard           fgeomN[0].invJ = &cgeom->suppInvJ[0][tp * dE * dE];
330*989fa639SBrad Aagaard           fgeomN[0].detJ = &cgeom->suppDetJ[0][tp];
331*989fa639SBrad Aagaard           fgeomN[1].J    = &cgeom->suppJ[1][tp * dE * dE];
332*989fa639SBrad Aagaard           fgeomN[1].invJ = &cgeom->suppInvJ[1][tp * dE * dE];
333*989fa639SBrad Aagaard           fgeomN[1].detJ = &cgeom->suppDetJ[1][tp];
334*989fa639SBrad Aagaard         }
3358c6c5593SMatthew G. Knepley       }
3368e6d238bSMatthew G. Knepley       if (coefficients) {
337*989fa639SBrad Aagaard         if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, fgeomN, coefficients, coefficients_t, u, u_x, u_t));
3388e6d238bSMatthew G. Knepley         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
3398e6d238bSMatthew G. Knepley       }
3409566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
3419566063dSJacob Faibussowitsch       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
342a6e0b375SMatthew 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]);
3431c531cf8SMatthew G. Knepley     }
3449566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
3459566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
346c330f8ffSToby Isaac     v += spDim;
34727f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
3485fedec97SMatthew G. Knepley     if (isCohesive && !cohesive) {
34927f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
35027f02ce8SMatthew G. Knepley     }
3518c6c5593SMatthew G. Knepley   }
352eb8f539aSJed Brown   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
3539566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
3548e6d238bSMatthew G. Knepley   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
3553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3568c6c5593SMatthew G. Knepley }
3578c6c5593SMatthew G. Knepley 
35879f2dbaeSMatthew 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[])
359d71ae5a4SJacob Faibussowitsch {
360b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
3614bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
362b18799e0SMatthew G. Knepley   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
363b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
364b18799e0SMatthew G. Knepley   const PetscScalar *constants;
365b18799e0SMatthew G. Knepley   PetscReal         *x;
36679f2dbaeSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
367*989fa639SBrad Aagaard   PetscFEGeom        fegeom, cgeom, fgeomN[2];
36879f2dbaeSMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed, *cone, *ornt;
36979f2dbaeSMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
37079f2dbaeSMatthew G. Knepley   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
37179f2dbaeSMatthew G. Knepley   DMPolytopeType     qct;
372b18799e0SMatthew G. Knepley 
373b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
3749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
37679f2dbaeSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
37779f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
37879f2dbaeSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
37979f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
38079f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
38179f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
38279f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
38379f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
38479f2dbaeSMatthew G. Knepley   PetscCall(DMHasBasisTransform(dmIn, &transform));
38579f2dbaeSMatthew G. Knepley   PetscCall(DMGetLocalSection(dmIn, &section));
38679f2dbaeSMatthew G. Knepley   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
38779f2dbaeSMatthew G. Knepley   // Get cohesive cell hanging off face
38879f2dbaeSMatthew G. Knepley   if (isCohesiveIn) {
38979f2dbaeSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
39079f2dbaeSMatthew 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)) {
39179f2dbaeSMatthew G. Knepley       DMPolytopeType  ct;
39279f2dbaeSMatthew G. Knepley       const PetscInt *support;
39379f2dbaeSMatthew G. Knepley       PetscInt        Ns, s;
39479f2dbaeSMatthew G. Knepley 
39579f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
39679f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
39779f2dbaeSMatthew G. Knepley       for (s = 0; s < Ns; ++s) {
39879f2dbaeSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
39979f2dbaeSMatthew 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)) {
40079f2dbaeSMatthew G. Knepley           inp = support[s];
40179f2dbaeSMatthew G. Knepley           break;
40279f2dbaeSMatthew G. Knepley         }
40379f2dbaeSMatthew G. Knepley       }
40479f2dbaeSMatthew G. Knepley       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
40579f2dbaeSMatthew G. Knepley       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
40679f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
40779f2dbaeSMatthew G. Knepley       face[0] = 0;
40879f2dbaeSMatthew G. Knepley       face[1] = 0;
40979f2dbaeSMatthew G. Knepley     }
41079f2dbaeSMatthew G. Knepley   }
41179f2dbaeSMatthew G. Knepley   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
412b18799e0SMatthew G. Knepley   if (dmAux) {
413a6e0b375SMatthew G. Knepley     PetscInt subp;
414b18799e0SMatthew G. Knepley 
4159566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
4169566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4179566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
4189566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4199566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4209566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
4219566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
422b18799e0SMatthew G. Knepley   }
423b18799e0SMatthew G. Knepley   /* Get values for closure */
4244bee2e38SMatthew G. Knepley   isAffine        = fgeom->isAffine;
425*989fa639SBrad Aagaard   fegeom.dim      = fgeom->dim;
426*989fa639SBrad Aagaard   fegeom.dimEmbed = fgeom->dimEmbed;
427ea78f98cSLisandro Dalcin   fegeom.n        = NULL;
428ea78f98cSLisandro Dalcin   fegeom.J        = NULL;
429ea78f98cSLisandro Dalcin   fegeom.v        = NULL;
430ea78f98cSLisandro Dalcin   fegeom.xi       = NULL;
431a6e0b375SMatthew G. Knepley   cgeom.dim       = fgeom->dim;
432a6e0b375SMatthew G. Knepley   cgeom.dimEmbed  = fgeom->dimEmbed;
433*989fa639SBrad Aagaard   if (isCohesiveIn) {
434*989fa639SBrad Aagaard     fgeomN[0].dim      = fgeom->dim;
435*989fa639SBrad Aagaard     fgeomN[0].dimEmbed = fgeom->dimEmbed;
436*989fa639SBrad Aagaard     fgeomN[1].dim      = fgeom->dim;
437*989fa639SBrad Aagaard     fgeomN[1].dimEmbed = fgeom->dimEmbed;
438*989fa639SBrad Aagaard   }
4394bee2e38SMatthew G. Knepley   if (isAffine) {
4404bee2e38SMatthew G. Knepley     fegeom.v    = x;
4414bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
4424bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
4434bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
4444bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
4454bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
4469f209ee4SMatthew G. Knepley 
4479f209ee4SMatthew G. Knepley     cgeom.J    = fgeom->suppJ[0];
4489f209ee4SMatthew G. Knepley     cgeom.invJ = fgeom->suppInvJ[0];
4499f209ee4SMatthew G. Knepley     cgeom.detJ = fgeom->suppDetJ[0];
450*989fa639SBrad Aagaard 
451*989fa639SBrad Aagaard     if (isCohesiveIn) {
452*989fa639SBrad Aagaard       fgeomN[0].J    = fgeom->suppJ[0];
453*989fa639SBrad Aagaard       fgeomN[0].invJ = fgeom->suppInvJ[0];
454*989fa639SBrad Aagaard       fgeomN[0].detJ = fgeom->suppDetJ[0];
455*989fa639SBrad Aagaard       fgeomN[1].J    = fgeom->suppJ[1];
456*989fa639SBrad Aagaard       fgeomN[1].invJ = fgeom->suppInvJ[1];
457*989fa639SBrad Aagaard       fgeomN[1].detJ = fgeom->suppDetJ[1];
458*989fa639SBrad Aagaard     }
4594bee2e38SMatthew G. Knepley   }
460b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
461b18799e0SMatthew G. Knepley     PetscQuadrature  allPoints;
462b18799e0SMatthew G. Knepley     PetscInt         q, dim, numPoints;
463b18799e0SMatthew G. Knepley     const PetscReal *points;
464b18799e0SMatthew G. Knepley     PetscScalar     *pointEval;
46579f2dbaeSMatthew G. Knepley     PetscBool        cohesive;
466b18799e0SMatthew G. Knepley     DM               dm;
467b18799e0SMatthew G. Knepley 
468b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
46979f2dbaeSMatthew G. Knepley     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
4709566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
471b18799e0SMatthew G. Knepley     if (!funcs[f]) {
472b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
47379f2dbaeSMatthew G. Knepley       if (isCohesive && !cohesive) {
47479f2dbaeSMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
47579f2dbaeSMatthew G. Knepley       }
476b18799e0SMatthew G. Knepley       continue;
477b18799e0SMatthew G. Knepley     }
4789566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
4799566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
4809566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
4819566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
482b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
48379f2dbaeSMatthew G. Knepley       PetscInt qpt[2];
48479f2dbaeSMatthew G. Knepley 
48579f2dbaeSMatthew G. Knepley       if (isCohesiveIn) {
48679f2dbaeSMatthew G. Knepley         // These points are not integration quadratures, but dual space quadratures
487d8b4a066SPierre Jolivet         // If they had multiple points we should match them from both sides, similar to hybrid residual eval
48879f2dbaeSMatthew G. Knepley         qpt[0] = qpt[1] = q;
48979f2dbaeSMatthew G. Knepley       }
490b18799e0SMatthew G. Knepley       if (isAffine) {
491ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x);
492b18799e0SMatthew G. Knepley       } else {
4933fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp * dE];
4949f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp * dE * dE];
4959f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp * dE * dE];
4964bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
4974bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp * dE];
4989f209ee4SMatthew G. Knepley 
4999f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][tp * dE * dE];
5009f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE];
5019f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][tp];
502*989fa639SBrad Aagaard         if (isCohesiveIn) {
503*989fa639SBrad Aagaard           fgeomN[0].J    = &fgeom->suppJ[0][tp * dE * dE];
504*989fa639SBrad Aagaard           fgeomN[0].invJ = &fgeom->suppInvJ[0][tp * dE * dE];
505*989fa639SBrad Aagaard           fgeomN[0].detJ = &fgeom->suppDetJ[0][tp];
506*989fa639SBrad Aagaard           fgeomN[1].J    = &fgeom->suppJ[1][tp * dE * dE];
507*989fa639SBrad Aagaard           fgeomN[1].invJ = &fgeom->suppInvJ[1][tp * dE * dE];
508*989fa639SBrad Aagaard           fgeomN[1].detJ = &fgeom->suppDetJ[1][tp];
509*989fa639SBrad Aagaard         }
510b18799e0SMatthew G. Knepley       }
511a6e0b375SMatthew 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 */
51279f2dbaeSMatthew G. Knepley       if (coefficients) {
513*989fa639SBrad Aagaard         if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, fgeomN, coefficients, coefficients_t, u, u_x, u_t));
51479f2dbaeSMatthew G. Knepley         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
51579f2dbaeSMatthew G. Knepley       }
5169566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
51779f2dbaeSMatthew G. Knepley       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
51879f2dbaeSMatthew 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]);
519b18799e0SMatthew G. Knepley     }
5209566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
5219566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
522b18799e0SMatthew G. Knepley     v += spDim;
52379f2dbaeSMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
52479f2dbaeSMatthew G. Knepley     if (isCohesive && !cohesive) {
52579f2dbaeSMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
526b18799e0SMatthew G. Knepley     }
52779f2dbaeSMatthew G. Knepley   }
52879f2dbaeSMatthew G. Knepley   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
5299566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
53079f2dbaeSMatthew G. Knepley   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
5313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
532b18799e0SMatthew G. Knepley }
533b18799e0SMatthew G. Knepley 
534d71ae5a4SJacob 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[])
535d71ae5a4SJacob Faibussowitsch {
5368c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
5378c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
5388c6c5593SMatthew G. Knepley 
5398c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
5409566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5419566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
5429566063dSJacob Faibussowitsch   if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
5438c6c5593SMatthew G. Knepley   switch (type) {
5448c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
545d71ae5a4SJacob Faibussowitsch   case DM_BC_NATURAL:
546d71ae5a4SJacob 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));
547d71ae5a4SJacob Faibussowitsch     break;
5488c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
5498c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
5509371c9d4SSatish 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));
551d0609cedSBarry Smith     break;
552b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
55379f2dbaeSMatthew 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));
554d0609cedSBarry Smith     break;
555d71ae5a4SJacob Faibussowitsch   default:
556d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type);
5578c6c5593SMatthew G. Knepley   }
5583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5598c6c5593SMatthew G. Knepley }
5608c6c5593SMatthew G. Knepley 
561d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
562d71ae5a4SJacob Faibussowitsch {
563df5c1128SToby Isaac   PetscReal *points;
564df5c1128SToby Isaac   PetscInt   f, numPoints;
565df5c1128SToby Isaac 
566df5c1128SToby Isaac   PetscFunctionBegin;
56719552267SMatthew G. Knepley   if (!dim) {
56819552267SMatthew G. Knepley     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5693ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
57019552267SMatthew G. Knepley   }
571df5c1128SToby Isaac   numPoints = 0;
572df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
573df5c1128SToby Isaac     if (funcs[f]) {
574df5c1128SToby Isaac       PetscQuadrature fAllPoints;
575df5c1128SToby Isaac       PetscInt        fNumPoints;
576df5c1128SToby Isaac 
5779566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5789566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
579df5c1128SToby Isaac       numPoints += fNumPoints;
580df5c1128SToby Isaac     }
581df5c1128SToby Isaac   }
5829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * numPoints, &points));
583df5c1128SToby Isaac   numPoints = 0;
584df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
585df5c1128SToby Isaac     if (funcs[f]) {
586df5c1128SToby Isaac       PetscQuadrature  fAllPoints;
58754ee0cceSMatthew G. Knepley       PetscInt         qdim, fNumPoints, q;
588df5c1128SToby Isaac       const PetscReal *fPoints;
589df5c1128SToby Isaac 
5909566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5919566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
59263a3b9bcSJacob 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);
593df5c1128SToby Isaac       for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q];
594df5c1128SToby Isaac       numPoints += fNumPoints;
595df5c1128SToby Isaac     }
596df5c1128SToby Isaac   }
5979566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5989566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL));
5993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
600df5c1128SToby Isaac }
601df5c1128SToby Isaac 
6025f15299fSJeremy L Thompson /*@C
60320f4b53cSBarry 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`.
6045f15299fSJeremy L Thompson 
6055f15299fSJeremy L Thompson   Input Parameters:
6062fe279fdSBarry Smith + dm     - the `DM`
6072fe279fdSBarry Smith . odm    - the enclosing `DM`
6082fe279fdSBarry Smith . label  - label for `DM` domain, or `NULL` for whole domain
6092fe279fdSBarry Smith . numIds - the number of `ids`
6102fe279fdSBarry Smith . ids    - An array of the label ids in sequence for the domain
6112fe279fdSBarry Smith - height - Height of target cells in `DMPLEX` topology
6125f15299fSJeremy L Thompson 
6135f15299fSJeremy L Thompson   Output Parameters:
6142fe279fdSBarry Smith + point - the first labeled point
615a3b724e8SBarry Smith - ds    - the `PetscDS` corresponding to the first labeled point
6165f15299fSJeremy L Thompson 
6175f15299fSJeremy L Thompson   Level: developer
618a1cb98faSBarry Smith 
6191cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS`
620817da375SSatish Balay @*/
621d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
622d71ae5a4SJacob Faibussowitsch {
623a6e0b375SMatthew G. Knepley   DM              plex;
624a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
625e9f0ba4eSJed Brown   PetscInt        ls = -1;
626e5e52638SMatthew G. Knepley 
627e5e52638SMatthew G. Knepley   PetscFunctionBegin;
6285f15299fSJeremy L Thompson   if (point) *point = -1;
6293ba16761SJacob Faibussowitsch   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6309566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
6319566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
632e9f0ba4eSJed Brown   for (PetscInt i = 0; i < numIds; ++i) {
633e9f0ba4eSJed Brown     IS       labelIS;
634e9f0ba4eSJed Brown     PetscInt num_points, pStart, pEnd;
6359566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
636e9f0ba4eSJed Brown     if (!labelIS) continue; /* No points with that id on this process */
6379566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
6389566063dSJacob Faibussowitsch     PetscCall(ISGetSize(labelIS, &num_points));
639e9f0ba4eSJed Brown     if (num_points) {
640e5e52638SMatthew G. Knepley       const PetscInt *points;
6419566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(labelIS, &points));
642e9f0ba4eSJed Brown       for (PetscInt i = 0; i < num_points; i++) {
643e9f0ba4eSJed Brown         PetscInt point;
6449566063dSJacob Faibussowitsch         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
645e9f0ba4eSJed Brown         if (pStart <= point && point < pEnd) {
646a6e0b375SMatthew G. Knepley           ls = point;
64779f2dbaeSMatthew G. Knepley           if (ds) {
64879f2dbaeSMatthew G. Knepley             // If this is a face of a cohesive cell, then prefer that DS
64979f2dbaeSMatthew G. Knepley             if (height == 1) {
65079f2dbaeSMatthew G. Knepley               const PetscInt *supp;
65179f2dbaeSMatthew G. Knepley               PetscInt        suppSize;
65279f2dbaeSMatthew G. Knepley               DMPolytopeType  ct;
65379f2dbaeSMatthew G. Knepley 
65479f2dbaeSMatthew G. Knepley               PetscCall(DMPlexGetSupport(dm, ls, &supp));
65579f2dbaeSMatthew G. Knepley               PetscCall(DMPlexGetSupportSize(dm, ls, &suppSize));
65679f2dbaeSMatthew G. Knepley               for (PetscInt s = 0; s < suppSize; ++s) {
65779f2dbaeSMatthew G. Knepley                 PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
65879f2dbaeSMatthew 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)) {
65979f2dbaeSMatthew G. Knepley                   ls = supp[s];
66079f2dbaeSMatthew G. Knepley                   break;
66179f2dbaeSMatthew G. Knepley                 }
66279f2dbaeSMatthew G. Knepley               }
66379f2dbaeSMatthew G. Knepley             }
66479f2dbaeSMatthew G. Knepley             PetscCall(DMGetCellDS(dm, ls, ds, NULL));
66579f2dbaeSMatthew G. Knepley           }
6668e6d238bSMatthew G. Knepley           if (ls >= 0) break;
667e5e52638SMatthew G. Knepley         }
668e9f0ba4eSJed Brown       }
6699566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(labelIS, &points));
670e9f0ba4eSJed Brown     }
6719566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&labelIS));
672e5e52638SMatthew G. Knepley     if (ls >= 0) break;
673e5e52638SMatthew G. Knepley   }
6745f15299fSJeremy L Thompson   if (point) *point = ls;
6759566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
6763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
677e5e52638SMatthew G. Knepley }
678e5e52638SMatthew G. Knepley 
6790de51b56SMatthew G. Knepley /*
6800de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
6810de51b56SMatthew G. Knepley 
6820de51b56SMatthew G. Knepley   There are several different scenarios:
6830de51b56SMatthew G. Knepley 
6840de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
6850de51b56SMatthew G. Knepley 
6860de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
6870de51b56SMatthew G. Knepley 
6880de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
6890de51b56SMatthew G. Knepley 
6900de51b56SMatthew 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.
6910de51b56SMatthew G. Knepley 
6920de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
6930de51b56SMatthew G. Knepley 
6940de51b56SMatthew 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.
6950de51b56SMatthew G. Knepley 
696636b9322SMatthew G. Knepley   4) Volumetric input mesh with boundary output mesh
697636b9322SMatthew G. Knepley 
698636b9322SMatthew G. Knepley      Here we must get a subspace for the input DS
699636b9322SMatthew G. Knepley 
7000de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
7010de51b56SMatthew G. Knepley 
7020de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
7030de51b56SMatthew G. Knepley 
7040de51b56SMatthew 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.
7050de51b56SMatthew 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
7065f790a90SMatthew 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
7070de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
7080de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
7090de51b56SMatthew G. Knepley 
7100de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
7110de51b56SMatthew G. Knepley 
7120de51b56SMatthew 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.
7130de51b56SMatthew G. Knepley */
714d71ae5a4SJacob 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)
715d71ae5a4SJacob Faibussowitsch {
716a6e0b375SMatthew G. Knepley   DM               plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
717a6e0b375SMatthew G. Knepley   DMEnclosureType  encIn, encAux;
718a6e0b375SMatthew G. Knepley   PetscDS          ds = NULL, dsIn = NULL, dsAux = NULL;
719ca3d3a14SMatthew G. Knepley   Vec              localA = NULL, tv;
720aa0eca99SMatthew G. Knepley   IS               fieldIS;
72147923291SMatthew G. Knepley   PetscSection     section;
722636b9322SMatthew G. Knepley   PetscDualSpace  *sp, *cellsp, *spIn, *cellspIn;
723ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
7248c6c5593SMatthew G. Knepley   PetscInt        *Nc;
725*989fa639SBrad Aagaard   PetscInt         dim, dimEmbed, depth, pStart, pEnd, lStart = PETSC_DETERMINE, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, minHeightIn, minHeightAux, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
7268e6d238bSMatthew G. Knepley   PetscBool       *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform;
727c330f8ffSToby Isaac   DMField          coordField;
728c330f8ffSToby Isaac   DMLabel          depthLabel;
729c330f8ffSToby Isaac   PetscQuadrature  allPoints = NULL;
73047923291SMatthew G. Knepley 
73147923291SMatthew G. Knepley   PetscFunctionBegin;
7329566063dSJacob Faibussowitsch   if (localU) PetscCall(VecGetDM(localU, &dmIn));
733ad540459SPierre Jolivet   else dmIn = dm;
7349566063dSJacob Faibussowitsch   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
7359371c9d4SSatish Balay   if (localA) PetscCall(VecGetDM(localA, &dmAux));
736ad540459SPierre Jolivet   else dmAux = NULL;
7379566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
7389566063dSJacob Faibussowitsch   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
7399566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
7409566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
7419566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
7439566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
7449566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
7459566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
7460de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
747a6e0b375SMatthew G. Knepley   if (dmAux) {
7489566063dSJacob Faibussowitsch     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
749ad540459SPierre 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");
750636b9322SMatthew G. Knepley   }
7515ee220baSJed Brown   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plexIn, PETSC_TRUE, localU, time, NULL, NULL, NULL));
752e04ae0b4SMatthew G. Knepley   PetscCall(DMGetCoordinateField(dm, &coordField));
753e44f6aebSMatthew G. Knepley   PetscCheck(coordField, PETSC_COMM_SELF, PETSC_ERR_USER, "DM must have a coordinate field");
754e04ae0b4SMatthew G. Knepley   /**** No collective calls below this point ****/
755636b9322SMatthew G. Knepley   /* Determine height for iteration of all meshes */
756636b9322SMatthew G. Knepley   {
757636b9322SMatthew G. Knepley     DMPolytopeType ct, ctIn, ctAux;
758*989fa639SBrad Aagaard     PetscInt       p, pStartIn, pStartAux, pEndAux;
75919552267SMatthew G. Knepley     PetscInt       dim = -1, dimIn = -1, dimAux = -1;
76088aca1feSMatthew G. Knepley 
7619566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
762636b9322SMatthew G. Knepley     if (pEnd > pStart) {
7639566063dSJacob Faibussowitsch       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
764636b9322SMatthew G. Knepley       p = lStart < 0 ? pStart : lStart;
7659566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plex, p, &ct));
766636b9322SMatthew G. Knepley       dim = DMPolytopeTypeGetDim(ct);
7679566063dSJacob Faibussowitsch       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
7689566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
7699566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
770636b9322SMatthew G. Knepley       dimIn = DMPolytopeTypeGetDim(ctIn);
771636b9322SMatthew G. Knepley       if (dmAux) {
7729566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
77319552267SMatthew G. Knepley         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
77419552267SMatthew G. Knepley         if (pStartAux < pEndAux) {
7759566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
776636b9322SMatthew G. Knepley           dimAux = DMPolytopeTypeGetDim(ctAux);
77719552267SMatthew G. Knepley         }
778636b9322SMatthew G. Knepley       } else dimAux = dim;
779e04ae0b4SMatthew G. Knepley     } else {
780e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plex));
781e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plexIn));
782e04ae0b4SMatthew G. Knepley       if (dmAux) PetscCall(DMDestroy(&plexAux));
7833ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
784e39fcb4eSStefano Zampini     }
785636b9322SMatthew G. Knepley     if (dim < 0) {
786636b9322SMatthew G. Knepley       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
787636b9322SMatthew G. Knepley 
788636b9322SMatthew G. Knepley       /* Fall back to determination based on being a submesh */
7899566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plex, &spmap));
7909566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
7919566063dSJacob Faibussowitsch       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
792636b9322SMatthew G. Knepley       dim    = spmap ? 1 : 0;
793636b9322SMatthew G. Knepley       dimIn  = spmapIn ? 1 : 0;
794636b9322SMatthew G. Knepley       dimAux = spmapAux ? 1 : 0;
79588aca1feSMatthew G. Knepley     }
796636b9322SMatthew G. Knepley     {
79757508eceSPierre Jolivet       PetscInt dimProj   = PetscMin(PetscMin(dim, dimIn), dimAux < 0 ? PETSC_INT_MAX : dimAux);
79819552267SMatthew G. Knepley       PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;
799636b9322SMatthew G. Knepley 
80019552267SMatthew 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");
801636b9322SMatthew G. Knepley       if (dimProj < dim) minHeight = 1;
802636b9322SMatthew G. Knepley       htInc    = dim - dimProj;
803636b9322SMatthew G. Knepley       htIncIn  = dimIn - dimProj;
80419552267SMatthew G. Knepley       htIncAux = dimAuxEff - dimProj;
8050de51b56SMatthew G. Knepley     }
806a6e0b375SMatthew G. Knepley   }
8079566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
8089566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
8099566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
8102af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
81163a3b9bcSJacob 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);
81279f2dbaeSMatthew G. Knepley   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, NULL, &ds));
8139566063dSJacob Faibussowitsch   if (!ds) PetscCall(DMGetDS(dm, &ds));
81479f2dbaeSMatthew G. Knepley   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, minHeight, NULL, &dsIn));
8159566063dSJacob Faibussowitsch   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
816*989fa639SBrad Aagaard   if (!dsIn) {
817*989fa639SBrad Aagaard     if (encIn == DM_ENC_SUPERMESH) {
818*989fa639SBrad Aagaard       PetscInt p = pStart, pIn;
819*989fa639SBrad Aagaard 
820*989fa639SBrad Aagaard       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? p : lStart, &pIn));
821*989fa639SBrad Aagaard       // If the input mesh is higher dimensional than the output mesh, get a cell from the output mesh
822*989fa639SBrad Aagaard       if (htIncIn) PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &p, NULL));
823*989fa639SBrad Aagaard       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? p : lStart, &pIn));
824*989fa639SBrad Aagaard       PetscCall(DMGetCellDS(dmIn, pIn, &dsIn, NULL));
825*989fa639SBrad Aagaard     } else PetscCall(DMGetDS(dmIn, &dsIn));
826*989fa639SBrad Aagaard   }
8279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
8289566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
8298e6d238bSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
8308e6d238bSMatthew G. Knepley   if (isCohesiveIn) --htIncIn; // Should be rearranged
8319566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &NfTot));
8329566063dSJacob Faibussowitsch   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
83307218a29SMatthew G. Knepley   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL));
8349566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
8359566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
8369566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
8370c898c18SMatthew G. Knepley   if (dmAux) {
838*989fa639SBrad Aagaard     if (encAux == DM_ENC_SUPERMESH) {
839*989fa639SBrad Aagaard       PetscInt p = pStart, pAux;
840*989fa639SBrad Aagaard 
841*989fa639SBrad Aagaard       // If the auxiliary mesh is higher dimensional than the output mesh, get a cell from the output mesh
842*989fa639SBrad Aagaard       if (htIncAux) PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &p, NULL));
843*989fa639SBrad Aagaard       PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, lStart < 0 ? p : lStart, &pAux));
844*989fa639SBrad Aagaard       PetscCall(DMGetCellDS(dmAux, pAux, &dsAux, NULL));
845*989fa639SBrad Aagaard     } else PetscCall(DMGetDS(dmAux, &dsAux));
8469566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
8470c898c18SMatthew G. Knepley   }
8489566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
8499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
8509566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
8519371c9d4SSatish Balay   else {
8529371c9d4SSatish Balay     cellsp   = sp;
8539371c9d4SSatish Balay     cellspIn = spIn;
8549371c9d4SSatish Balay   }
8558c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
85647923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
857665f567fSMatthew G. Knepley     PetscDiscType disctype;
85847923291SMatthew G. Knepley 
8599566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
860665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
861665f567fSMatthew G. Knepley       PetscFE fe;
86247923291SMatthew G. Knepley 
86347923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
864665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
8659566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
8669566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
867665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
868665f567fSMatthew G. Knepley       PetscFV fv;
86947923291SMatthew G. Knepley 
87047923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
871665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
8729566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv));
8739566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
874665f567fSMatthew G. Knepley     } else {
875665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
876665f567fSMatthew G. Knepley       cellsp[f] = NULL;
877665f567fSMatthew G. Knepley     }
87847923291SMatthew G. Knepley   }
879636b9322SMatthew G. Knepley   for (f = 0; f < NfIn; ++f) {
880636b9322SMatthew G. Knepley     PetscDiscType disctype;
881636b9322SMatthew G. Knepley 
8829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
883636b9322SMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
884636b9322SMatthew G. Knepley       PetscFE fe;
885636b9322SMatthew G. Knepley 
8869566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
8879566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
888636b9322SMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
889636b9322SMatthew G. Knepley       PetscFV fv;
890636b9322SMatthew G. Knepley 
8919566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv));
8929566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
893636b9322SMatthew G. Knepley     } else {
894636b9322SMatthew G. Knepley       cellspIn[f] = NULL;
895636b9322SMatthew G. Knepley     }
896636b9322SMatthew G. Knepley   }
897636b9322SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
8989371c9d4SSatish Balay     if (!htInc) {
8999371c9d4SSatish Balay       sp[f] = cellsp[f];
9009371c9d4SSatish Balay     } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
901636b9322SMatthew G. Knepley   }
902ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
90374fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
904665f567fSMatthew G. Knepley     PetscDiscType    disctype;
9054a3e9fdbSToby Isaac     const PetscReal *points;
906409bb69aSDarsh Nathawani     PetscInt         numPoints, k;
9070c898c18SMatthew G. Knepley 
90808401ef6SPierre Jolivet     PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
9099566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints));
9109566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
9119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
912a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
9139371c9d4SSatish Balay       if (!htIncIn) {
9149371c9d4SSatish Balay         spIn[f] = cellspIn[f];
9159371c9d4SSatish Balay       } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
916636b9322SMatthew G. Knepley 
9179566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
918665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
9199566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem));
9209371c9d4SSatish Balay       if (!htIncIn) {
9219371c9d4SSatish Balay         subfem = fem;
9229371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
923409bb69aSDarsh Nathawani       PetscCall(PetscDSGetJetDegree(dsIn, f, &k));
924409bb69aSDarsh Nathawani       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &T[f]));
9250c898c18SMatthew G. Knepley     }
9260c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
9279566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
928665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
9299566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem));
9309371c9d4SSatish Balay       if (!htIncAux) {
9319371c9d4SSatish Balay         subfem = fem;
9329371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
933409bb69aSDarsh Nathawani       PetscCall(PetscDSGetJetDegree(dsAux, f, &k));
934409bb69aSDarsh Nathawani       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &TAux[f]));
9350c898c18SMatthew G. Knepley     }
9360c898c18SMatthew G. Knepley   }
93747923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
9382af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
939636b9322SMatthew G. Knepley     PetscInt     hEff     = h - minHeight + htInc;
940636b9322SMatthew G. Knepley     PetscInt     hEffIn   = h - minHeight + htIncIn;
941636b9322SMatthew G. Knepley     PetscInt     hEffAux  = h - minHeight + htIncAux;
942a6e0b375SMatthew G. Knepley     PetscDS      dsEff    = ds;
943636b9322SMatthew G. Knepley     PetscDS      dsEffIn  = dsIn;
944636b9322SMatthew G. Knepley     PetscDS      dsEffAux = dsAux;
9458c6c5593SMatthew G. Knepley     PetscScalar *values;
946b7260050SToby Isaac     PetscBool   *fieldActive;
947b7260050SToby Isaac     PetscInt     maxDegree;
948*989fa639SBrad Aagaard     PetscInt     p, spDim, totDim, numValues;
949c330f8ffSToby Isaac     IS           heightIS;
9508c6c5593SMatthew G. Knepley 
951636b9322SMatthew G. Knepley     if (h > minHeight) {
9529566063dSJacob Faibussowitsch       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
953636b9322SMatthew G. Knepley     }
9549566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
9559566063dSJacob Faibussowitsch     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
9569566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
957c330f8ffSToby Isaac     if (pEnd <= pStart) {
9589566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
959c330f8ffSToby Isaac       continue;
960c330f8ffSToby Isaac     }
96147923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
96247923291SMatthew G. Knepley     totDim = 0;
96347923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
9645fedec97SMatthew G. Knepley       PetscBool cohesive;
9655fedec97SMatthew G. Knepley 
966665f567fSMatthew G. Knepley       if (!sp[f]) continue;
9679566063dSJacob Faibussowitsch       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
9689566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
9699c3cf19fSMatthew G. Knepley       totDim += spDim;
9705fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) totDim += spDim;
97147923291SMatthew G. Knepley     }
972636b9322SMatthew G. Knepley     p = lStart < 0 ? pStart : lStart;
9739566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
97463a3b9bcSJacob 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);
975c330f8ffSToby Isaac     if (!totDim) {
9769566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
977c330f8ffSToby Isaac       continue;
978c330f8ffSToby Isaac     }
9799566063dSJacob Faibussowitsch     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
980636b9322SMatthew G. Knepley     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
981636b9322SMatthew G. Knepley     if (localU) {
982636b9322SMatthew G. Knepley       PetscInt totDimIn, pIn, numValuesIn;
983636b9322SMatthew G. Knepley 
984636b9322SMatthew G. Knepley       totDimIn = 0;
985636b9322SMatthew G. Knepley       for (f = 0; f < NfIn; ++f) {
9865fedec97SMatthew G. Knepley         PetscBool cohesive;
9875fedec97SMatthew G. Knepley 
988636b9322SMatthew G. Knepley         if (!spIn[f]) continue;
9899566063dSJacob Faibussowitsch         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
9909566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
991636b9322SMatthew G. Knepley         totDimIn += spDim;
9928e6d238bSMatthew G. Knepley         if (isCohesiveIn && !cohesive) totDimIn += spDim;
993636b9322SMatthew G. Knepley       }
9949566063dSJacob Faibussowitsch       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
9959566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
9968e6d238bSMatthew G. Knepley       // TODO We could check that pIn is a cohesive cell for this check
9978e6d238bSMatthew 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);
9989566063dSJacob Faibussowitsch       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
999636b9322SMatthew G. Knepley     }
10009566063dSJacob Faibussowitsch     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
100147923291SMatthew G. Knepley     /* Loop over points at this height */
10029566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
10039566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
1004aa0eca99SMatthew G. Knepley     {
1005aa0eca99SMatthew G. Knepley       const PetscInt *fields;
1006aa0eca99SMatthew G. Knepley 
10079566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
1008ad540459SPierre Jolivet       for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE;
1009ad540459SPierre Jolivet       for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
10109566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
1011aa0eca99SMatthew G. Knepley     }
10128c6c5593SMatthew G. Knepley     if (label) {
10138c6c5593SMatthew G. Knepley       PetscInt i;
101447923291SMatthew G. Knepley 
101547923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
1016c330f8ffSToby Isaac         IS              pointIS, isectIS;
101747923291SMatthew G. Knepley         const PetscInt *points;
10188c6c5593SMatthew G. Knepley         PetscInt        n;
1019c330f8ffSToby Isaac         PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
1020c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
102147923291SMatthew G. Knepley 
10229566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
102347923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
10249566063dSJacob Faibussowitsch         PetscCall(ISIntersect(pointIS, heightIS, &isectIS));
10259566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&pointIS));
1026c330f8ffSToby Isaac         if (!isectIS) continue;
10279566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(isectIS, &n));
10289566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(isectIS, &points));
10299566063dSJacob Faibussowitsch         PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree));
103048a46eb9SPierre Jolivet         if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad));
1031c330f8ffSToby Isaac         if (!quad) {
1032c330f8ffSToby Isaac           if (!h && allPoints) {
1033c330f8ffSToby Isaac             quad      = allPoints;
1034c330f8ffSToby Isaac             allPoints = NULL;
1035c330f8ffSToby Isaac           } else {
10369566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad));
1037c330f8ffSToby Isaac           }
1038c330f8ffSToby Isaac         }
1039ac9d17c7SMatthew G. Knepley         PetscFEGeomMode geommode = htInc && h == minHeight ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC;
104079f2dbaeSMatthew G. Knepley 
104179f2dbaeSMatthew G. Knepley         if (n) {
104279f2dbaeSMatthew G. Knepley           PetscInt depth, dep;
104379f2dbaeSMatthew G. Knepley 
104479f2dbaeSMatthew G. Knepley           PetscCall(DMPlexGetDepth(dm, &depth));
104579f2dbaeSMatthew G. Knepley           PetscCall(DMPlexGetPointDepth(dm, points[0], &dep));
1046ac9d17c7SMatthew G. Knepley           if (dep < depth && h == minHeight) geommode = PETSC_FEGEOM_BOUNDARY;
104779f2dbaeSMatthew G. Knepley         }
1048ac9d17c7SMatthew G. Knepley         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, geommode, &fegeom));
104947923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
105047923291SMatthew G. Knepley           const PetscInt point = points[p];
105147923291SMatthew G. Knepley 
10529566063dSJacob Faibussowitsch           PetscCall(PetscArrayzero(values, numValues));
10539566063dSJacob Faibussowitsch           PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom));
10549566063dSJacob Faibussowitsch           PetscCall(DMPlexSetActivePoint(dm, point));
1055d0609cedSBarry 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));
10569566063dSJacob Faibussowitsch           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
10579566063dSJacob Faibussowitsch           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
105847923291SMatthew G. Knepley         }
10599566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom));
10609566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&fegeom));
10619566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureDestroy(&quad));
10629566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(isectIS, &points));
10639566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isectIS));
106447923291SMatthew G. Knepley       }
10658c6c5593SMatthew G. Knepley     } else {
1066c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
1067c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
1068c330f8ffSToby Isaac       IS              pointIS;
1069c330f8ffSToby Isaac 
10709566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS));
10719566063dSJacob Faibussowitsch       PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
107248a46eb9SPierre Jolivet       if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
1073c330f8ffSToby Isaac       if (!quad) {
1074c330f8ffSToby Isaac         if (!h && allPoints) {
1075c330f8ffSToby Isaac           quad      = allPoints;
1076c330f8ffSToby Isaac           allPoints = NULL;
1077c330f8ffSToby Isaac         } else {
10789566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad));
1079c330f8ffSToby Isaac         }
1080c330f8ffSToby Isaac       }
1081ac9d17c7SMatthew G. Knepley       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC, &fegeom));
10828c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
10839566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(values, numValues));
10849566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom));
10859566063dSJacob Faibussowitsch         PetscCall(DMPlexSetActivePoint(dm, p));
1086d0609cedSBarry 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));
10879566063dSJacob Faibussowitsch         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
10889566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
10898c6c5593SMatthew G. Knepley       }
10909566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom));
10919566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomDestroy(&fegeom));
10929566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
10939566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
10948c6c5593SMatthew G. Knepley     }
10959566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&heightIS));
10969566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
10979566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
109847923291SMatthew G. Knepley   }
10998c6c5593SMatthew G. Knepley   /* Cleanup */
1100ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
11019566063dSJacob Faibussowitsch     for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f]));
11029566063dSJacob Faibussowitsch     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
11039566063dSJacob Faibussowitsch     PetscCall(PetscFree2(T, TAux));
11040c898c18SMatthew G. Knepley   }
11059566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&allPoints));
11069566063dSJacob Faibussowitsch   PetscCall(PetscFree3(isFE, sp, spIn));
11079566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
11089566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
11099566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plexIn));
11109566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMDestroy(&plexAux));
11113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
111247923291SMatthew G. Knepley }
11138c6c5593SMatthew G. Knepley 
1114d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
1115d71ae5a4SJacob Faibussowitsch {
11168c6c5593SMatthew G. Knepley   PetscFunctionBegin;
11179566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
11183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11198c6c5593SMatthew G. Knepley }
11208c6c5593SMatthew G. Knepley 
1121d71ae5a4SJacob 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)
1122d71ae5a4SJacob Faibussowitsch {
11238c6c5593SMatthew G. Knepley   PetscFunctionBegin;
11249566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
11253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112647923291SMatthew G. Knepley }
112747923291SMatthew G. Knepley 
1128d71ae5a4SJacob 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)
1129d71ae5a4SJacob Faibussowitsch {
113047923291SMatthew G. Knepley   PetscFunctionBegin;
11319566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
11323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113347923291SMatthew G. Knepley }
113447923291SMatthew G. Knepley 
1135d71ae5a4SJacob 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)
1136d71ae5a4SJacob Faibussowitsch {
11378c6c5593SMatthew G. Knepley   PetscFunctionBegin;
11389566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
11393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114047923291SMatthew G. Knepley }
1141ece3a9fcSMatthew G. Knepley 
1142d71ae5a4SJacob 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)
1143d71ae5a4SJacob Faibussowitsch {
1144ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
11459566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX));
11463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1147ece3a9fcSMatthew G. Knepley }
1148