xref: /petsc/src/dm/impls/plex/plexproject.c (revision ac9d17c75cb5e991458e9675ca6ec07d361dad47)
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];
2084bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
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;
2684bee2e38SMatthew G. Knepley   if (isAffine) {
2694bee2e38SMatthew G. Knepley     fegeom.v    = x;
2704bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
2714bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
2724bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
2734bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
2744bee2e38SMatthew G. Knepley   }
2758c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
276c330f8ffSToby Isaac     PetscQuadrature  allPoints;
277c330f8ffSToby Isaac     PetscInt         q, dim, numPoints;
278c330f8ffSToby Isaac     const PetscReal *points;
279c330f8ffSToby Isaac     PetscScalar     *pointEval;
2805fedec97SMatthew G. Knepley     PetscBool        cohesive;
281c330f8ffSToby Isaac     DM               dm;
282c330f8ffSToby Isaac 
2838c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
2859566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
286c330f8ffSToby Isaac     if (!funcs[f]) {
287be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
2885fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
28927f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
29027f02ce8SMatthew G. Knepley       }
291c330f8ffSToby Isaac       continue;
292c330f8ffSToby Isaac     }
2936f0e67eaSMatthew G. Knepley     const PetscInt ***perms;
2949566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
2956f0e67eaSMatthew G. Knepley     PetscCall(PetscDualSpaceGetSymmetries(sp[f], &perms, NULL));
2969566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
2979566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
2989566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
2990c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
3008e6d238bSMatthew G. Knepley       PetscInt qpt[2];
3018e6d238bSMatthew G. Knepley 
3028e6d238bSMatthew G. Knepley       if (isCohesiveIn) {
3036f0e67eaSMatthew G. Knepley         qpt[0] = perms ? perms[0][ornt[0]][q] : q;
3046f0e67eaSMatthew G. Knepley         qpt[1] = perms ? perms[0][DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0)][q] : q;
3058e6d238bSMatthew G. Knepley       }
306c330f8ffSToby Isaac       if (isAffine) {
307ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x);
3081c531cf8SMatthew G. Knepley       } else {
3094bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp * dE];
3104bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp * dE * dE];
3114bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp * dE * dE];
3124bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
3138c6c5593SMatthew G. Knepley       }
3148e6d238bSMatthew G. Knepley       if (coefficients) {
3158e6d238bSMatthew 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));
3168e6d238bSMatthew G. Knepley         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
3178e6d238bSMatthew G. Knepley       }
3189566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
3199566063dSJacob Faibussowitsch       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
320a6e0b375SMatthew 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]);
3211c531cf8SMatthew G. Knepley     }
3229566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
3239566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
324c330f8ffSToby Isaac     v += spDim;
32527f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
3265fedec97SMatthew G. Knepley     if (isCohesive && !cohesive) {
32727f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
32827f02ce8SMatthew G. Knepley     }
3298c6c5593SMatthew G. Knepley   }
330eb8f539aSJed Brown   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
3319566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
3328e6d238bSMatthew G. Knepley   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3348c6c5593SMatthew G. Knepley }
3358c6c5593SMatthew G. Knepley 
33679f2dbaeSMatthew 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[])
337d71ae5a4SJacob Faibussowitsch {
338b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
3394bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
340b18799e0SMatthew G. Knepley   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
341b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
342b18799e0SMatthew G. Knepley   const PetscScalar *constants;
343b18799e0SMatthew G. Knepley   PetscReal         *x;
34479f2dbaeSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
3459f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
34679f2dbaeSMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed, *cone, *ornt;
34779f2dbaeSMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
34879f2dbaeSMatthew G. Knepley   PetscBool          isAffine, isCohesive, isCohesiveIn, transform;
34979f2dbaeSMatthew G. Knepley   DMPolytopeType     qct;
350b18799e0SMatthew G. Knepley 
351b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
3529566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3539566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
35479f2dbaeSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
35579f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
35679f2dbaeSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
35779f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
35879f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
35979f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
36079f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
36179f2dbaeSMatthew G. Knepley   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
36279f2dbaeSMatthew G. Knepley   PetscCall(DMHasBasisTransform(dmIn, &transform));
36379f2dbaeSMatthew G. Knepley   PetscCall(DMGetLocalSection(dmIn, &section));
36479f2dbaeSMatthew G. Knepley   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
36579f2dbaeSMatthew G. Knepley   // Get cohesive cell hanging off face
36679f2dbaeSMatthew G. Knepley   if (isCohesiveIn) {
36779f2dbaeSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
36879f2dbaeSMatthew 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)) {
36979f2dbaeSMatthew G. Knepley       DMPolytopeType  ct;
37079f2dbaeSMatthew G. Knepley       const PetscInt *support;
37179f2dbaeSMatthew G. Knepley       PetscInt        Ns, s;
37279f2dbaeSMatthew G. Knepley 
37379f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetSupport(dmIn, inp, &support));
37479f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
37579f2dbaeSMatthew G. Knepley       for (s = 0; s < Ns; ++s) {
37679f2dbaeSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
37779f2dbaeSMatthew 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)) {
37879f2dbaeSMatthew G. Knepley           inp = support[s];
37979f2dbaeSMatthew G. Knepley           break;
38079f2dbaeSMatthew G. Knepley         }
38179f2dbaeSMatthew G. Knepley       }
38279f2dbaeSMatthew G. Knepley       PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
38379f2dbaeSMatthew G. Knepley       PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
38479f2dbaeSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
38579f2dbaeSMatthew G. Knepley       face[0] = 0;
38679f2dbaeSMatthew G. Knepley       face[1] = 0;
38779f2dbaeSMatthew G. Knepley     }
38879f2dbaeSMatthew G. Knepley   }
38979f2dbaeSMatthew G. Knepley   if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
390b18799e0SMatthew G. Knepley   if (dmAux) {
391a6e0b375SMatthew G. Knepley     PetscInt subp;
392b18799e0SMatthew G. Knepley 
3939566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
3949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3959566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
3969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
3999566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
400b18799e0SMatthew G. Knepley   }
401b18799e0SMatthew G. Knepley   /* Get values for closure */
4024bee2e38SMatthew G. Knepley   isAffine       = fgeom->isAffine;
403ea78f98cSLisandro Dalcin   fegeom.n       = NULL;
404ea78f98cSLisandro Dalcin   fegeom.J       = NULL;
405ea78f98cSLisandro Dalcin   fegeom.v       = NULL;
406ea78f98cSLisandro Dalcin   fegeom.xi      = NULL;
407a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
408a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
4094bee2e38SMatthew G. Knepley   if (isAffine) {
4104bee2e38SMatthew G. Knepley     fegeom.v    = x;
4114bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
4124bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
4134bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
4144bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
4154bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
4169f209ee4SMatthew G. Knepley 
4179f209ee4SMatthew G. Knepley     cgeom.J    = fgeom->suppJ[0];
4189f209ee4SMatthew G. Knepley     cgeom.invJ = fgeom->suppInvJ[0];
4199f209ee4SMatthew G. Knepley     cgeom.detJ = fgeom->suppDetJ[0];
4204bee2e38SMatthew G. Knepley   }
421b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
422b18799e0SMatthew G. Knepley     PetscQuadrature  allPoints;
423b18799e0SMatthew G. Knepley     PetscInt         q, dim, numPoints;
424b18799e0SMatthew G. Knepley     const PetscReal *points;
425b18799e0SMatthew G. Knepley     PetscScalar     *pointEval;
42679f2dbaeSMatthew G. Knepley     PetscBool        cohesive;
427b18799e0SMatthew G. Knepley     DM               dm;
428b18799e0SMatthew G. Knepley 
429b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
43079f2dbaeSMatthew G. Knepley     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
4319566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
432b18799e0SMatthew G. Knepley     if (!funcs[f]) {
433b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
43479f2dbaeSMatthew G. Knepley       if (isCohesive && !cohesive) {
43579f2dbaeSMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
43679f2dbaeSMatthew G. Knepley       }
437b18799e0SMatthew G. Knepley       continue;
438b18799e0SMatthew G. Knepley     }
4399566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
4409566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
4419566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
4429566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
443b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
44479f2dbaeSMatthew G. Knepley       PetscInt qpt[2];
44579f2dbaeSMatthew G. Knepley 
44679f2dbaeSMatthew G. Knepley       if (isCohesiveIn) {
44779f2dbaeSMatthew G. Knepley         // These points are not integration quadratures, but dual space quadratures
448d8b4a066SPierre Jolivet         // If they had multiple points we should match them from both sides, similar to hybrid residual eval
44979f2dbaeSMatthew G. Knepley         qpt[0] = qpt[1] = q;
45079f2dbaeSMatthew G. Knepley       }
451b18799e0SMatthew G. Knepley       if (isAffine) {
452ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x);
453b18799e0SMatthew G. Knepley       } else {
4543fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp * dE];
4559f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp * dE * dE];
4569f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp * dE * dE];
4574bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
4584bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp * dE];
4599f209ee4SMatthew G. Knepley 
4609f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][tp * dE * dE];
4619f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE];
4629f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][tp];
463b18799e0SMatthew G. Knepley       }
464a6e0b375SMatthew 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 */
46579f2dbaeSMatthew G. Knepley       if (coefficients) {
46679f2dbaeSMatthew 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));
46779f2dbaeSMatthew G. Knepley         else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
46879f2dbaeSMatthew G. Knepley       }
4699566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
47079f2dbaeSMatthew G. Knepley       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
47179f2dbaeSMatthew 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]);
472b18799e0SMatthew G. Knepley     }
4739566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
4749566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
475b18799e0SMatthew G. Knepley     v += spDim;
47679f2dbaeSMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
47779f2dbaeSMatthew G. Knepley     if (isCohesive && !cohesive) {
47879f2dbaeSMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
479b18799e0SMatthew G. Knepley     }
48079f2dbaeSMatthew G. Knepley   }
48179f2dbaeSMatthew G. Knepley   if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
4829566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
48379f2dbaeSMatthew G. Knepley   if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
485b18799e0SMatthew G. Knepley }
486b18799e0SMatthew G. Knepley 
487d71ae5a4SJacob 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[])
488d71ae5a4SJacob Faibussowitsch {
4898c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
4908c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
4918c6c5593SMatthew G. Knepley 
4928c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
4939566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4949566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
4959566063dSJacob Faibussowitsch   if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
4968c6c5593SMatthew G. Knepley   switch (type) {
4978c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
498d71ae5a4SJacob Faibussowitsch   case DM_BC_NATURAL:
499d71ae5a4SJacob 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));
500d71ae5a4SJacob Faibussowitsch     break;
5018c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
5028c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
5039371c9d4SSatish 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));
504d0609cedSBarry Smith     break;
505b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
50679f2dbaeSMatthew 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));
507d0609cedSBarry Smith     break;
508d71ae5a4SJacob Faibussowitsch   default:
509d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type);
5108c6c5593SMatthew G. Knepley   }
5113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5128c6c5593SMatthew G. Knepley }
5138c6c5593SMatthew G. Knepley 
514d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
515d71ae5a4SJacob Faibussowitsch {
516df5c1128SToby Isaac   PetscReal *points;
517df5c1128SToby Isaac   PetscInt   f, numPoints;
518df5c1128SToby Isaac 
519df5c1128SToby Isaac   PetscFunctionBegin;
52019552267SMatthew G. Knepley   if (!dim) {
52119552267SMatthew G. Knepley     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5223ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
52319552267SMatthew G. Knepley   }
524df5c1128SToby Isaac   numPoints = 0;
525df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
526df5c1128SToby Isaac     if (funcs[f]) {
527df5c1128SToby Isaac       PetscQuadrature fAllPoints;
528df5c1128SToby Isaac       PetscInt        fNumPoints;
529df5c1128SToby Isaac 
5309566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5319566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
532df5c1128SToby Isaac       numPoints += fNumPoints;
533df5c1128SToby Isaac     }
534df5c1128SToby Isaac   }
5359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * numPoints, &points));
536df5c1128SToby Isaac   numPoints = 0;
537df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
538df5c1128SToby Isaac     if (funcs[f]) {
539df5c1128SToby Isaac       PetscQuadrature  fAllPoints;
54054ee0cceSMatthew G. Knepley       PetscInt         qdim, fNumPoints, q;
541df5c1128SToby Isaac       const PetscReal *fPoints;
542df5c1128SToby Isaac 
5439566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5449566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
54563a3b9bcSJacob 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);
546df5c1128SToby Isaac       for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q];
547df5c1128SToby Isaac       numPoints += fNumPoints;
548df5c1128SToby Isaac     }
549df5c1128SToby Isaac   }
5509566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5519566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL));
5523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
553df5c1128SToby Isaac }
554df5c1128SToby Isaac 
5555f15299fSJeremy L Thompson /*@C
55620f4b53cSBarry 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`.
5575f15299fSJeremy L Thompson 
5585f15299fSJeremy L Thompson   Input Parameters:
5592fe279fdSBarry Smith + dm     - the `DM`
5602fe279fdSBarry Smith . odm    - the enclosing `DM`
5612fe279fdSBarry Smith . label  - label for `DM` domain, or `NULL` for whole domain
5622fe279fdSBarry Smith . numIds - the number of `ids`
5632fe279fdSBarry Smith . ids    - An array of the label ids in sequence for the domain
5642fe279fdSBarry Smith - height - Height of target cells in `DMPLEX` topology
5655f15299fSJeremy L Thompson 
5665f15299fSJeremy L Thompson   Output Parameters:
5672fe279fdSBarry Smith + point - the first labeled point
568a3b724e8SBarry Smith - ds    - the `PetscDS` corresponding to the first labeled point
5695f15299fSJeremy L Thompson 
5705f15299fSJeremy L Thompson   Level: developer
571a1cb98faSBarry Smith 
5721cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS`
573817da375SSatish Balay @*/
574d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
575d71ae5a4SJacob Faibussowitsch {
576a6e0b375SMatthew G. Knepley   DM              plex;
577a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
578e9f0ba4eSJed Brown   PetscInt        ls = -1;
579e5e52638SMatthew G. Knepley 
580e5e52638SMatthew G. Knepley   PetscFunctionBegin;
5815f15299fSJeremy L Thompson   if (point) *point = -1;
5823ba16761SJacob Faibussowitsch   if (!label) PetscFunctionReturn(PETSC_SUCCESS);
5839566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
5849566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
585e9f0ba4eSJed Brown   for (PetscInt i = 0; i < numIds; ++i) {
586e9f0ba4eSJed Brown     IS       labelIS;
587e9f0ba4eSJed Brown     PetscInt num_points, pStart, pEnd;
5889566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
589e9f0ba4eSJed Brown     if (!labelIS) continue; /* No points with that id on this process */
5909566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
5919566063dSJacob Faibussowitsch     PetscCall(ISGetSize(labelIS, &num_points));
592e9f0ba4eSJed Brown     if (num_points) {
593e5e52638SMatthew G. Knepley       const PetscInt *points;
5949566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(labelIS, &points));
595e9f0ba4eSJed Brown       for (PetscInt i = 0; i < num_points; i++) {
596e9f0ba4eSJed Brown         PetscInt point;
5979566063dSJacob Faibussowitsch         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
598e9f0ba4eSJed Brown         if (pStart <= point && point < pEnd) {
599a6e0b375SMatthew G. Knepley           ls = point;
60079f2dbaeSMatthew G. Knepley           if (ds) {
60179f2dbaeSMatthew G. Knepley             // If this is a face of a cohesive cell, then prefer that DS
60279f2dbaeSMatthew G. Knepley             if (height == 1) {
60379f2dbaeSMatthew G. Knepley               const PetscInt *supp;
60479f2dbaeSMatthew G. Knepley               PetscInt        suppSize;
60579f2dbaeSMatthew G. Knepley               DMPolytopeType  ct;
60679f2dbaeSMatthew G. Knepley 
60779f2dbaeSMatthew G. Knepley               PetscCall(DMPlexGetSupport(dm, ls, &supp));
60879f2dbaeSMatthew G. Knepley               PetscCall(DMPlexGetSupportSize(dm, ls, &suppSize));
60979f2dbaeSMatthew G. Knepley               for (PetscInt s = 0; s < suppSize; ++s) {
61079f2dbaeSMatthew G. Knepley                 PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
61179f2dbaeSMatthew 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)) {
61279f2dbaeSMatthew G. Knepley                   ls = supp[s];
61379f2dbaeSMatthew G. Knepley                   break;
61479f2dbaeSMatthew G. Knepley                 }
61579f2dbaeSMatthew G. Knepley               }
61679f2dbaeSMatthew G. Knepley             }
61779f2dbaeSMatthew G. Knepley             PetscCall(DMGetCellDS(dm, ls, ds, NULL));
61879f2dbaeSMatthew G. Knepley           }
6198e6d238bSMatthew G. Knepley           if (ls >= 0) break;
620e5e52638SMatthew G. Knepley         }
621e9f0ba4eSJed Brown       }
6229566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(labelIS, &points));
623e9f0ba4eSJed Brown     }
6249566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&labelIS));
625e5e52638SMatthew G. Knepley     if (ls >= 0) break;
626e5e52638SMatthew G. Knepley   }
6275f15299fSJeremy L Thompson   if (point) *point = ls;
6289566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
6293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
630e5e52638SMatthew G. Knepley }
631e5e52638SMatthew G. Knepley 
6320de51b56SMatthew G. Knepley /*
6330de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
6340de51b56SMatthew G. Knepley 
6350de51b56SMatthew G. Knepley   There are several different scenarios:
6360de51b56SMatthew G. Knepley 
6370de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
6380de51b56SMatthew G. Knepley 
6390de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
6400de51b56SMatthew G. Knepley 
6410de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
6420de51b56SMatthew G. Knepley 
6430de51b56SMatthew 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.
6440de51b56SMatthew G. Knepley 
6450de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
6460de51b56SMatthew G. Knepley 
6470de51b56SMatthew 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.
6480de51b56SMatthew G. Knepley 
649636b9322SMatthew G. Knepley   4) Volumetric input mesh with boundary output mesh
650636b9322SMatthew G. Knepley 
651636b9322SMatthew G. Knepley      Here we must get a subspace for the input DS
652636b9322SMatthew G. Knepley 
6530de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
6540de51b56SMatthew G. Knepley 
6550de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
6560de51b56SMatthew G. Knepley 
6570de51b56SMatthew 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.
6580de51b56SMatthew 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
6595f790a90SMatthew 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
6600de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
6610de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
6620de51b56SMatthew G. Knepley 
6630de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
6640de51b56SMatthew G. Knepley 
6650de51b56SMatthew 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.
6660de51b56SMatthew G. Knepley */
667d71ae5a4SJacob 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)
668d71ae5a4SJacob Faibussowitsch {
669a6e0b375SMatthew G. Knepley   DM               plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
670a6e0b375SMatthew G. Knepley   DMEnclosureType  encIn, encAux;
671a6e0b375SMatthew G. Knepley   PetscDS          ds = NULL, dsIn = NULL, dsAux = NULL;
672ca3d3a14SMatthew G. Knepley   Vec              localA = NULL, tv;
673aa0eca99SMatthew G. Knepley   IS               fieldIS;
67447923291SMatthew G. Knepley   PetscSection     section;
675636b9322SMatthew G. Knepley   PetscDualSpace  *sp, *cellsp, *spIn, *cellspIn;
676ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
6778c6c5593SMatthew G. Knepley   PetscInt        *Nc;
67879f2dbaeSMatthew G. Knepley   PetscInt         dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, minHeightIn, minHeightAux, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
6798e6d238bSMatthew G. Knepley   PetscBool       *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform;
680c330f8ffSToby Isaac   DMField          coordField;
681c330f8ffSToby Isaac   DMLabel          depthLabel;
682c330f8ffSToby Isaac   PetscQuadrature  allPoints = NULL;
68347923291SMatthew G. Knepley 
68447923291SMatthew G. Knepley   PetscFunctionBegin;
6859566063dSJacob Faibussowitsch   if (localU) PetscCall(VecGetDM(localU, &dmIn));
686ad540459SPierre Jolivet   else dmIn = dm;
6879566063dSJacob Faibussowitsch   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
6889371c9d4SSatish Balay   if (localA) PetscCall(VecGetDM(localA, &dmAux));
689ad540459SPierre Jolivet   else dmAux = NULL;
6909566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
6919566063dSJacob Faibussowitsch   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
6929566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
6939566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
6949566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
6959566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
6969566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
6979566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
6989566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
6990de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
700a6e0b375SMatthew G. Knepley   if (dmAux) {
7019566063dSJacob Faibussowitsch     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
702ad540459SPierre 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");
703636b9322SMatthew G. Knepley   }
7045ee220baSJed Brown   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plexIn, PETSC_TRUE, localU, time, NULL, NULL, NULL));
705e04ae0b4SMatthew G. Knepley   PetscCall(DMGetCoordinateField(dm, &coordField));
706e44f6aebSMatthew G. Knepley   PetscCheck(coordField, PETSC_COMM_SELF, PETSC_ERR_USER, "DM must have a coordinate field");
707e04ae0b4SMatthew G. Knepley   /**** No collective calls below this point ****/
708636b9322SMatthew G. Knepley   /* Determine height for iteration of all meshes */
709636b9322SMatthew G. Knepley   {
710636b9322SMatthew G. Knepley     DMPolytopeType ct, ctIn, ctAux;
71179f2dbaeSMatthew G. Knepley     PetscInt       lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux;
71219552267SMatthew G. Knepley     PetscInt       dim = -1, dimIn = -1, dimAux = -1;
71388aca1feSMatthew G. Knepley 
7149566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
715636b9322SMatthew G. Knepley     if (pEnd > pStart) {
7169566063dSJacob Faibussowitsch       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
717636b9322SMatthew G. Knepley       p = lStart < 0 ? pStart : lStart;
7189566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plex, p, &ct));
719636b9322SMatthew G. Knepley       dim = DMPolytopeTypeGetDim(ct);
7209566063dSJacob Faibussowitsch       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
7219566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
7229566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
723636b9322SMatthew G. Knepley       dimIn = DMPolytopeTypeGetDim(ctIn);
724636b9322SMatthew G. Knepley       if (dmAux) {
7259566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
72619552267SMatthew G. Knepley         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
72719552267SMatthew G. Knepley         if (pStartAux < pEndAux) {
7289566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
729636b9322SMatthew G. Knepley           dimAux = DMPolytopeTypeGetDim(ctAux);
73019552267SMatthew G. Knepley         }
731636b9322SMatthew G. Knepley       } else dimAux = dim;
732e04ae0b4SMatthew G. Knepley     } else {
733e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plex));
734e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plexIn));
735e04ae0b4SMatthew G. Knepley       if (dmAux) PetscCall(DMDestroy(&plexAux));
7363ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
737e39fcb4eSStefano Zampini     }
738636b9322SMatthew G. Knepley     if (dim < 0) {
739636b9322SMatthew G. Knepley       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
740636b9322SMatthew G. Knepley 
741636b9322SMatthew G. Knepley       /* Fall back to determination based on being a submesh */
7429566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plex, &spmap));
7439566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
7449566063dSJacob Faibussowitsch       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
745636b9322SMatthew G. Knepley       dim    = spmap ? 1 : 0;
746636b9322SMatthew G. Knepley       dimIn  = spmapIn ? 1 : 0;
747636b9322SMatthew G. Knepley       dimAux = spmapAux ? 1 : 0;
74888aca1feSMatthew G. Knepley     }
749636b9322SMatthew G. Knepley     {
75057508eceSPierre Jolivet       PetscInt dimProj   = PetscMin(PetscMin(dim, dimIn), dimAux < 0 ? PETSC_INT_MAX : dimAux);
75119552267SMatthew G. Knepley       PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;
752636b9322SMatthew G. Knepley 
75319552267SMatthew 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");
754636b9322SMatthew G. Knepley       if (dimProj < dim) minHeight = 1;
755636b9322SMatthew G. Knepley       htInc    = dim - dimProj;
756636b9322SMatthew G. Knepley       htIncIn  = dimIn - dimProj;
75719552267SMatthew G. Knepley       htIncAux = dimAuxEff - dimProj;
7580de51b56SMatthew G. Knepley     }
759a6e0b375SMatthew G. Knepley   }
7609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
7619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
7629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
7632af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
76463a3b9bcSJacob 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);
76579f2dbaeSMatthew G. Knepley   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, NULL, &ds));
7669566063dSJacob Faibussowitsch   if (!ds) PetscCall(DMGetDS(dm, &ds));
76779f2dbaeSMatthew G. Knepley   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, minHeight, NULL, &dsIn));
7689566063dSJacob Faibussowitsch   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
7699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
7709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
7718e6d238bSMatthew G. Knepley   PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
7728e6d238bSMatthew G. Knepley   if (isCohesiveIn) --htIncIn; // Should be rearranged
7739566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &NfTot));
7749566063dSJacob Faibussowitsch   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
77507218a29SMatthew G. Knepley   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL));
7769566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
7779566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
7789566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
7790c898c18SMatthew G. Knepley   if (dmAux) {
7809566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmAux, &dsAux));
7819566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
7820c898c18SMatthew G. Knepley   }
7839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
7849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
7859566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
7869371c9d4SSatish Balay   else {
7879371c9d4SSatish Balay     cellsp   = sp;
7889371c9d4SSatish Balay     cellspIn = spIn;
7899371c9d4SSatish Balay   }
7908c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
79147923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
792665f567fSMatthew G. Knepley     PetscDiscType disctype;
79347923291SMatthew G. Knepley 
7949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
795665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
796665f567fSMatthew G. Knepley       PetscFE fe;
79747923291SMatthew G. Knepley 
79847923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
799665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
8009566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
8019566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
802665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
803665f567fSMatthew G. Knepley       PetscFV fv;
80447923291SMatthew G. Knepley 
80547923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
806665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
8079566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv));
8089566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
809665f567fSMatthew G. Knepley     } else {
810665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
811665f567fSMatthew G. Knepley       cellsp[f] = NULL;
812665f567fSMatthew G. Knepley     }
81347923291SMatthew G. Knepley   }
814636b9322SMatthew G. Knepley   for (f = 0; f < NfIn; ++f) {
815636b9322SMatthew G. Knepley     PetscDiscType disctype;
816636b9322SMatthew G. Knepley 
8179566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
818636b9322SMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
819636b9322SMatthew G. Knepley       PetscFE fe;
820636b9322SMatthew G. Knepley 
8219566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
8229566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
823636b9322SMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
824636b9322SMatthew G. Knepley       PetscFV fv;
825636b9322SMatthew G. Knepley 
8269566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv));
8279566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
828636b9322SMatthew G. Knepley     } else {
829636b9322SMatthew G. Knepley       cellspIn[f] = NULL;
830636b9322SMatthew G. Knepley     }
831636b9322SMatthew G. Knepley   }
832636b9322SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
8339371c9d4SSatish Balay     if (!htInc) {
8349371c9d4SSatish Balay       sp[f] = cellsp[f];
8359371c9d4SSatish Balay     } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
836636b9322SMatthew G. Knepley   }
837ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
83874fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
839665f567fSMatthew G. Knepley     PetscDiscType    disctype;
8404a3e9fdbSToby Isaac     const PetscReal *points;
841409bb69aSDarsh Nathawani     PetscInt         numPoints, k;
8420c898c18SMatthew G. Knepley 
84308401ef6SPierre Jolivet     PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
8449566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints));
8459566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
8469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
847a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
8489371c9d4SSatish Balay       if (!htIncIn) {
8499371c9d4SSatish Balay         spIn[f] = cellspIn[f];
8509371c9d4SSatish Balay       } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
851636b9322SMatthew G. Knepley 
8529566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
853665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
8549566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem));
8559371c9d4SSatish Balay       if (!htIncIn) {
8569371c9d4SSatish Balay         subfem = fem;
8579371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
858409bb69aSDarsh Nathawani       PetscCall(PetscDSGetJetDegree(dsIn, f, &k));
859409bb69aSDarsh Nathawani       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &T[f]));
8600c898c18SMatthew G. Knepley     }
8610c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
8629566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
863665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
8649566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem));
8659371c9d4SSatish Balay       if (!htIncAux) {
8669371c9d4SSatish Balay         subfem = fem;
8679371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
868409bb69aSDarsh Nathawani       PetscCall(PetscDSGetJetDegree(dsAux, f, &k));
869409bb69aSDarsh Nathawani       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &TAux[f]));
8700c898c18SMatthew G. Knepley     }
8710c898c18SMatthew G. Knepley   }
87247923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
8732af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
874636b9322SMatthew G. Knepley     PetscInt     hEff     = h - minHeight + htInc;
875636b9322SMatthew G. Knepley     PetscInt     hEffIn   = h - minHeight + htIncIn;
876636b9322SMatthew G. Knepley     PetscInt     hEffAux  = h - minHeight + htIncAux;
877a6e0b375SMatthew G. Knepley     PetscDS      dsEff    = ds;
878636b9322SMatthew G. Knepley     PetscDS      dsEffIn  = dsIn;
879636b9322SMatthew G. Knepley     PetscDS      dsEffAux = dsAux;
8808c6c5593SMatthew G. Knepley     PetscScalar *values;
881b7260050SToby Isaac     PetscBool   *fieldActive;
882b7260050SToby Isaac     PetscInt     maxDegree;
883e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
884c330f8ffSToby Isaac     IS           heightIS;
8858c6c5593SMatthew G. Knepley 
886636b9322SMatthew G. Knepley     if (h > minHeight) {
8879566063dSJacob Faibussowitsch       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
888636b9322SMatthew G. Knepley     }
8899566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
8909566063dSJacob Faibussowitsch     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
8919566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
892c330f8ffSToby Isaac     if (pEnd <= pStart) {
8939566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
894c330f8ffSToby Isaac       continue;
895c330f8ffSToby Isaac     }
89647923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
89747923291SMatthew G. Knepley     totDim = 0;
89847923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
8995fedec97SMatthew G. Knepley       PetscBool cohesive;
9005fedec97SMatthew G. Knepley 
901665f567fSMatthew G. Knepley       if (!sp[f]) continue;
9029566063dSJacob Faibussowitsch       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
9039566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
9049c3cf19fSMatthew G. Knepley       totDim += spDim;
9055fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) totDim += spDim;
90647923291SMatthew G. Knepley     }
907636b9322SMatthew G. Knepley     p = lStart < 0 ? pStart : lStart;
9089566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
90963a3b9bcSJacob 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);
910c330f8ffSToby Isaac     if (!totDim) {
9119566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
912c330f8ffSToby Isaac       continue;
913c330f8ffSToby Isaac     }
9149566063dSJacob Faibussowitsch     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
915636b9322SMatthew G. Knepley     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
916636b9322SMatthew G. Knepley     if (localU) {
917636b9322SMatthew G. Knepley       PetscInt totDimIn, pIn, numValuesIn;
918636b9322SMatthew G. Knepley 
919636b9322SMatthew G. Knepley       totDimIn = 0;
920636b9322SMatthew G. Knepley       for (f = 0; f < NfIn; ++f) {
9215fedec97SMatthew G. Knepley         PetscBool cohesive;
9225fedec97SMatthew G. Knepley 
923636b9322SMatthew G. Knepley         if (!spIn[f]) continue;
9249566063dSJacob Faibussowitsch         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
9259566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
926636b9322SMatthew G. Knepley         totDimIn += spDim;
9278e6d238bSMatthew G. Knepley         if (isCohesiveIn && !cohesive) totDimIn += spDim;
928636b9322SMatthew G. Knepley       }
9299566063dSJacob Faibussowitsch       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
9309566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
9318e6d238bSMatthew G. Knepley       // TODO We could check that pIn is a cohesive cell for this check
9328e6d238bSMatthew 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);
9339566063dSJacob Faibussowitsch       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
934636b9322SMatthew G. Knepley     }
9359566063dSJacob Faibussowitsch     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
93647923291SMatthew G. Knepley     /* Loop over points at this height */
9379566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
9389566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
939aa0eca99SMatthew G. Knepley     {
940aa0eca99SMatthew G. Knepley       const PetscInt *fields;
941aa0eca99SMatthew G. Knepley 
9429566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
943ad540459SPierre Jolivet       for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE;
944ad540459SPierre Jolivet       for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
9459566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
946aa0eca99SMatthew G. Knepley     }
9478c6c5593SMatthew G. Knepley     if (label) {
9488c6c5593SMatthew G. Knepley       PetscInt i;
94947923291SMatthew G. Knepley 
95047923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
951c330f8ffSToby Isaac         IS              pointIS, isectIS;
95247923291SMatthew G. Knepley         const PetscInt *points;
9538c6c5593SMatthew G. Knepley         PetscInt        n;
954c330f8ffSToby Isaac         PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
955c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
95647923291SMatthew G. Knepley 
9579566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
95847923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
9599566063dSJacob Faibussowitsch         PetscCall(ISIntersect(pointIS, heightIS, &isectIS));
9609566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&pointIS));
961c330f8ffSToby Isaac         if (!isectIS) continue;
9629566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(isectIS, &n));
9639566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(isectIS, &points));
9649566063dSJacob Faibussowitsch         PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree));
96548a46eb9SPierre Jolivet         if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad));
966c330f8ffSToby Isaac         if (!quad) {
967c330f8ffSToby Isaac           if (!h && allPoints) {
968c330f8ffSToby Isaac             quad      = allPoints;
969c330f8ffSToby Isaac             allPoints = NULL;
970c330f8ffSToby Isaac           } else {
9719566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad));
972c330f8ffSToby Isaac           }
973c330f8ffSToby Isaac         }
974*ac9d17c7SMatthew G. Knepley         PetscFEGeomMode geommode = htInc && h == minHeight ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC;
97579f2dbaeSMatthew G. Knepley 
97679f2dbaeSMatthew G. Knepley         if (n) {
97779f2dbaeSMatthew G. Knepley           PetscInt depth, dep;
97879f2dbaeSMatthew G. Knepley 
97979f2dbaeSMatthew G. Knepley           PetscCall(DMPlexGetDepth(dm, &depth));
98079f2dbaeSMatthew G. Knepley           PetscCall(DMPlexGetPointDepth(dm, points[0], &dep));
981*ac9d17c7SMatthew G. Knepley           if (dep < depth && h == minHeight) geommode = PETSC_FEGEOM_BOUNDARY;
98279f2dbaeSMatthew G. Knepley         }
983*ac9d17c7SMatthew G. Knepley         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, geommode, &fegeom));
98447923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
98547923291SMatthew G. Knepley           const PetscInt point = points[p];
98647923291SMatthew G. Knepley 
9879566063dSJacob Faibussowitsch           PetscCall(PetscArrayzero(values, numValues));
9889566063dSJacob Faibussowitsch           PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom));
9899566063dSJacob Faibussowitsch           PetscCall(DMPlexSetActivePoint(dm, point));
990d0609cedSBarry 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));
9919566063dSJacob Faibussowitsch           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
9929566063dSJacob Faibussowitsch           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
99347923291SMatthew G. Knepley         }
9949566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom));
9959566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&fegeom));
9969566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureDestroy(&quad));
9979566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(isectIS, &points));
9989566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isectIS));
99947923291SMatthew G. Knepley       }
10008c6c5593SMatthew G. Knepley     } else {
1001c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
1002c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
1003c330f8ffSToby Isaac       IS              pointIS;
1004c330f8ffSToby Isaac 
10059566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS));
10069566063dSJacob Faibussowitsch       PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
100748a46eb9SPierre Jolivet       if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
1008c330f8ffSToby Isaac       if (!quad) {
1009c330f8ffSToby Isaac         if (!h && allPoints) {
1010c330f8ffSToby Isaac           quad      = allPoints;
1011c330f8ffSToby Isaac           allPoints = NULL;
1012c330f8ffSToby Isaac         } else {
10139566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad));
1014c330f8ffSToby Isaac         }
1015c330f8ffSToby Isaac       }
1016*ac9d17c7SMatthew G. Knepley       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC, &fegeom));
10178c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
10189566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(values, numValues));
10199566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom));
10209566063dSJacob Faibussowitsch         PetscCall(DMPlexSetActivePoint(dm, p));
1021d0609cedSBarry 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));
10229566063dSJacob Faibussowitsch         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
10239566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
10248c6c5593SMatthew G. Knepley       }
10259566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom));
10269566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomDestroy(&fegeom));
10279566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
10289566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
10298c6c5593SMatthew G. Knepley     }
10309566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&heightIS));
10319566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
10329566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
103347923291SMatthew G. Knepley   }
10348c6c5593SMatthew G. Knepley   /* Cleanup */
1035ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
10369566063dSJacob Faibussowitsch     for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f]));
10379566063dSJacob Faibussowitsch     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
10389566063dSJacob Faibussowitsch     PetscCall(PetscFree2(T, TAux));
10390c898c18SMatthew G. Knepley   }
10409566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&allPoints));
10419566063dSJacob Faibussowitsch   PetscCall(PetscFree3(isFE, sp, spIn));
10429566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
10439566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
10449566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plexIn));
10459566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMDestroy(&plexAux));
10463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104747923291SMatthew G. Knepley }
10488c6c5593SMatthew G. Knepley 
1049d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
1050d71ae5a4SJacob Faibussowitsch {
10518c6c5593SMatthew G. Knepley   PetscFunctionBegin;
10529566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
10533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10548c6c5593SMatthew G. Knepley }
10558c6c5593SMatthew G. Knepley 
1056d71ae5a4SJacob 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)
1057d71ae5a4SJacob Faibussowitsch {
10588c6c5593SMatthew G. Knepley   PetscFunctionBegin;
10599566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
10603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106147923291SMatthew G. Knepley }
106247923291SMatthew G. Knepley 
1063d71ae5a4SJacob 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)
1064d71ae5a4SJacob Faibussowitsch {
106547923291SMatthew G. Knepley   PetscFunctionBegin;
10669566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
10673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106847923291SMatthew G. Knepley }
106947923291SMatthew G. Knepley 
1070d71ae5a4SJacob 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)
1071d71ae5a4SJacob Faibussowitsch {
10728c6c5593SMatthew G. Knepley   PetscFunctionBegin;
10739566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
10743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
107547923291SMatthew G. Knepley }
1076ece3a9fcSMatthew G. Knepley 
1077d71ae5a4SJacob 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)
1078d71ae5a4SJacob Faibussowitsch {
1079ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
10809566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX));
10813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1082ece3a9fcSMatthew G. Knepley }
1083