xref: /petsc/src/dm/impls/plex/plexproject.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 
181b32699bSMatthew G. Knepley .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 
381b32699bSMatthew G. Knepley .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 
68a6e0b375SMatthew G. Knepley .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));
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;
11898921bdaSJacob Faibussowitsch               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %D != %D 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           }
1259566063dSJacob Faibussowitsch           if (transform) {PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx)); v0 = x;}
1269566063dSJacob Faibussowitsch           PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx));
127c330f8ffSToby Isaac         }
1284bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
1299566063dSJacob Faibussowitsch         PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
1309566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
1319566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x));
1329566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
133c330f8ffSToby Isaac         v += spDim;
1345fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) {
13527f02ce8SMatthew G. Knepley           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
13627f02ce8SMatthew G. Knepley         }
137c330f8ffSToby Isaac       } else {
138c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
1399566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
140c330f8ffSToby Isaac         }
141c330f8ffSToby Isaac       }
142c330f8ffSToby Isaac     } else {
14327f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
1445fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
14527f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
14627f02ce8SMatthew G. Knepley       }
1478c6c5593SMatthew G. Knepley     }
1489c3cf19fSMatthew G. Knepley   }
1498c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1508c6c5593SMatthew G. Knepley }
1518c6c5593SMatthew G. Knepley 
152a6e0b375SMatthew G. Knepley /*
153a6e0b375SMatthew G. Knepley   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
154a6e0b375SMatthew G. Knepley 
155a6e0b375SMatthew G. Knepley   Input Parameters:
156a6e0b375SMatthew G. Knepley + dm             - The output DM
157a6e0b375SMatthew G. Knepley . ds             - The output DS
158a6e0b375SMatthew G. Knepley . dmIn           - The input DM
159a6e0b375SMatthew G. Knepley . dsIn           - The input DS
160a6e0b375SMatthew G. Knepley . dmAux          - The auxiliary DM, which is always for the input space
161a6e0b375SMatthew G. Knepley . dsAux          - The auxiliary DS, which is always for the input space
162a6e0b375SMatthew G. Knepley . time           - The time for this evaluation
163a6e0b375SMatthew G. Knepley . localU         - The local solution
164a6e0b375SMatthew G. Knepley . localA         - The local auziliary fields
165a6e0b375SMatthew G. Knepley . cgeom          - The FE geometry for this point
166a6e0b375SMatthew G. Knepley . sp             - The output PetscDualSpace for each field
167a6e0b375SMatthew G. Knepley . p              - The point in the output DM
168a5b23f4aSJose E. Roman . T              - Input basis and derivatives for each field tabulated on the quadrature points
169ef0bb6c7SMatthew G. Knepley . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
170a6e0b375SMatthew G. Knepley . funcs          - The evaluation function for each field
171a6e0b375SMatthew G. Knepley - ctxs           - The user context for each field
172a6e0b375SMatthew G. Knepley 
173a6e0b375SMatthew G. Knepley   Output Parameter:
174a6e0b375SMatthew G. Knepley . values         - The value for each dual basis vector in the output dual space
175a6e0b375SMatthew G. Knepley 
176a6e0b375SMatthew G. Knepley   Note: Not supported for FV
177a6e0b375SMatthew G. Knepley 
178a6e0b375SMatthew G. Knepley   Level: developer
179a6e0b375SMatthew G. Knepley 
180a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private()
181a6e0b375SMatthew G. Knepley */
182a6e0b375SMatthew 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,
183ef0bb6c7SMatthew G. Knepley                                                    PetscTabulation *T, PetscTabulation *TAux,
1848c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
1858c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1868c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
187191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
1888c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
1898c6c5593SMatthew G. Knepley {
1908c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
1914bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
1928c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
1938c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
194191494d9SMatthew G. Knepley   const PetscScalar *constants;
1958c6c5593SMatthew G. Knepley   PetscReal         *x;
196ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
1974bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
1984bee2e38SMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed;
199ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
2005fedec97SMatthew G. Knepley   PetscBool          isAffine, isCohesive, transform;
2018c6c5593SMatthew G. Knepley 
2028c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
2039566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
2059566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
2069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
2079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
2089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
2099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
2109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
2119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
2129566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
2139566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmIn, &section));
2149566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
2159566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
2168c6c5593SMatthew G. Knepley   if (dmAux) {
21744171101SMatthew G. Knepley     PetscInt subp;
2181ba34526SMatthew G. Knepley 
2199566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
2209566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2219566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
2229566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2239566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2249566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
2259566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
2268c6c5593SMatthew G. Knepley   }
2278c6c5593SMatthew G. Knepley   /* Get values for closure */
2284bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
22927f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
23027f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
2314bee2e38SMatthew G. Knepley   if (isAffine) {
2324bee2e38SMatthew G. Knepley     fegeom.v    = x;
2334bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
2344bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
2354bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
2364bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
2374bee2e38SMatthew G. Knepley   }
2388c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
239c330f8ffSToby Isaac     PetscQuadrature  allPoints;
240c330f8ffSToby Isaac     PetscInt         q, dim, numPoints;
241c330f8ffSToby Isaac     const PetscReal *points;
242c330f8ffSToby Isaac     PetscScalar     *pointEval;
2435fedec97SMatthew G. Knepley     PetscBool        cohesive;
244c330f8ffSToby Isaac     DM               dm;
245c330f8ffSToby Isaac 
2468c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2479566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
2489566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
249c330f8ffSToby Isaac     if (!funcs[f]) {
250be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
2515fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
25227f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
25327f02ce8SMatthew G. Knepley       }
254c330f8ffSToby Isaac       continue;
255c330f8ffSToby Isaac     }
2569566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f],&dm));
2579566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
2589566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL));
2599566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
2600c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
261c330f8ffSToby Isaac       if (isAffine) {
262ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
2631c531cf8SMatthew G. Knepley       } else {
2644bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp*dE];
2654bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp*dE*dE];
2664bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
2674bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
2688c6c5593SMatthew G. Knepley       }
2699566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
2709566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
2719566063dSJacob Faibussowitsch       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
272a6e0b375SMatthew 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]);
2731c531cf8SMatthew G. Knepley     }
2749566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
2759566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
276c330f8ffSToby Isaac     v += spDim;
27727f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
2785fedec97SMatthew G. Knepley     if (isCohesive && !cohesive) {
27927f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
28027f02ce8SMatthew G. Knepley     }
2818c6c5593SMatthew G. Knepley   }
2829566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
2839566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
2848c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2858c6c5593SMatthew G. Knepley }
2868c6c5593SMatthew G. Knepley 
287a6e0b375SMatthew 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,
288ef0bb6c7SMatthew G. Knepley                                                      PetscTabulation *T, PetscTabulation *TAux,
289b18799e0SMatthew G. Knepley                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
290b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
291b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
292b18799e0SMatthew G. Knepley                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
293b18799e0SMatthew G. Knepley                                                      PetscScalar values[])
294b18799e0SMatthew G. Knepley {
295b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
2964bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
297b18799e0SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
298b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
299b18799e0SMatthew G. Knepley   const PetscScalar *constants;
300b18799e0SMatthew G. Knepley   PetscReal         *x;
301ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
3029f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
3034bee2e38SMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed;
304ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
305b18799e0SMatthew G. Knepley   PetscBool          isAffine;
306b18799e0SMatthew G. Knepley 
307b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
308*08401ef6SPierre Jolivet   PetscCheck(dm == dmIn,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
3099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
3119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
3129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
3139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x));
3149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
3159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
3169566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
3179566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients));
318b18799e0SMatthew G. Knepley   if (dmAux) {
319a6e0b375SMatthew G. Knepley     PetscInt subp;
320b18799e0SMatthew G. Knepley 
3219566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
3229566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3239566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
3249566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3259566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3269566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
3279566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
328b18799e0SMatthew G. Knepley   }
329b18799e0SMatthew G. Knepley   /* Get values for closure */
3304bee2e38SMatthew G. Knepley   isAffine = fgeom->isAffine;
331ea78f98cSLisandro Dalcin   fegeom.n  = NULL;
332ea78f98cSLisandro Dalcin   fegeom.J  = NULL;
333ea78f98cSLisandro Dalcin   fegeom.v  = NULL;
334ea78f98cSLisandro Dalcin   fegeom.xi = NULL;
335a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
336a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
3374bee2e38SMatthew G. Knepley   if (isAffine) {
3384bee2e38SMatthew G. Knepley     fegeom.v    = x;
3394bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
3404bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
3414bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
3424bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
3434bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
3449f209ee4SMatthew G. Knepley 
3459f209ee4SMatthew G. Knepley     cgeom.J     = fgeom->suppJ[0];
3469f209ee4SMatthew G. Knepley     cgeom.invJ  = fgeom->suppInvJ[0];
3479f209ee4SMatthew G. Knepley     cgeom.detJ  = fgeom->suppDetJ[0];
3484bee2e38SMatthew G. Knepley   }
349b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
350b18799e0SMatthew G. Knepley     PetscQuadrature   allPoints;
351b18799e0SMatthew G. Knepley     PetscInt          q, dim, numPoints;
352b18799e0SMatthew G. Knepley     const PetscReal   *points;
353b18799e0SMatthew G. Knepley     PetscScalar       *pointEval;
354b18799e0SMatthew G. Knepley     DM                dm;
355b18799e0SMatthew G. Knepley 
356b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
3579566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
358b18799e0SMatthew G. Knepley     if (!funcs[f]) {
359b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
360b18799e0SMatthew G. Knepley       continue;
361b18799e0SMatthew G. Knepley     }
3629566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f],&dm));
3639566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
3649566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL));
3659566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
366b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
367b18799e0SMatthew G. Knepley       if (isAffine) {
368ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
369b18799e0SMatthew G. Knepley       } else {
3703fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp*dE];
3719f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp*dE*dE];
3729f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp*dE*dE];
3734bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
3744bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp*dE];
3759f209ee4SMatthew G. Knepley 
3769f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][tp*dE*dE];
3779f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][tp*dE*dE];
3789f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][tp];
379b18799e0SMatthew G. Knepley       }
380a6e0b375SMatthew 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 */
3819566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
3829566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
3834bee2e38SMatthew 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]);
384b18799e0SMatthew G. Knepley     }
3859566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
3869566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval));
387b18799e0SMatthew G. Knepley     v += spDim;
388b18799e0SMatthew G. Knepley   }
3899566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients));
3909566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
391b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
392b18799e0SMatthew G. Knepley }
393b18799e0SMatthew G. Knepley 
394a6e0b375SMatthew 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[],
395ef0bb6c7SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux,
3968c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
3978c6c5593SMatthew G. Knepley {
3988c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
3998c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
4008c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
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:
412ef0bb6c7SMatthew G. Knepley     ierr = 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[],
4169566063dSJacob Faibussowitsch                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);PetscCall(ierr);break;
417b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
418ef0bb6c7SMatthew G. Knepley     ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
419b18799e0SMatthew G. Knepley                                           (void (**)(PetscInt, PetscInt, PetscInt,
420b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
421b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
4229566063dSJacob Faibussowitsch                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);PetscCall(ierr);break;
42398921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
4248c6c5593SMatthew G. Knepley   }
4258c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
4268c6c5593SMatthew G. Knepley }
4278c6c5593SMatthew G. Knepley 
428df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
429df5c1128SToby Isaac {
430df5c1128SToby Isaac   PetscReal      *points;
431df5c1128SToby Isaac   PetscInt       f, numPoints;
432df5c1128SToby Isaac 
433df5c1128SToby Isaac   PetscFunctionBegin;
434df5c1128SToby Isaac   numPoints = 0;
435df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
436df5c1128SToby Isaac     if (funcs[f]) {
437df5c1128SToby Isaac       PetscQuadrature fAllPoints;
438df5c1128SToby Isaac       PetscInt        fNumPoints;
439df5c1128SToby Isaac 
4409566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL));
4419566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
442df5c1128SToby Isaac       numPoints += fNumPoints;
443df5c1128SToby Isaac     }
444df5c1128SToby Isaac   }
4459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim*numPoints,&points));
446df5c1128SToby Isaac   numPoints = 0;
447df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
448df5c1128SToby Isaac     if (funcs[f]) {
449df5c1128SToby Isaac       PetscQuadrature fAllPoints;
45054ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
451df5c1128SToby Isaac       const PetscReal *fPoints;
452df5c1128SToby Isaac 
4539566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL));
4549566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
455*08401ef6SPierre Jolivet       PetscCheck(qdim == dim,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %D for dual basis does not match input dimension %D", qdim, dim);
456df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
457df5c1128SToby Isaac       numPoints += fNumPoints;
458df5c1128SToby Isaac     }
459df5c1128SToby Isaac   }
4609566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,allPoints));
4619566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL));
462df5c1128SToby Isaac   PetscFunctionReturn(0);
463df5c1128SToby Isaac }
464df5c1128SToby Isaac 
4655f15299fSJeremy L Thompson /*@C
4665f15299fSJeremy 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.
4675f15299fSJeremy L Thompson 
4685f15299fSJeremy L Thompson   Input Parameters:
4695f15299fSJeremy L Thompson   dm - the DM
4705f15299fSJeremy L Thompson   odm - the enclosing DM
4715f15299fSJeremy L Thompson   label - label for DM domain, or NULL for whole domain
4725f15299fSJeremy L Thompson   numIds - the number of ids
4735f15299fSJeremy L Thompson   ids - An array of the label ids in sequence for the domain
4745f15299fSJeremy L Thompson   height - Height of target cells in DMPlex topology
4755f15299fSJeremy L Thompson 
4765f15299fSJeremy L Thompson   Output Parameters:
4775f15299fSJeremy L Thompson   point - the first labeled point
4785f15299fSJeremy L Thompson   ds - the ds corresponding to the first labeled point
4795f15299fSJeremy L Thompson 
4805f15299fSJeremy L Thompson   Level: developer
481817da375SSatish Balay @*/
4825f15299fSJeremy L Thompson PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
483e5e52638SMatthew G. Knepley {
484a6e0b375SMatthew G. Knepley   DM              plex;
485a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
486e9f0ba4eSJed Brown   PetscInt        ls = -1;
487e5e52638SMatthew G. Knepley 
488e5e52638SMatthew G. Knepley   PetscFunctionBegin;
4895f15299fSJeremy L Thompson   if (point) *point = -1;
490e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
4919566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
4929566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
493e9f0ba4eSJed Brown   for (PetscInt i = 0; i < numIds; ++i) {
494e9f0ba4eSJed Brown     IS       labelIS;
495e9f0ba4eSJed Brown     PetscInt num_points, pStart, pEnd;
4969566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
497e9f0ba4eSJed Brown     if (!labelIS) continue; /* No points with that id on this process */
4989566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
4999566063dSJacob Faibussowitsch     PetscCall(ISGetSize(labelIS, &num_points));
500e9f0ba4eSJed Brown     if (num_points) {
501e5e52638SMatthew G. Knepley       const PetscInt *points;
5029566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(labelIS, &points));
503e9f0ba4eSJed Brown       for (PetscInt i=0; i<num_points; i++) {
504e9f0ba4eSJed Brown         PetscInt point;
5059566063dSJacob Faibussowitsch         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
506e9f0ba4eSJed Brown         if (pStart <= point && point < pEnd) {
507a6e0b375SMatthew G. Knepley           ls = point;
5089566063dSJacob Faibussowitsch           if (ds) PetscCall(DMGetCellDS(dm, ls, ds));
509e5e52638SMatthew G. Knepley         }
510e9f0ba4eSJed Brown       }
5119566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(labelIS, &points));
512e9f0ba4eSJed Brown     }
5139566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&labelIS));
514e5e52638SMatthew G. Knepley     if (ls >= 0) break;
515e5e52638SMatthew G. Knepley   }
5165f15299fSJeremy L Thompson   if (point) *point = ls;
5179566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
518e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
519e5e52638SMatthew G. Knepley }
520e5e52638SMatthew G. Knepley 
5210de51b56SMatthew G. Knepley /*
5220de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
5230de51b56SMatthew G. Knepley 
5240de51b56SMatthew G. Knepley   There are several different scenarios:
5250de51b56SMatthew G. Knepley 
5260de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
5270de51b56SMatthew G. Knepley 
5280de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
5290de51b56SMatthew G. Knepley 
5300de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
5310de51b56SMatthew G. Knepley 
5320de51b56SMatthew 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.
5330de51b56SMatthew G. Knepley 
5340de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
5350de51b56SMatthew G. Knepley 
5360de51b56SMatthew 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.
5370de51b56SMatthew G. Knepley 
538636b9322SMatthew G. Knepley   4) Volumetric input mesh with boundary output mesh
539636b9322SMatthew G. Knepley 
540636b9322SMatthew G. Knepley      Here we must get a subspace for the input DS
541636b9322SMatthew G. Knepley 
5420de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
5430de51b56SMatthew G. Knepley 
5440de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
5450de51b56SMatthew G. Knepley 
5460de51b56SMatthew 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.
5470de51b56SMatthew 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
5485f790a90SMatthew 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
5490de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
5500de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
5510de51b56SMatthew G. Knepley 
5520de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
5530de51b56SMatthew G. Knepley 
5540de51b56SMatthew 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.
5550de51b56SMatthew G. Knepley */
55646fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
5571c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
5588c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
5598c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
5608c6c5593SMatthew G. Knepley {
561a6e0b375SMatthew G. Knepley   DM                 plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
562a6e0b375SMatthew G. Knepley   DMEnclosureType    encIn, encAux;
563a6e0b375SMatthew G. Knepley   PetscDS            ds = NULL, dsIn = NULL, dsAux = NULL;
564ca3d3a14SMatthew G. Knepley   Vec                localA = NULL, tv;
565aa0eca99SMatthew G. Knepley   IS                 fieldIS;
56647923291SMatthew G. Knepley   PetscSection       section;
567636b9322SMatthew G. Knepley   PetscDualSpace    *sp, *cellsp, *spIn, *cellspIn;
568ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
5698c6c5593SMatthew G. Knepley   PetscInt          *Nc;
570636b9322SMatthew G. Knepley   PetscInt           dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
5715fedec97SMatthew G. Knepley   PetscBool         *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, transform;
572c330f8ffSToby Isaac   DMField            coordField;
573c330f8ffSToby Isaac   DMLabel            depthLabel;
574c330f8ffSToby Isaac   PetscQuadrature    allPoints = NULL;
57547923291SMatthew G. Knepley   PetscErrorCode     ierr;
57647923291SMatthew G. Knepley 
57747923291SMatthew G. Knepley   PetscFunctionBegin;
5789566063dSJacob Faibussowitsch   if (localU) PetscCall(VecGetDM(localU, &dmIn));
579a6e0b375SMatthew G. Knepley   else        {dmIn = dm;}
5809566063dSJacob Faibussowitsch   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
5819566063dSJacob Faibussowitsch   if (localA) PetscCall(VecGetDM(localA, &dmAux)); else {dmAux = NULL;}
5829566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
5839566063dSJacob Faibussowitsch   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
5849566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
5859566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
5869566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5879566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
5889566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
5899566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
5909566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
5910de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
592a6e0b375SMatthew G. Knepley   if (dmAux) {
5939566063dSJacob Faibussowitsch     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
594a6e0b375SMatthew G. Knepley     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
59528b400f6SJacob Faibussowitsch       PetscCheck(localA,PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector");
596636b9322SMatthew G. Knepley     }
597636b9322SMatthew G. Knepley   }
598636b9322SMatthew G. Knepley   /* Determine height for iteration of all meshes */
599636b9322SMatthew G. Knepley   {
600636b9322SMatthew G. Knepley     DMPolytopeType ct, ctIn, ctAux;
601636b9322SMatthew G. Knepley     PetscInt       minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux;
602636b9322SMatthew G. Knepley     PetscInt       dim = -1, dimIn, dimAux;
60388aca1feSMatthew G. Knepley 
6049566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
605636b9322SMatthew G. Knepley     if (pEnd > pStart) {
6069566063dSJacob Faibussowitsch       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
607636b9322SMatthew G. Knepley       p    = lStart < 0 ? pStart : lStart;
6089566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plex, p, &ct));
609636b9322SMatthew G. Knepley       dim  = DMPolytopeTypeGetDim(ct);
6109566063dSJacob Faibussowitsch       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
6119566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
6129566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
613636b9322SMatthew G. Knepley       dimIn = DMPolytopeTypeGetDim(ctIn);
614636b9322SMatthew G. Knepley       if (dmAux) {
6159566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
6169566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, NULL));
6179566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
618636b9322SMatthew G. Knepley         dimAux = DMPolytopeTypeGetDim(ctAux);
619636b9322SMatthew G. Knepley       } else dimAux = dim;
620e39fcb4eSStefano Zampini     }
621636b9322SMatthew G. Knepley     if (dim < 0) {
622636b9322SMatthew G. Knepley       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
623636b9322SMatthew G. Knepley 
624636b9322SMatthew G. Knepley       /* Fall back to determination based on being a submesh */
6259566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plex,   &spmap));
6269566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
6279566063dSJacob Faibussowitsch       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
628636b9322SMatthew G. Knepley       dim    = spmap    ? 1 : 0;
629636b9322SMatthew G. Knepley       dimIn  = spmapIn  ? 1 : 0;
630636b9322SMatthew G. Knepley       dimAux = spmapAux ? 1 : 0;
63188aca1feSMatthew G. Knepley     }
632636b9322SMatthew G. Knepley     {
633636b9322SMatthew G. Knepley       PetscInt dimProj = PetscMin(PetscMin(dim, dimIn), dimAux);
634636b9322SMatthew G. Knepley 
6352c71b3e2SJacob Faibussowitsch       PetscCheckFalse(PetscAbsInt(dimProj - dim) > 1 || PetscAbsInt(dimProj - dimIn) > 1 || PetscAbsInt(dimProj - dimAux) > 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not currently support differences of more than 1 in dimension");
636636b9322SMatthew G. Knepley       if (dimProj < dim) minHeight = 1;
637636b9322SMatthew G. Knepley       htInc    =  dim    - dimProj;
638636b9322SMatthew G. Knepley       htIncIn  =  dimIn  - dimProj;
639636b9322SMatthew G. Knepley       htIncAux =  dimAux - dimProj;
6400de51b56SMatthew G. Knepley     }
641a6e0b375SMatthew G. Knepley   }
6429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
6439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
6449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
6452af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
6462c71b3e2SJacob Faibussowitsch   PetscCheck(maxHeight >= 0 && maxHeight <= dim,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
6479566063dSJacob Faibussowitsch   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, 0, NULL, &ds));
6489566063dSJacob Faibussowitsch   if (!ds) PetscCall(DMGetDS(dm, &ds));
6499566063dSJacob Faibussowitsch   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn));
6509566063dSJacob Faibussowitsch   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
6519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6529566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
6539566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &NfTot));
6549566063dSJacob Faibussowitsch   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
6559566063dSJacob Faibussowitsch   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL));
6569566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6579566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
6589566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
6590c898c18SMatthew G. Knepley   if (dmAux) {
6609566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmAux, &dsAux));
6619566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6620c898c18SMatthew G. Knepley   }
6639566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
6649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
6659566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
666636b9322SMatthew G. Knepley   else               {cellsp = sp; cellspIn = spIn;}
6679566063dSJacob Faibussowitsch   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL));
6688c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
66947923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
670665f567fSMatthew G. Knepley     PetscDiscType disctype;
67147923291SMatthew G. Knepley 
6729566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
673665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
674665f567fSMatthew G. Knepley       PetscFE fe;
67547923291SMatthew G. Knepley 
67647923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
677665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
6789566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *) &fe));
6799566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
680665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
681665f567fSMatthew G. Knepley       PetscFV fv;
68247923291SMatthew G. Knepley 
68347923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
684665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
6859566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *) &fv));
6869566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
687665f567fSMatthew G. Knepley     } else {
688665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
689665f567fSMatthew G. Knepley       cellsp[f] = NULL;
690665f567fSMatthew G. Knepley     }
69147923291SMatthew G. Knepley   }
692636b9322SMatthew G. Knepley   for (f = 0; f < NfIn; ++f) {
693636b9322SMatthew G. Knepley     PetscDiscType disctype;
694636b9322SMatthew G. Knepley 
6959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
696636b9322SMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
697636b9322SMatthew G. Knepley       PetscFE fe;
698636b9322SMatthew G. Knepley 
6999566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fe));
7009566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
701636b9322SMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
702636b9322SMatthew G. Knepley       PetscFV fv;
703636b9322SMatthew G. Knepley 
7049566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fv));
7059566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
706636b9322SMatthew G. Knepley     } else {
707636b9322SMatthew G. Knepley       cellspIn[f] = NULL;
708636b9322SMatthew G. Knepley     }
709636b9322SMatthew G. Knepley   }
7109566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm,&coordField));
711636b9322SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
712636b9322SMatthew G. Knepley     if (!htInc) {sp[f] = cellsp[f];}
7139566063dSJacob Faibussowitsch     else        PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
714636b9322SMatthew G. Knepley   }
715ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
71674fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
717665f567fSMatthew G. Knepley     PetscDiscType    disctype;
7184a3e9fdbSToby Isaac     const PetscReal *points;
7194a3e9fdbSToby Isaac     PetscInt         numPoints;
7200c898c18SMatthew G. Knepley 
721*08401ef6SPierre Jolivet     PetscCheck(maxHeight <= minHeight,PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
7229566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &allPoints));
7239566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
7249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
725a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
726636b9322SMatthew G. Knepley       if (!htIncIn) {spIn[f] = cellspIn[f];}
7279566063dSJacob Faibussowitsch       else          PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
728636b9322SMatthew G. Knepley 
7299566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
730665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
7319566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem));
732636b9322SMatthew G. Knepley       if (!htIncIn) {subfem = fem;}
7339566063dSJacob Faibussowitsch       else          PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
7349566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]));
7350c898c18SMatthew G. Knepley     }
7360c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
7379566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
738665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
7399566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem));
740636b9322SMatthew G. Knepley       if (!htIncAux) {subfem = fem;}
7419566063dSJacob Faibussowitsch       else           PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
7429566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]));
7430c898c18SMatthew G. Knepley     }
7440c898c18SMatthew G. Knepley   }
74547923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
7462af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
747636b9322SMatthew G. Knepley     PetscInt     hEff     = h - minHeight + htInc;
748636b9322SMatthew G. Knepley     PetscInt     hEffIn   = h - minHeight + htIncIn;
749636b9322SMatthew G. Knepley     PetscInt     hEffAux  = h - minHeight + htIncAux;
750a6e0b375SMatthew G. Knepley     PetscDS      dsEff    = ds;
751636b9322SMatthew G. Knepley     PetscDS      dsEffIn  = dsIn;
752636b9322SMatthew G. Knepley     PetscDS      dsEffAux = dsAux;
7538c6c5593SMatthew G. Knepley     PetscScalar *values;
754b7260050SToby Isaac     PetscBool   *fieldActive;
755b7260050SToby Isaac     PetscInt     maxDegree;
756e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
757c330f8ffSToby Isaac     IS           heightIS;
7588c6c5593SMatthew G. Knepley 
759636b9322SMatthew G. Knepley     if (h > minHeight) {
7609566063dSJacob Faibussowitsch       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
761636b9322SMatthew G. Knepley     }
7629566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
7639566063dSJacob Faibussowitsch     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
7649566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
765c330f8ffSToby Isaac     if (pEnd <= pStart) {
7669566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
767c330f8ffSToby Isaac       continue;
768c330f8ffSToby Isaac     }
76947923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
77047923291SMatthew G. Knepley     totDim = 0;
77147923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
7725fedec97SMatthew G. Knepley       PetscBool cohesive;
7735fedec97SMatthew G. Knepley 
774665f567fSMatthew G. Knepley       if (!sp[f]) continue;
7759566063dSJacob Faibussowitsch       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
7769566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
7779c3cf19fSMatthew G. Knepley       totDim += spDim;
7785fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) totDim += spDim;
77947923291SMatthew G. Knepley     }
780636b9322SMatthew G. Knepley     p    = lStart < 0 ? pStart : lStart;
7819566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
782*08401ef6SPierre Jolivet     PetscCheck(numValues == totDim,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The output section point (%D) closure size %D != dual space dimension %D at height %D in [%D, %D]", p, numValues, totDim, h, minHeight, maxHeight);
783c330f8ffSToby Isaac     if (!totDim) {
7849566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
785c330f8ffSToby Isaac       continue;
786c330f8ffSToby Isaac     }
7879566063dSJacob Faibussowitsch     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
788636b9322SMatthew G. Knepley     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
789636b9322SMatthew G. Knepley     if (localU) {
790636b9322SMatthew G. Knepley       PetscInt totDimIn, pIn, numValuesIn;
791636b9322SMatthew G. Knepley 
792636b9322SMatthew G. Knepley       totDimIn = 0;
793636b9322SMatthew G. Knepley       for (f = 0; f < NfIn; ++f) {
7945fedec97SMatthew G. Knepley         PetscBool cohesive;
7955fedec97SMatthew G. Knepley 
796636b9322SMatthew G. Knepley         if (!spIn[f]) continue;
7979566063dSJacob Faibussowitsch         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
7989566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
799636b9322SMatthew G. Knepley         totDimIn += spDim;
8005fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) totDimIn += spDim;
801636b9322SMatthew G. Knepley       }
8029566063dSJacob Faibussowitsch       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
8039566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
804*08401ef6SPierre Jolivet       PetscCheck(numValuesIn == totDimIn,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The input section point (%D) closure size %D != dual space dimension %D at height %D", pIn, numValuesIn, totDimIn, htIncIn);
8059566063dSJacob Faibussowitsch       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
806636b9322SMatthew G. Knepley     }
8079566063dSJacob Faibussowitsch     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
80847923291SMatthew G. Knepley     /* Loop over points at this height */
8099566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
8109566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
811aa0eca99SMatthew G. Knepley     {
812aa0eca99SMatthew G. Knepley       const PetscInt *fields;
813aa0eca99SMatthew G. Knepley 
8149566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
815aa0eca99SMatthew G. Knepley       for (f = 0; f < NfTot; ++f) {fieldActive[f] = PETSC_FALSE;}
816aa0eca99SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;}
8179566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
818aa0eca99SMatthew G. Knepley     }
8198c6c5593SMatthew G. Knepley     if (label) {
8208c6c5593SMatthew G. Knepley       PetscInt i;
82147923291SMatthew G. Knepley 
82247923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
823c330f8ffSToby Isaac         IS              pointIS, isectIS;
82447923291SMatthew G. Knepley         const PetscInt *points;
8258c6c5593SMatthew G. Knepley         PetscInt        n;
826c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
827c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
82847923291SMatthew G. Knepley 
8299566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
83047923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
8319566063dSJacob Faibussowitsch         PetscCall(ISIntersect(pointIS,heightIS,&isectIS));
8329566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&pointIS));
833c330f8ffSToby Isaac         if (!isectIS) continue;
8349566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(isectIS, &n));
8359566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(isectIS, &points));
8369566063dSJacob Faibussowitsch         PetscCall(DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree));
837b7260050SToby Isaac         if (maxDegree <= 1) {
8389566063dSJacob Faibussowitsch           PetscCall(DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad));
839c330f8ffSToby Isaac         }
840c330f8ffSToby Isaac         if (!quad) {
841c330f8ffSToby Isaac           if (!h && allPoints) {
842c330f8ffSToby Isaac             quad = allPoints;
843c330f8ffSToby Isaac             allPoints = NULL;
844c330f8ffSToby Isaac           } else {
8459566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf,sp,isCohesive ? dim-htInc-1 : dim-htInc,funcs,&quad));
846c330f8ffSToby Isaac           }
847c330f8ffSToby Isaac         }
8489566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
84947923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
85047923291SMatthew G. Knepley           const PetscInt  point = points[p];
85147923291SMatthew G. Knepley 
8529566063dSJacob Faibussowitsch           PetscCall(PetscArrayzero(values, numValues));
8539566063dSJacob Faibussowitsch           PetscCall(PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom));
8549566063dSJacob Faibussowitsch           PetscCall(DMPlexSetActivePoint(dm, point));
855636b9322SMatthew G. Knepley           ierr = 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);
85647923291SMatthew G. Knepley           if (ierr) {
8579566063dSJacob Faibussowitsch             PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
8589566063dSJacob Faibussowitsch             PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
8599566063dSJacob Faibussowitsch             PetscCall(ierr);
86047923291SMatthew G. Knepley           }
8619566063dSJacob Faibussowitsch           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
8629566063dSJacob Faibussowitsch           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
86347923291SMatthew G. Knepley         }
8649566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom));
8659566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&fegeom));
8669566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureDestroy(&quad));
8679566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(isectIS, &points));
8689566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isectIS));
86947923291SMatthew G. Knepley       }
8708c6c5593SMatthew G. Knepley     } else {
871c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
872c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
873c330f8ffSToby Isaac       IS              pointIS;
874c330f8ffSToby Isaac 
8759566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS));
8769566063dSJacob Faibussowitsch       PetscCall(DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree));
877b7260050SToby Isaac       if (maxDegree <= 1) {
8789566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad));
879c330f8ffSToby Isaac       }
880c330f8ffSToby Isaac       if (!quad) {
881c330f8ffSToby Isaac         if (!h && allPoints) {
882c330f8ffSToby Isaac           quad = allPoints;
883c330f8ffSToby Isaac           allPoints = NULL;
884c330f8ffSToby Isaac         } else {
8859566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &quad));
886c330f8ffSToby Isaac         }
887c330f8ffSToby Isaac       }
8889566063dSJacob Faibussowitsch       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
8898c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
8909566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(values, numValues));
8919566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom));
8929566063dSJacob Faibussowitsch         PetscCall(DMPlexSetActivePoint(dm, p));
893636b9322SMatthew G. Knepley         ierr = 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);
8948c6c5593SMatthew G. Knepley         if (ierr) {
8959566063dSJacob Faibussowitsch           PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
8969566063dSJacob Faibussowitsch           PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
8979566063dSJacob Faibussowitsch           PetscCall(ierr);
8988c6c5593SMatthew G. Knepley         }
8999566063dSJacob Faibussowitsch         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
9009566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
9018c6c5593SMatthew G. Knepley       }
9029566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom));
9039566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomDestroy(&fegeom));
9049566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
9059566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
9068c6c5593SMatthew G. Knepley     }
9079566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&heightIS));
9089566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
9099566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
91047923291SMatthew G. Knepley   }
9118c6c5593SMatthew G. Knepley   /* Cleanup */
912ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
9139566063dSJacob Faibussowitsch     for (f = 0; f < NfIn;  ++f) PetscCall(PetscTabulationDestroy(&T[f]));
9149566063dSJacob Faibussowitsch     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
9159566063dSJacob Faibussowitsch     PetscCall(PetscFree2(T, TAux));
9160c898c18SMatthew G. Knepley   }
9179566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&allPoints));
9189566063dSJacob Faibussowitsch   PetscCall(PetscFree3(isFE, sp, spIn));
9199566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
9209566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
9219566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plexIn));
9229566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMDestroy(&plexAux));
9238c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
92447923291SMatthew G. Knepley }
9258c6c5593SMatthew G. Knepley 
9268c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
9278c6c5593SMatthew G. Knepley {
9288c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9299566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX));
9308c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
9318c6c5593SMatthew G. Knepley }
9328c6c5593SMatthew G. Knepley 
9331c531cf8SMatthew 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)
9348c6c5593SMatthew G. Knepley {
9358c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9369566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX));
93747923291SMatthew G. Knepley   PetscFunctionReturn(0);
93847923291SMatthew G. Knepley }
93947923291SMatthew G. Knepley 
9408c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
94147923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
94247923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
94347923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
944191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
94547923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
94647923291SMatthew G. Knepley {
94747923291SMatthew G. Knepley   PetscFunctionBegin;
9489566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX));
9498c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
95047923291SMatthew G. Knepley }
95147923291SMatthew G. Knepley 
9521c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
9538c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
9548c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
9558c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
956191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
9578c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
9588c6c5593SMatthew G. Knepley {
9598c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9609566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX));
96147923291SMatthew G. Knepley   PetscFunctionReturn(0);
96247923291SMatthew G. Knepley }
963ece3a9fcSMatthew G. Knepley 
964ece3a9fcSMatthew G. Knepley PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
965ece3a9fcSMatthew G. Knepley                                                void (**funcs)(PetscInt, PetscInt, PetscInt,
966ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
967ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
968ece3a9fcSMatthew G. Knepley                                                               PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
969ece3a9fcSMatthew G. Knepley                                                InsertMode mode, Vec localX)
970ece3a9fcSMatthew G. Knepley {
971ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
9729566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void)) funcs, NULL, mode, localX));
973ece3a9fcSMatthew G. Knepley   PetscFunctionReturn(0);
974ece3a9fcSMatthew G. Knepley }
975