xref: /petsc/src/dm/impls/plex/plexproject.c (revision e04ae0b46864b7ce274fa2f3362ed93b23f5a1f7)
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 
81b32699bSMatthew G. Knepley   Not collective
91b32699bSMatthew G. Knepley 
104165533cSJose E. Roman   Input Parameter:
111b32699bSMatthew G. Knepley . 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 
18db781477SPatrick Sanan .seealso: `DMPlexSetActivePoint()`
191b32699bSMatthew G. Knepley @*/
20bdb10af2SPierre Jolivet PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point)
21bdb10af2SPierre Jolivet {
221b32699bSMatthew G. Knepley   PetscFunctionBeginHot;
231b32699bSMatthew G. Knepley   *point = ((DM_Plex *) dm->data)->activePoint;
241b32699bSMatthew G. Knepley   PetscFunctionReturn(0);
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 
301b32699bSMatthew G. Knepley   Not collective
311b32699bSMatthew G. Knepley 
324165533cSJose E. Roman   Input Parameters:
331b32699bSMatthew G. Knepley + 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 
38db781477SPatrick Sanan .seealso: `DMPlexGetActivePoint()`
391b32699bSMatthew G. Knepley @*/
40bdb10af2SPierre Jolivet PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point)
41bdb10af2SPierre Jolivet {
421b32699bSMatthew G. Knepley   PetscFunctionBeginHot;
431b32699bSMatthew G. Knepley   ((DM_Plex *) dm->data)->activePoint = point;
441b32699bSMatthew G. Knepley   PetscFunctionReturn(0);
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:
51a6e0b375SMatthew G. Knepley + dm     - The output DM
52a6e0b375SMatthew G. Knepley . ds     - The output DS
53a6e0b375SMatthew G. Knepley . dmIn   - The input DM
54a6e0b375SMatthew G. Knepley . 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
59a6e0b375SMatthew G. Knepley . 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 
68db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()`
69a6e0b375SMatthew G. Knepley */
70a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
718c6c5593SMatthew G. Knepley                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
728c6c5593SMatthew G. Knepley                                                   PetscScalar values[])
7347923291SMatthew G. Knepley {
74a6e0b375SMatthew G. Knepley   PetscInt       coordDim, Nf, *Nc, f, spDim, d, v, tp;
755fedec97SMatthew G. Knepley   PetscBool      isAffine, isCohesive, transform;
768c6c5593SMatthew G. Knepley 
778c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
789566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmIn, &coordDim));
799566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
829566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
838c6c5593SMatthew G. Knepley   /* Get values for closure */
84c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
85c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
868c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
875fedec97SMatthew G. Knepley     PetscBool    cohesive;
888c6c5593SMatthew G. Knepley 
898c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
909566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
919566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
928c6c5593SMatthew G. Knepley     if (funcs[f]) {
93c330f8ffSToby Isaac       if (isFE[f]) {
94c330f8ffSToby Isaac         PetscQuadrature   allPoints;
95c330f8ffSToby Isaac         PetscInt          q, dim, numPoints;
96c330f8ffSToby Isaac         const PetscReal   *points;
97c330f8ffSToby Isaac         PetscScalar       *pointEval;
98c330f8ffSToby Isaac         PetscReal         *x;
99ca3d3a14SMatthew G. Knepley         DM                rdm;
100c330f8ffSToby Isaac 
1019566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDM(sp[f],&rdm));
1029566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
1039566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL));
1049566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
1059566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x));
106737e23dcSMatthew G. Knepley         PetscCall(PetscArrayzero(pointEval, numPoints*Nc[f]));
107c330f8ffSToby Isaac         for (q = 0; q < numPoints; q++, tp++) {
108c330f8ffSToby Isaac           const PetscReal *v0;
109c330f8ffSToby Isaac 
110c330f8ffSToby Isaac           if (isAffine) {
111665f567fSMatthew G. Knepley             const PetscReal *refpoint = &points[q*dim];
112665f567fSMatthew G. Knepley             PetscReal        injpoint[3] = {0., 0., 0.};
113665f567fSMatthew G. Knepley 
114665f567fSMatthew G. Knepley             if (dim != fegeom->dim) {
1155fedec97SMatthew G. Knepley               if (isCohesive) {
116665f567fSMatthew G. Knepley                 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
117665f567fSMatthew G. Knepley                 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
118665f567fSMatthew G. Knepley                 refpoint = injpoint;
11963a3b9bcSJacob Faibussowitsch               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim);
120665f567fSMatthew G. Knepley             }
121665f567fSMatthew G. Knepley             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x);
122c330f8ffSToby Isaac             v0 = x;
1238c6c5593SMatthew G. Knepley           } else {
124c330f8ffSToby Isaac             v0 = &fegeom->v[tp*coordDim];
1258c6c5593SMatthew G. Knepley           }
1269566063dSJacob Faibussowitsch           if (transform) {PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx)); v0 = x;}
1279566063dSJacob Faibussowitsch           PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx));
128c330f8ffSToby Isaac         }
1294bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
1309566063dSJacob Faibussowitsch         PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
1319566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
1329566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x));
1339566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
134c330f8ffSToby Isaac         v += spDim;
1355fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) {
13627f02ce8SMatthew G. Knepley           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
13727f02ce8SMatthew G. Knepley         }
138c330f8ffSToby Isaac       } else {
139c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
1409566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
141c330f8ffSToby Isaac         }
142c330f8ffSToby Isaac       }
143c330f8ffSToby Isaac     } else {
14427f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
1455fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
14627f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
14727f02ce8SMatthew G. Knepley       }
1488c6c5593SMatthew G. Knepley     }
1499c3cf19fSMatthew G. Knepley   }
1508c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1518c6c5593SMatthew G. Knepley }
1528c6c5593SMatthew G. Knepley 
153a6e0b375SMatthew G. Knepley /*
154a6e0b375SMatthew G. Knepley   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
155a6e0b375SMatthew G. Knepley 
156a6e0b375SMatthew G. Knepley   Input Parameters:
157a6e0b375SMatthew G. Knepley + dm             - The output DM
158a6e0b375SMatthew G. Knepley . ds             - The output DS
159a6e0b375SMatthew G. Knepley . dmIn           - The input DM
160a6e0b375SMatthew G. Knepley . dsIn           - The input DS
161a6e0b375SMatthew G. Knepley . dmAux          - The auxiliary DM, which is always for the input space
162a6e0b375SMatthew G. Knepley . dsAux          - The auxiliary DS, which is always for the input space
163a6e0b375SMatthew G. Knepley . time           - The time for this evaluation
164a6e0b375SMatthew G. Knepley . localU         - The local solution
165a6e0b375SMatthew G. Knepley . localA         - The local auziliary fields
166a6e0b375SMatthew G. Knepley . cgeom          - The FE geometry for this point
167a6e0b375SMatthew G. Knepley . sp             - The output PetscDualSpace for each field
168a6e0b375SMatthew G. Knepley . p              - The point in the output DM
169a5b23f4aSJose E. Roman . T              - Input basis and derivatives for each field tabulated on the quadrature points
170ef0bb6c7SMatthew G. Knepley . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
171a6e0b375SMatthew G. Knepley . funcs          - The evaluation function for each field
172a6e0b375SMatthew G. Knepley - ctxs           - The user context for each field
173a6e0b375SMatthew G. Knepley 
174a6e0b375SMatthew G. Knepley   Output Parameter:
175a6e0b375SMatthew G. Knepley . values         - The value for each dual basis vector in the output dual space
176a6e0b375SMatthew G. Knepley 
177a6e0b375SMatthew G. Knepley   Note: Not supported for FV
178a6e0b375SMatthew G. Knepley 
179a6e0b375SMatthew G. Knepley   Level: developer
180a6e0b375SMatthew G. Knepley 
181db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()`
182a6e0b375SMatthew G. Knepley */
183a6e0b375SMatthew G. Knepley 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,
184ef0bb6c7SMatthew G. Knepley                                                    PetscTabulation *T, PetscTabulation *TAux,
1858c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
1868c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1878c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
188191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
1898c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
1908c6c5593SMatthew G. Knepley {
1918c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
1924bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
1938c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
1948c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
195191494d9SMatthew G. Knepley   const PetscScalar *constants;
1968c6c5593SMatthew G. Knepley   PetscReal         *x;
197ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
1984bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
1994bee2e38SMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed;
200ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
2015fedec97SMatthew G. Knepley   PetscBool          isAffine, isCohesive, transform;
2028c6c5593SMatthew G. Knepley 
2038c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
2049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
2069566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
2079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
2089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
2099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
2109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
2119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
2129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
2139566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
2149566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmIn, &section));
2159566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
2169566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
2178c6c5593SMatthew G. Knepley   if (dmAux) {
21844171101SMatthew G. Knepley     PetscInt subp;
2191ba34526SMatthew G. Knepley 
2209566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
2219566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2229566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
2239566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2249566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2259566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
2269566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
2278c6c5593SMatthew G. Knepley   }
2288c6c5593SMatthew G. Knepley   /* Get values for closure */
2294bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
23027f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
23127f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
2324bee2e38SMatthew G. Knepley   if (isAffine) {
2334bee2e38SMatthew G. Knepley     fegeom.v    = x;
2344bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
2354bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
2364bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
2374bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
2384bee2e38SMatthew G. Knepley   }
2398c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
240c330f8ffSToby Isaac     PetscQuadrature  allPoints;
241c330f8ffSToby Isaac     PetscInt         q, dim, numPoints;
242c330f8ffSToby Isaac     const PetscReal *points;
243c330f8ffSToby Isaac     PetscScalar     *pointEval;
2445fedec97SMatthew G. Knepley     PetscBool        cohesive;
245c330f8ffSToby Isaac     DM               dm;
246c330f8ffSToby Isaac 
2478c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2489566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
2499566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
250c330f8ffSToby Isaac     if (!funcs[f]) {
251be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
2525fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
25327f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
25427f02ce8SMatthew G. Knepley       }
255c330f8ffSToby Isaac       continue;
256c330f8ffSToby Isaac     }
2579566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f],&dm));
2589566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
2599566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL));
2609566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
2610c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
262c330f8ffSToby Isaac       if (isAffine) {
263ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
2641c531cf8SMatthew G. Knepley       } else {
2654bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp*dE];
2664bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp*dE*dE];
2674bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
2684bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
2698c6c5593SMatthew G. Knepley       }
2709566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
2719566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
2729566063dSJacob Faibussowitsch       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
273a6e0b375SMatthew 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]);
2741c531cf8SMatthew G. Knepley     }
2759566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
2769566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
277c330f8ffSToby Isaac     v += spDim;
27827f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
2795fedec97SMatthew G. Knepley     if (isCohesive && !cohesive) {
28027f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
28127f02ce8SMatthew G. Knepley     }
2828c6c5593SMatthew G. Knepley   }
2839566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
2849566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
2858c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2868c6c5593SMatthew G. Knepley }
2878c6c5593SMatthew G. Knepley 
288a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p,
289ef0bb6c7SMatthew G. Knepley                                                      PetscTabulation *T, PetscTabulation *TAux,
290b18799e0SMatthew G. Knepley                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
291b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
292b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
293b18799e0SMatthew G. Knepley                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
294b18799e0SMatthew G. Knepley                                                      PetscScalar values[])
295b18799e0SMatthew G. Knepley {
296b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
2974bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
298b18799e0SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
299b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
300b18799e0SMatthew G. Knepley   const PetscScalar *constants;
301b18799e0SMatthew G. Knepley   PetscReal         *x;
302ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
3039f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
3044bee2e38SMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed;
305ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
306b18799e0SMatthew G. Knepley   PetscBool          isAffine;
307b18799e0SMatthew G. Knepley 
308b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
30908401ef6SPierre Jolivet   PetscCheck(dm == dmIn,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
3109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
3129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
3139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
3149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x));
3159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
3169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
3179566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
3189566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients));
319b18799e0SMatthew G. Knepley   if (dmAux) {
320a6e0b375SMatthew G. Knepley     PetscInt subp;
321b18799e0SMatthew G. Knepley 
3229566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
3239566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3249566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
3259566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3269566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3279566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
3289566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
329b18799e0SMatthew G. Knepley   }
330b18799e0SMatthew G. Knepley   /* Get values for closure */
3314bee2e38SMatthew G. Knepley   isAffine = fgeom->isAffine;
332ea78f98cSLisandro Dalcin   fegeom.n  = NULL;
333ea78f98cSLisandro Dalcin   fegeom.J  = NULL;
334ea78f98cSLisandro Dalcin   fegeom.v  = NULL;
335ea78f98cSLisandro Dalcin   fegeom.xi = NULL;
336a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
337a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
3384bee2e38SMatthew G. Knepley   if (isAffine) {
3394bee2e38SMatthew G. Knepley     fegeom.v    = x;
3404bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
3414bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
3424bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
3434bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
3444bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
3459f209ee4SMatthew G. Knepley 
3469f209ee4SMatthew G. Knepley     cgeom.J     = fgeom->suppJ[0];
3479f209ee4SMatthew G. Knepley     cgeom.invJ  = fgeom->suppInvJ[0];
3489f209ee4SMatthew G. Knepley     cgeom.detJ  = fgeom->suppDetJ[0];
3494bee2e38SMatthew G. Knepley   }
350b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
351b18799e0SMatthew G. Knepley     PetscQuadrature   allPoints;
352b18799e0SMatthew G. Knepley     PetscInt          q, dim, numPoints;
353b18799e0SMatthew G. Knepley     const PetscReal   *points;
354b18799e0SMatthew G. Knepley     PetscScalar       *pointEval;
355b18799e0SMatthew G. Knepley     DM                dm;
356b18799e0SMatthew G. Knepley 
357b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
3589566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
359b18799e0SMatthew G. Knepley     if (!funcs[f]) {
360b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
361b18799e0SMatthew G. Knepley       continue;
362b18799e0SMatthew G. Knepley     }
3639566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f],&dm));
3649566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
3659566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL));
3669566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
367b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
368b18799e0SMatthew G. Knepley       if (isAffine) {
369ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
370b18799e0SMatthew G. Knepley       } else {
3713fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp*dE];
3729f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp*dE*dE];
3739f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp*dE*dE];
3744bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
3754bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp*dE];
3769f209ee4SMatthew G. Knepley 
3779f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][tp*dE*dE];
3789f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][tp*dE*dE];
3799f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][tp];
380b18799e0SMatthew G. Knepley       }
381a6e0b375SMatthew 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 */
3829566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
3839566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
3844bee2e38SMatthew G. Knepley       (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, fegeom.n, numConstants, constants, &pointEval[Nc[f]*q]);
385b18799e0SMatthew G. Knepley     }
3869566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
3879566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
388b18799e0SMatthew G. Knepley     v += spDim;
389b18799e0SMatthew G. Knepley   }
3909566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients));
3919566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
392b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
393b18799e0SMatthew G. Knepley }
394b18799e0SMatthew G. Knepley 
395a6e0b375SMatthew G. Knepley 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[],
396ef0bb6c7SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux,
3978c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
3988c6c5593SMatthew G. Knepley {
3998c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
4008c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
4018c6c5593SMatthew G. Knepley 
4028c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
4039566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4049566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
4059566063dSJacob Faibussowitsch   if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
4068c6c5593SMatthew G. Knepley   switch (type) {
4078c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
4088c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
4099566063dSJacob 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));break;
4108c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
4118c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
412d0609cedSBarry Smith     PetscCall(DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
4130c898c18SMatthew G. Knepley                                            (void (**)(PetscInt, PetscInt, PetscInt,
4140c898c18SMatthew G. Knepley                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
4150c898c18SMatthew G. Knepley                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
416d0609cedSBarry Smith                                                       PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values));
417d0609cedSBarry Smith     break;
418b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
419d0609cedSBarry Smith     PetscCall(DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
420b18799e0SMatthew G. Knepley                                              (void (**)(PetscInt, PetscInt, PetscInt,
421b18799e0SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
422b18799e0SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
423d0609cedSBarry Smith                                                         PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values));
424d0609cedSBarry Smith     break;
42598921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
4268c6c5593SMatthew G. Knepley   }
4278c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
4288c6c5593SMatthew G. Knepley }
4298c6c5593SMatthew G. Knepley 
430df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
431df5c1128SToby Isaac {
432df5c1128SToby Isaac   PetscReal      *points;
433df5c1128SToby Isaac   PetscInt       f, numPoints;
434df5c1128SToby Isaac 
435df5c1128SToby Isaac   PetscFunctionBegin;
43619552267SMatthew G. Knepley   if (!dim) {
43719552267SMatthew G. Knepley     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
43819552267SMatthew G. Knepley     PetscFunctionReturn(0);
43919552267SMatthew G. Knepley   }
440df5c1128SToby Isaac   numPoints = 0;
441df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
442df5c1128SToby Isaac     if (funcs[f]) {
443df5c1128SToby Isaac       PetscQuadrature fAllPoints;
444df5c1128SToby Isaac       PetscInt        fNumPoints;
445df5c1128SToby Isaac 
4469566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL));
4479566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
448df5c1128SToby Isaac       numPoints += fNumPoints;
449df5c1128SToby Isaac     }
450df5c1128SToby Isaac   }
4519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim*numPoints,&points));
452df5c1128SToby Isaac   numPoints = 0;
453df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
454df5c1128SToby Isaac     if (funcs[f]) {
455df5c1128SToby Isaac       PetscQuadrature fAllPoints;
45654ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
457df5c1128SToby Isaac       const PetscReal *fPoints;
458df5c1128SToby Isaac 
4599566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL));
4609566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
46163a3b9bcSJacob 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);
462df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
463df5c1128SToby Isaac       numPoints += fNumPoints;
464df5c1128SToby Isaac     }
465df5c1128SToby Isaac   }
4669566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,allPoints));
4679566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL));
468df5c1128SToby Isaac   PetscFunctionReturn(0);
469df5c1128SToby Isaac }
470df5c1128SToby Isaac 
4715f15299fSJeremy L Thompson /*@C
4725f15299fSJeremy L Thompson   DMGetFirstLabeledPoint - Find first labeled point p_o in odm such that the corresponding point p in dm has the specified height. Return p and the corresponding ds.
4735f15299fSJeremy L Thompson 
4745f15299fSJeremy L Thompson   Input Parameters:
4755f15299fSJeremy L Thompson   dm - the DM
4765f15299fSJeremy L Thompson   odm - the enclosing DM
4775f15299fSJeremy L Thompson   label - label for DM domain, or NULL for whole domain
4785f15299fSJeremy L Thompson   numIds - the number of ids
4795f15299fSJeremy L Thompson   ids - An array of the label ids in sequence for the domain
4805f15299fSJeremy L Thompson   height - Height of target cells in DMPlex topology
4815f15299fSJeremy L Thompson 
4825f15299fSJeremy L Thompson   Output Parameters:
4835f15299fSJeremy L Thompson   point - the first labeled point
4845f15299fSJeremy L Thompson   ds - the ds corresponding to the first labeled point
4855f15299fSJeremy L Thompson 
4865f15299fSJeremy L Thompson   Level: developer
487817da375SSatish Balay @*/
4885f15299fSJeremy L Thompson PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
489e5e52638SMatthew G. Knepley {
490a6e0b375SMatthew G. Knepley   DM              plex;
491a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
492e9f0ba4eSJed Brown   PetscInt        ls = -1;
493e5e52638SMatthew G. Knepley 
494e5e52638SMatthew G. Knepley   PetscFunctionBegin;
4955f15299fSJeremy L Thompson   if (point) *point = -1;
496e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
4979566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
4989566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
499e9f0ba4eSJed Brown   for (PetscInt i = 0; i < numIds; ++i) {
500e9f0ba4eSJed Brown     IS       labelIS;
501e9f0ba4eSJed Brown     PetscInt num_points, pStart, pEnd;
5029566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
503e9f0ba4eSJed Brown     if (!labelIS) continue; /* No points with that id on this process */
5049566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
5059566063dSJacob Faibussowitsch     PetscCall(ISGetSize(labelIS, &num_points));
506e9f0ba4eSJed Brown     if (num_points) {
507e5e52638SMatthew G. Knepley       const PetscInt *points;
5089566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(labelIS, &points));
509e9f0ba4eSJed Brown       for (PetscInt i=0; i<num_points; i++) {
510e9f0ba4eSJed Brown         PetscInt point;
5119566063dSJacob Faibussowitsch         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
512e9f0ba4eSJed Brown         if (pStart <= point && point < pEnd) {
513a6e0b375SMatthew G. Knepley           ls = point;
5149566063dSJacob Faibussowitsch           if (ds) PetscCall(DMGetCellDS(dm, ls, ds));
515e5e52638SMatthew G. Knepley         }
516e9f0ba4eSJed Brown       }
5179566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(labelIS, &points));
518e9f0ba4eSJed Brown     }
5199566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&labelIS));
520e5e52638SMatthew G. Knepley     if (ls >= 0) break;
521e5e52638SMatthew G. Knepley   }
5225f15299fSJeremy L Thompson   if (point) *point = ls;
5239566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
524e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
525e5e52638SMatthew G. Knepley }
526e5e52638SMatthew G. Knepley 
5270de51b56SMatthew G. Knepley /*
5280de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
5290de51b56SMatthew G. Knepley 
5300de51b56SMatthew G. Knepley   There are several different scenarios:
5310de51b56SMatthew G. Knepley 
5320de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
5330de51b56SMatthew G. Knepley 
5340de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
5350de51b56SMatthew G. Knepley 
5360de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
5370de51b56SMatthew G. Knepley 
5380de51b56SMatthew 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.
5390de51b56SMatthew G. Knepley 
5400de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
5410de51b56SMatthew G. Knepley 
5420de51b56SMatthew 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.
5430de51b56SMatthew G. Knepley 
544636b9322SMatthew G. Knepley   4) Volumetric input mesh with boundary output mesh
545636b9322SMatthew G. Knepley 
546636b9322SMatthew G. Knepley      Here we must get a subspace for the input DS
547636b9322SMatthew G. Knepley 
5480de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
5490de51b56SMatthew G. Knepley 
5500de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
5510de51b56SMatthew G. Knepley 
5520de51b56SMatthew 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.
5530de51b56SMatthew 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
5545f790a90SMatthew 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
5550de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
5560de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
5570de51b56SMatthew G. Knepley 
5580de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
5590de51b56SMatthew G. Knepley 
5600de51b56SMatthew 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.
5610de51b56SMatthew G. Knepley */
56246fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
5631c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
5648c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
5658c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
5668c6c5593SMatthew G. Knepley {
567a6e0b375SMatthew G. Knepley   DM                 plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
568a6e0b375SMatthew G. Knepley   DMEnclosureType    encIn, encAux;
569a6e0b375SMatthew G. Knepley   PetscDS            ds = NULL, dsIn = NULL, dsAux = NULL;
570ca3d3a14SMatthew G. Knepley   Vec                localA = NULL, tv;
571aa0eca99SMatthew G. Knepley   IS                 fieldIS;
57247923291SMatthew G. Knepley   PetscSection       section;
573636b9322SMatthew G. Knepley   PetscDualSpace    *sp, *cellsp, *spIn, *cellspIn;
574ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
5758c6c5593SMatthew G. Knepley   PetscInt          *Nc;
576636b9322SMatthew G. Knepley   PetscInt           dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
5775fedec97SMatthew G. Knepley   PetscBool         *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, transform;
578c330f8ffSToby Isaac   DMField            coordField;
579c330f8ffSToby Isaac   DMLabel            depthLabel;
580c330f8ffSToby Isaac   PetscQuadrature    allPoints = NULL;
58147923291SMatthew G. Knepley 
58247923291SMatthew G. Knepley   PetscFunctionBegin;
5839566063dSJacob Faibussowitsch   if (localU) PetscCall(VecGetDM(localU, &dmIn));
584a6e0b375SMatthew G. Knepley   else        {dmIn = dm;}
5859566063dSJacob Faibussowitsch   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
5869566063dSJacob Faibussowitsch   if (localA) PetscCall(VecGetDM(localA, &dmAux)); else {dmAux = NULL;}
5879566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
5889566063dSJacob Faibussowitsch   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
5899566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
5909566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
5919566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5929566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
5939566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
5949566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
5959566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
5960de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
597a6e0b375SMatthew G. Knepley   if (dmAux) {
5989566063dSJacob Faibussowitsch     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
599a6e0b375SMatthew G. Knepley     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
60028b400f6SJacob Faibussowitsch       PetscCheck(localA,PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector");
601636b9322SMatthew G. Knepley     }
602636b9322SMatthew G. Knepley   }
603*e04ae0b4SMatthew G. Knepley   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL));
604*e04ae0b4SMatthew G. Knepley   PetscCall(DMGetCoordinateField(dm,&coordField));
605*e04ae0b4SMatthew G. Knepley   /**** No collective calls below this point ****/
606636b9322SMatthew G. Knepley   /* Determine height for iteration of all meshes */
607636b9322SMatthew G. Knepley   {
608636b9322SMatthew G. Knepley     DMPolytopeType ct, ctIn, ctAux;
60919552267SMatthew G. Knepley     PetscInt       minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux;
61019552267SMatthew G. Knepley     PetscInt       dim = -1, dimIn = -1, dimAux = -1;
61188aca1feSMatthew G. Knepley 
6129566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
613636b9322SMatthew G. Knepley     if (pEnd > pStart) {
6149566063dSJacob Faibussowitsch       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
615636b9322SMatthew G. Knepley       p    = lStart < 0 ? pStart : lStart;
6169566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plex, p, &ct));
617636b9322SMatthew G. Knepley       dim  = DMPolytopeTypeGetDim(ct);
6189566063dSJacob Faibussowitsch       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
6199566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
6209566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
621636b9322SMatthew G. Knepley       dimIn = DMPolytopeTypeGetDim(ctIn);
622636b9322SMatthew G. Knepley       if (dmAux) {
6239566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
62419552267SMatthew G. Knepley         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
62519552267SMatthew G. Knepley         if (pStartAux < pEndAux) {
6269566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
627636b9322SMatthew G. Knepley           dimAux = DMPolytopeTypeGetDim(ctAux);
62819552267SMatthew G. Knepley         }
629636b9322SMatthew G. Knepley       } else dimAux = dim;
630*e04ae0b4SMatthew G. Knepley     } else {
631*e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plex));
632*e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plexIn));
633*e04ae0b4SMatthew G. Knepley       if (dmAux) PetscCall(DMDestroy(&plexAux));
634*e04ae0b4SMatthew G. Knepley       PetscFunctionReturn(0);
635e39fcb4eSStefano Zampini     }
636636b9322SMatthew G. Knepley     if (dim < 0) {
637636b9322SMatthew G. Knepley       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
638636b9322SMatthew G. Knepley 
639636b9322SMatthew G. Knepley       /* Fall back to determination based on being a submesh */
6409566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plex,   &spmap));
6419566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
6429566063dSJacob Faibussowitsch       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
643636b9322SMatthew G. Knepley       dim    = spmap    ? 1 : 0;
644636b9322SMatthew G. Knepley       dimIn  = spmapIn  ? 1 : 0;
645636b9322SMatthew G. Knepley       dimAux = spmapAux ? 1 : 0;
64688aca1feSMatthew G. Knepley     }
647636b9322SMatthew G. Knepley     {
64819552267SMatthew G. Knepley       PetscInt dimProj   = PetscMin(PetscMin(dim, dimIn), (dimAux < 0 ? PETSC_MAX_INT : dimAux));
64919552267SMatthew G. Knepley       PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;
650636b9322SMatthew G. Knepley 
65119552267SMatthew 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");
652636b9322SMatthew G. Knepley       if (dimProj < dim) minHeight = 1;
653636b9322SMatthew G. Knepley       htInc    =  dim       - dimProj;
654636b9322SMatthew G. Knepley       htIncIn  =  dimIn     - dimProj;
65519552267SMatthew G. Knepley       htIncAux =  dimAuxEff - dimProj;
6560de51b56SMatthew G. Knepley     }
657a6e0b375SMatthew G. Knepley   }
6589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
6599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
6609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
6612af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
66263a3b9bcSJacob 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);
6639566063dSJacob Faibussowitsch   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, 0, NULL, &ds));
6649566063dSJacob Faibussowitsch   if (!ds) PetscCall(DMGetDS(dm, &ds));
6659566063dSJacob Faibussowitsch   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn));
6669566063dSJacob Faibussowitsch   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
6679566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
6699566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &NfTot));
6709566063dSJacob Faibussowitsch   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
6719566063dSJacob Faibussowitsch   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL));
6729566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6739566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
6749566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
6750c898c18SMatthew G. Knepley   if (dmAux) {
6769566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmAux, &dsAux));
6779566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6780c898c18SMatthew G. Knepley   }
6799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
6809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
6819566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
682636b9322SMatthew G. Knepley   else               {cellsp = sp; cellspIn = spIn;}
6838c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
68447923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
685665f567fSMatthew G. Knepley     PetscDiscType disctype;
68647923291SMatthew G. Knepley 
6879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
688665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
689665f567fSMatthew G. Knepley       PetscFE fe;
69047923291SMatthew G. Knepley 
69147923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
692665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
6939566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *) &fe));
6949566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
695665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
696665f567fSMatthew G. Knepley       PetscFV fv;
69747923291SMatthew G. Knepley 
69847923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
699665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
7009566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *) &fv));
7019566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
702665f567fSMatthew G. Knepley     } else {
703665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
704665f567fSMatthew G. Knepley       cellsp[f] = NULL;
705665f567fSMatthew G. Knepley     }
70647923291SMatthew G. Knepley   }
707636b9322SMatthew G. Knepley   for (f = 0; f < NfIn; ++f) {
708636b9322SMatthew G. Knepley     PetscDiscType disctype;
709636b9322SMatthew G. Knepley 
7109566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
711636b9322SMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
712636b9322SMatthew G. Knepley       PetscFE fe;
713636b9322SMatthew G. Knepley 
7149566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fe));
7159566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
716636b9322SMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
717636b9322SMatthew G. Knepley       PetscFV fv;
718636b9322SMatthew G. Knepley 
7199566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fv));
7209566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
721636b9322SMatthew G. Knepley     } else {
722636b9322SMatthew G. Knepley       cellspIn[f] = NULL;
723636b9322SMatthew G. Knepley     }
724636b9322SMatthew G. Knepley   }
725636b9322SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
726636b9322SMatthew G. Knepley     if (!htInc) {sp[f] = cellsp[f];}
7279566063dSJacob Faibussowitsch     else        PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
728636b9322SMatthew G. Knepley   }
729ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
73074fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
731665f567fSMatthew G. Knepley     PetscDiscType    disctype;
7324a3e9fdbSToby Isaac     const PetscReal *points;
7334a3e9fdbSToby Isaac     PetscInt         numPoints;
7340c898c18SMatthew G. Knepley 
73508401ef6SPierre Jolivet     PetscCheck(maxHeight <= minHeight,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
7369566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &allPoints));
7379566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
7389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
739a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
740636b9322SMatthew G. Knepley       if (!htIncIn) {spIn[f] = cellspIn[f];}
7419566063dSJacob Faibussowitsch       else          PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
742636b9322SMatthew G. Knepley 
7439566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
744665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
7459566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem));
746636b9322SMatthew G. Knepley       if (!htIncIn) {subfem = fem;}
7479566063dSJacob Faibussowitsch       else          PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
7489566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]));
7490c898c18SMatthew G. Knepley     }
7500c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
7519566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
752665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
7539566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem));
754636b9322SMatthew G. Knepley       if (!htIncAux) {subfem = fem;}
7559566063dSJacob Faibussowitsch       else           PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
7569566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]));
7570c898c18SMatthew G. Knepley     }
7580c898c18SMatthew G. Knepley   }
75947923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
7602af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
761636b9322SMatthew G. Knepley     PetscInt     hEff     = h - minHeight + htInc;
762636b9322SMatthew G. Knepley     PetscInt     hEffIn   = h - minHeight + htIncIn;
763636b9322SMatthew G. Knepley     PetscInt     hEffAux  = h - minHeight + htIncAux;
764a6e0b375SMatthew G. Knepley     PetscDS      dsEff    = ds;
765636b9322SMatthew G. Knepley     PetscDS      dsEffIn  = dsIn;
766636b9322SMatthew G. Knepley     PetscDS      dsEffAux = dsAux;
7678c6c5593SMatthew G. Knepley     PetscScalar *values;
768b7260050SToby Isaac     PetscBool   *fieldActive;
769b7260050SToby Isaac     PetscInt     maxDegree;
770e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
771c330f8ffSToby Isaac     IS           heightIS;
7728c6c5593SMatthew G. Knepley 
773636b9322SMatthew G. Knepley     if (h > minHeight) {
7749566063dSJacob Faibussowitsch       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
775636b9322SMatthew G. Knepley     }
7769566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
7779566063dSJacob Faibussowitsch     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
7789566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
779c330f8ffSToby Isaac     if (pEnd <= pStart) {
7809566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
781c330f8ffSToby Isaac       continue;
782c330f8ffSToby Isaac     }
78347923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
78447923291SMatthew G. Knepley     totDim = 0;
78547923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
7865fedec97SMatthew G. Knepley       PetscBool cohesive;
7875fedec97SMatthew G. Knepley 
788665f567fSMatthew G. Knepley       if (!sp[f]) continue;
7899566063dSJacob Faibussowitsch       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
7909566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
7919c3cf19fSMatthew G. Knepley       totDim += spDim;
7925fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) totDim += spDim;
79347923291SMatthew G. Knepley     }
794636b9322SMatthew G. Knepley     p    = lStart < 0 ? pStart : lStart;
7959566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
79663a3b9bcSJacob 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);
797c330f8ffSToby Isaac     if (!totDim) {
7989566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
799c330f8ffSToby Isaac       continue;
800c330f8ffSToby Isaac     }
8019566063dSJacob Faibussowitsch     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
802636b9322SMatthew G. Knepley     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
803636b9322SMatthew G. Knepley     if (localU) {
804636b9322SMatthew G. Knepley       PetscInt totDimIn, pIn, numValuesIn;
805636b9322SMatthew G. Knepley 
806636b9322SMatthew G. Knepley       totDimIn = 0;
807636b9322SMatthew G. Knepley       for (f = 0; f < NfIn; ++f) {
8085fedec97SMatthew G. Knepley         PetscBool cohesive;
8095fedec97SMatthew G. Knepley 
810636b9322SMatthew G. Knepley         if (!spIn[f]) continue;
8119566063dSJacob Faibussowitsch         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
8129566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
813636b9322SMatthew G. Knepley         totDimIn += spDim;
8145fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) totDimIn += spDim;
815636b9322SMatthew G. Knepley       }
8169566063dSJacob Faibussowitsch       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
8179566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
81863a3b9bcSJacob Faibussowitsch       PetscCheck(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);
8199566063dSJacob Faibussowitsch       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
820636b9322SMatthew G. Knepley     }
8219566063dSJacob Faibussowitsch     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
82247923291SMatthew G. Knepley     /* Loop over points at this height */
8239566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
8249566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
825aa0eca99SMatthew G. Knepley     {
826aa0eca99SMatthew G. Knepley       const PetscInt *fields;
827aa0eca99SMatthew G. Knepley 
8289566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
829aa0eca99SMatthew G. Knepley       for (f = 0; f < NfTot; ++f) {fieldActive[f] = PETSC_FALSE;}
830aa0eca99SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;}
8319566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
832aa0eca99SMatthew G. Knepley     }
8338c6c5593SMatthew G. Knepley     if (label) {
8348c6c5593SMatthew G. Knepley       PetscInt i;
83547923291SMatthew G. Knepley 
83647923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
837c330f8ffSToby Isaac         IS              pointIS, isectIS;
83847923291SMatthew G. Knepley         const PetscInt *points;
8398c6c5593SMatthew G. Knepley         PetscInt        n;
840c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
841c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
84247923291SMatthew G. Knepley 
8439566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
84447923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
8459566063dSJacob Faibussowitsch         PetscCall(ISIntersect(pointIS,heightIS,&isectIS));
8469566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&pointIS));
847c330f8ffSToby Isaac         if (!isectIS) continue;
8489566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(isectIS, &n));
8499566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(isectIS, &points));
8509566063dSJacob Faibussowitsch         PetscCall(DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree));
851b7260050SToby Isaac         if (maxDegree <= 1) {
8529566063dSJacob Faibussowitsch           PetscCall(DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad));
853c330f8ffSToby Isaac         }
854c330f8ffSToby Isaac         if (!quad) {
855c330f8ffSToby Isaac           if (!h && allPoints) {
856c330f8ffSToby Isaac             quad = allPoints;
857c330f8ffSToby Isaac             allPoints = NULL;
858c330f8ffSToby Isaac           } else {
8599566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf,sp,isCohesive ? dim-htInc-1 : dim-htInc,funcs,&quad));
860c330f8ffSToby Isaac           }
861c330f8ffSToby Isaac         }
8629566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
86347923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
86447923291SMatthew G. Knepley           const PetscInt  point = points[p];
86547923291SMatthew G. Knepley 
8669566063dSJacob Faibussowitsch           PetscCall(PetscArrayzero(values, numValues));
8679566063dSJacob Faibussowitsch           PetscCall(PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom));
8689566063dSJacob Faibussowitsch           PetscCall(DMPlexSetActivePoint(dm, point));
869d0609cedSBarry 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));
8709566063dSJacob Faibussowitsch           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
8719566063dSJacob Faibussowitsch           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
87247923291SMatthew G. Knepley         }
8739566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom));
8749566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&fegeom));
8759566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureDestroy(&quad));
8769566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(isectIS, &points));
8779566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isectIS));
87847923291SMatthew G. Knepley       }
8798c6c5593SMatthew G. Knepley     } else {
880c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
881c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
882c330f8ffSToby Isaac       IS              pointIS;
883c330f8ffSToby Isaac 
8849566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS));
8859566063dSJacob Faibussowitsch       PetscCall(DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree));
886b7260050SToby Isaac       if (maxDegree <= 1) {
8879566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad));
888c330f8ffSToby Isaac       }
889c330f8ffSToby Isaac       if (!quad) {
890c330f8ffSToby Isaac         if (!h && allPoints) {
891c330f8ffSToby Isaac           quad = allPoints;
892c330f8ffSToby Isaac           allPoints = NULL;
893c330f8ffSToby Isaac         } else {
8949566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &quad));
895c330f8ffSToby Isaac         }
896c330f8ffSToby Isaac       }
8979566063dSJacob Faibussowitsch       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
8988c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
8999566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(values, numValues));
9009566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom));
9019566063dSJacob Faibussowitsch         PetscCall(DMPlexSetActivePoint(dm, p));
902d0609cedSBarry 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));
9039566063dSJacob Faibussowitsch         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
9049566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
9058c6c5593SMatthew G. Knepley       }
9069566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom));
9079566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomDestroy(&fegeom));
9089566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
9099566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
9108c6c5593SMatthew G. Knepley     }
9119566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&heightIS));
9129566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
9139566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
91447923291SMatthew G. Knepley   }
9158c6c5593SMatthew G. Knepley   /* Cleanup */
916ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
9179566063dSJacob Faibussowitsch     for (f = 0; f < NfIn;  ++f) PetscCall(PetscTabulationDestroy(&T[f]));
9189566063dSJacob Faibussowitsch     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
9199566063dSJacob Faibussowitsch     PetscCall(PetscFree2(T, TAux));
9200c898c18SMatthew G. Knepley   }
9219566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&allPoints));
9229566063dSJacob Faibussowitsch   PetscCall(PetscFree3(isFE, sp, spIn));
9239566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
9249566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
9259566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plexIn));
9269566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMDestroy(&plexAux));
9278c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
92847923291SMatthew G. Knepley }
9298c6c5593SMatthew G. Knepley 
9308c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
9318c6c5593SMatthew G. Knepley {
9328c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9339566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX));
9348c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
9358c6c5593SMatthew G. Knepley }
9368c6c5593SMatthew G. Knepley 
9371c531cf8SMatthew G. Knepley 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)
9388c6c5593SMatthew G. Knepley {
9398c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9409566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX));
94147923291SMatthew G. Knepley   PetscFunctionReturn(0);
94247923291SMatthew G. Knepley }
94347923291SMatthew G. Knepley 
9448c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
94547923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
94647923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
94747923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
948191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
94947923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
95047923291SMatthew G. Knepley {
95147923291SMatthew G. Knepley   PetscFunctionBegin;
9529566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX));
9538c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
95447923291SMatthew G. Knepley }
95547923291SMatthew G. Knepley 
9561c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
9578c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
9588c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
9598c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
960191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
9618c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
9628c6c5593SMatthew G. Knepley {
9638c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9649566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX));
96547923291SMatthew G. Knepley   PetscFunctionReturn(0);
96647923291SMatthew G. Knepley }
967ece3a9fcSMatthew G. Knepley 
968ece3a9fcSMatthew G. Knepley PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
969ece3a9fcSMatthew G. Knepley                                                void (**funcs)(PetscInt, PetscInt, PetscInt,
970ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
971ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
972ece3a9fcSMatthew G. Knepley                                                               PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
973ece3a9fcSMatthew G. Knepley                                                InsertMode mode, Vec localX)
974ece3a9fcSMatthew G. Knepley {
975ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
9769566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void)) funcs, NULL, mode, localX));
977ece3a9fcSMatthew G. Knepley   PetscFunctionReturn(0);
978ece3a9fcSMatthew G. Knepley }
979