xref: /petsc/src/dm/impls/plex/plexproject.c (revision e9f0ba4e0790d84194ea0735b7d9de705f1a5366)
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   PetscErrorCode ierr;
778c6c5593SMatthew G. Knepley 
788c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
79a6e0b375SMatthew G. Knepley   ierr = DMGetCoordinateDim(dmIn, &coordDim);CHKERRQ(ierr);
80a6e0b375SMatthew G. Knepley   ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr);
81a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
82a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
835fedec97SMatthew G. Knepley   ierr = PetscDSIsCohesive(ds, &isCohesive);CHKERRQ(ierr);
848c6c5593SMatthew G. Knepley   /* Get values for closure */
85c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
86c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
878c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
885fedec97SMatthew G. Knepley     PetscBool    cohesive;
898c6c5593SMatthew G. Knepley 
908c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
915fedec97SMatthew G. Knepley     ierr = PetscDSGetCohesive(ds, f, &cohesive);CHKERRQ(ierr);
928c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
938c6c5593SMatthew G. Knepley     if (funcs[f]) {
94c330f8ffSToby Isaac       if (isFE[f]) {
95c330f8ffSToby Isaac         PetscQuadrature   allPoints;
96c330f8ffSToby Isaac         PetscInt          q, dim, numPoints;
97c330f8ffSToby Isaac         const PetscReal   *points;
98c330f8ffSToby Isaac         PetscScalar       *pointEval;
99c330f8ffSToby Isaac         PetscReal         *x;
100ca3d3a14SMatthew G. Knepley         DM                rdm;
101c330f8ffSToby Isaac 
102ca3d3a14SMatthew G. Knepley         ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr);
103b4457527SToby Isaac         ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
104c330f8ffSToby Isaac         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
105ca3d3a14SMatthew G. Knepley         ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
106ca3d3a14SMatthew G. Knepley         ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
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;
119665f567fSMatthew G. Knepley               } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %D != %D 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           }
126a6e0b375SMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;}
127c330f8ffSToby Isaac           ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
128c330f8ffSToby Isaac         }
1294bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
1302edcad52SToby Isaac         ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr);
131c330f8ffSToby Isaac         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
132ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
133ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
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) {
140c330f8ffSToby Isaac           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
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 
181a6e0b375SMatthew G. Knepley .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   PetscErrorCode     ierr;
2038c6c5593SMatthew G. Knepley 
2048c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
205a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
206a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
2075fedec97SMatthew G. Knepley   ierr = PetscDSIsCohesive(ds, &isCohesive);CHKERRQ(ierr);
208a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
209a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(dsIn, &uOff);CHKERRQ(ierr);
210a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x);CHKERRQ(ierr);
211a6e0b375SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
212a6e0b375SMatthew G. Knepley   ierr = PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
213a6e0b375SMatthew G. Knepley   ierr = PetscDSGetConstants(dsIn, &numConstants, &constants);CHKERRQ(ierr);
214a6e0b375SMatthew G. Knepley   ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr);
215a6e0b375SMatthew G. Knepley   ierr = DMGetLocalSection(dmIn, &section);CHKERRQ(ierr);
216a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp);CHKERRQ(ierr);
217a6e0b375SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
2188c6c5593SMatthew G. Knepley   if (dmAux) {
21944171101SMatthew G. Knepley     PetscInt subp;
2201ba34526SMatthew G. Knepley 
221a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
222a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
22392fd8e1eSJed Brown     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
224a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
225a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
226a6e0b375SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
2271ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
2288c6c5593SMatthew G. Knepley   }
2298c6c5593SMatthew G. Knepley   /* Get values for closure */
2304bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
23127f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
23227f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
2334bee2e38SMatthew G. Knepley   if (isAffine) {
2344bee2e38SMatthew G. Knepley     fegeom.v    = x;
2354bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
2364bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
2374bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
2384bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
2394bee2e38SMatthew G. Knepley   }
2408c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
241c330f8ffSToby Isaac     PetscQuadrature  allPoints;
242c330f8ffSToby Isaac     PetscInt         q, dim, numPoints;
243c330f8ffSToby Isaac     const PetscReal *points;
244c330f8ffSToby Isaac     PetscScalar     *pointEval;
2455fedec97SMatthew G. Knepley     PetscBool        cohesive;
246c330f8ffSToby Isaac     DM               dm;
247c330f8ffSToby Isaac 
2488c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2495fedec97SMatthew G. Knepley     ierr = PetscDSGetCohesive(ds, f, &cohesive);CHKERRQ(ierr);
2508c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
251c330f8ffSToby Isaac     if (!funcs[f]) {
252be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
2535fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
25427f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
25527f02ce8SMatthew G. Knepley       }
256c330f8ffSToby Isaac       continue;
257c330f8ffSToby Isaac     }
258c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
259b4457527SToby Isaac     ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
260c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
261c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
2620c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
263c330f8ffSToby Isaac       if (isAffine) {
264ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
2651c531cf8SMatthew G. Knepley       } else {
2664bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp*dE];
2674bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp*dE*dE];
2684bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
2694bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
2708c6c5593SMatthew G. Knepley       }
271ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
272ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
273a6e0b375SMatthew G. Knepley       if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);}
274a6e0b375SMatthew 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]);
2751c531cf8SMatthew G. Knepley     }
276c330f8ffSToby Isaac     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
277c330f8ffSToby Isaac     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
278c330f8ffSToby Isaac     v += spDim;
27927f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
2805fedec97SMatthew G. Knepley     if (isCohesive && !cohesive) {
28127f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
28227f02ce8SMatthew G. Knepley     }
2838c6c5593SMatthew G. Knepley   }
284a6e0b375SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
2858c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
2868c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2878c6c5593SMatthew G. Knepley }
2888c6c5593SMatthew G. Knepley 
289a6e0b375SMatthew 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,
290ef0bb6c7SMatthew G. Knepley                                                      PetscTabulation *T, PetscTabulation *TAux,
291b18799e0SMatthew G. Knepley                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
292b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
293b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
294b18799e0SMatthew G. Knepley                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
295b18799e0SMatthew G. Knepley                                                      PetscScalar values[])
296b18799e0SMatthew G. Knepley {
297b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
2984bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
299b18799e0SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
300b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
301b18799e0SMatthew G. Knepley   const PetscScalar *constants;
302b18799e0SMatthew G. Knepley   PetscReal         *x;
303ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
3049f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
3054bee2e38SMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed;
306ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
307b18799e0SMatthew G. Knepley   PetscBool          isAffine;
308b18799e0SMatthew G. Knepley   PetscErrorCode     ierr;
309b18799e0SMatthew G. Knepley 
310b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
311a6e0b375SMatthew G. Knepley   if (dm != dmIn) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
312a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
313a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
314a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
315a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
316a6e0b375SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
317a6e0b375SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
318a6e0b375SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
31992fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
320a6e0b375SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
321b18799e0SMatthew G. Knepley   if (dmAux) {
322a6e0b375SMatthew G. Knepley     PetscInt subp;
323b18799e0SMatthew G. Knepley 
324a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
325a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
32692fd8e1eSJed Brown     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
327a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
328a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
329a6e0b375SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
330b18799e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
331b18799e0SMatthew G. Knepley   }
332b18799e0SMatthew G. Knepley   /* Get values for closure */
3334bee2e38SMatthew G. Knepley   isAffine = fgeom->isAffine;
334ea78f98cSLisandro Dalcin   fegeom.n  = NULL;
335ea78f98cSLisandro Dalcin   fegeom.J  = NULL;
336ea78f98cSLisandro Dalcin   fegeom.v  = NULL;
337ea78f98cSLisandro Dalcin   fegeom.xi = NULL;
338a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
339a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
3404bee2e38SMatthew G. Knepley   if (isAffine) {
3414bee2e38SMatthew G. Knepley     fegeom.v    = x;
3424bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
3434bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
3444bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
3454bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
3464bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
3479f209ee4SMatthew G. Knepley 
3489f209ee4SMatthew G. Knepley     cgeom.J     = fgeom->suppJ[0];
3499f209ee4SMatthew G. Knepley     cgeom.invJ  = fgeom->suppInvJ[0];
3509f209ee4SMatthew G. Knepley     cgeom.detJ  = fgeom->suppDetJ[0];
3514bee2e38SMatthew G. Knepley   }
352b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
353b18799e0SMatthew G. Knepley     PetscQuadrature   allPoints;
354b18799e0SMatthew G. Knepley     PetscInt          q, dim, numPoints;
355b18799e0SMatthew G. Knepley     const PetscReal   *points;
356b18799e0SMatthew G. Knepley     PetscScalar       *pointEval;
357b18799e0SMatthew G. Knepley     DM                dm;
358b18799e0SMatthew G. Knepley 
359b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
360b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
361b18799e0SMatthew G. Knepley     if (!funcs[f]) {
362b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
363b18799e0SMatthew G. Knepley       continue;
364b18799e0SMatthew G. Knepley     }
365b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
366b4457527SToby Isaac     ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
367b18799e0SMatthew G. Knepley     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
368b18799e0SMatthew G. Knepley     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
369b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
370b18799e0SMatthew G. Knepley       if (isAffine) {
371ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
372b18799e0SMatthew G. Knepley       } else {
3733fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp*dE];
3749f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp*dE*dE];
3759f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp*dE*dE];
3764bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
3774bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp*dE];
3789f209ee4SMatthew G. Knepley 
3799f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][tp*dE*dE];
3809f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][tp*dE*dE];
3819f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][tp];
382b18799e0SMatthew G. Knepley       }
383a6e0b375SMatthew 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 */
384ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
385ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
3864bee2e38SMatthew 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]);
387b18799e0SMatthew G. Knepley     }
388b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
389b18799e0SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
390b18799e0SMatthew G. Knepley     v += spDim;
391b18799e0SMatthew G. Knepley   }
392a6e0b375SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
393b18799e0SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
394b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
395b18799e0SMatthew G. Knepley }
396b18799e0SMatthew G. Knepley 
397a6e0b375SMatthew 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[],
398ef0bb6c7SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux,
3998c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
4008c6c5593SMatthew G. Knepley {
4018c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
4028c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
4038c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
4048c6c5593SMatthew G. Knepley 
4058c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
4068c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4078c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
4088c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
4098c6c5593SMatthew G. Knepley   switch (type) {
4108c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
4118c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
412a6e0b375SMatthew G. Knepley     ierr = DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
4138c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
4148c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
415ef0bb6c7SMatthew G. Knepley     ierr = DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
4160c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
4170c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
4180c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
419191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
420b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
421ef0bb6c7SMatthew G. Knepley     ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
422b18799e0SMatthew G. Knepley                                           (void (**)(PetscInt, PetscInt, PetscInt,
423b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
424b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
425b18799e0SMatthew G. Knepley                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
4268c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
4278c6c5593SMatthew G. Knepley   }
4288c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
4298c6c5593SMatthew G. Knepley }
4308c6c5593SMatthew G. Knepley 
431df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
432df5c1128SToby Isaac {
433df5c1128SToby Isaac   PetscReal      *points;
434df5c1128SToby Isaac   PetscInt       f, numPoints;
435df5c1128SToby Isaac   PetscErrorCode ierr;
436df5c1128SToby Isaac 
437df5c1128SToby Isaac   PetscFunctionBegin;
438df5c1128SToby Isaac   numPoints = 0;
439df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
440df5c1128SToby Isaac     if (funcs[f]) {
441df5c1128SToby Isaac       PetscQuadrature fAllPoints;
442df5c1128SToby Isaac       PetscInt        fNumPoints;
443df5c1128SToby Isaac 
444b4457527SToby Isaac       ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr);
445df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
446df5c1128SToby Isaac       numPoints += fNumPoints;
447df5c1128SToby Isaac     }
448df5c1128SToby Isaac   }
449df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
450df5c1128SToby Isaac   numPoints = 0;
451df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
452df5c1128SToby Isaac     if (funcs[f]) {
453df5c1128SToby Isaac       PetscQuadrature fAllPoints;
45454ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
455df5c1128SToby Isaac       const PetscReal *fPoints;
456df5c1128SToby Isaac 
457b4457527SToby Isaac       ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr);
45854ee0cceSMatthew G. Knepley       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
45954ee0cceSMatthew G. Knepley       if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %D for dual basis does not match input dimension %D", qdim, dim);
460df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
461df5c1128SToby Isaac       numPoints += fNumPoints;
462df5c1128SToby Isaac     }
463df5c1128SToby Isaac   }
464df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
465df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
466df5c1128SToby Isaac   PetscFunctionReturn(0);
467df5c1128SToby Isaac }
468df5c1128SToby Isaac 
469*e9f0ba4eSJed Brown /*
470*e9f0ba4eSJed Brown  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.
471*e9f0ba4eSJed Brown */
472a2c9b50fSJeremy L Thompson PetscErrorCode DMGetFirstLabelEntry_Internal(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *ds)
473e5e52638SMatthew G. Knepley {
474a6e0b375SMatthew G. Knepley   DM              plex;
475a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
476*e9f0ba4eSJed Brown   PetscInt        ls = -1;
477e5e52638SMatthew G. Knepley   PetscErrorCode  ierr;
478e5e52638SMatthew G. Knepley 
479e5e52638SMatthew G. Knepley   PetscFunctionBegin;
480e5e52638SMatthew G. Knepley   if (lStart) *lStart = -1;
481e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
482a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dm, odm, &enc);CHKERRQ(ierr);
483a6e0b375SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
484*e9f0ba4eSJed Brown   for (PetscInt i = 0; i < numIds; ++i) {
485*e9f0ba4eSJed Brown     IS       labelIS;
486*e9f0ba4eSJed Brown     PetscInt num_points, pStart, pEnd;
487*e9f0ba4eSJed Brown     ierr = DMLabelGetStratumIS(label, ids[i], &labelIS);CHKERRQ(ierr);
488*e9f0ba4eSJed Brown     if (!labelIS) continue; /* No points with that id on this process */
489*e9f0ba4eSJed Brown     ierr = DMPlexGetHeightStratum(plex, height, &pStart, &pEnd);CHKERRQ(ierr);
490*e9f0ba4eSJed Brown     ierr = ISGetSize(labelIS, &num_points);CHKERRQ(ierr);
491*e9f0ba4eSJed Brown     if (num_points) {
492e5e52638SMatthew G. Knepley       const PetscInt *points;
493*e9f0ba4eSJed Brown       ierr = ISGetIndices(labelIS, &points);CHKERRQ(ierr);
494*e9f0ba4eSJed Brown       for (PetscInt i=0; i<num_points; i++) {
495*e9f0ba4eSJed Brown         PetscInt point;
496*e9f0ba4eSJed Brown         ierr = DMGetEnclosurePoint(dm, odm, enc, points[i], &point);CHKERRQ(ierr);
497*e9f0ba4eSJed Brown         if (pStart <= point && point < pEnd) {
498a6e0b375SMatthew G. Knepley           ls = point;
499a6e0b375SMatthew G. Knepley           if (ds) {ierr = DMGetCellDS(dm, ls, ds);CHKERRQ(ierr);}
500e5e52638SMatthew G. Knepley         }
501*e9f0ba4eSJed Brown       }
502*e9f0ba4eSJed Brown       ierr = ISRestoreIndices(labelIS, &points);CHKERRQ(ierr);
503*e9f0ba4eSJed Brown     }
504*e9f0ba4eSJed Brown     ierr = ISDestroy(&labelIS);CHKERRQ(ierr);
505e5e52638SMatthew G. Knepley     if (ls >= 0) break;
506e5e52638SMatthew G. Knepley   }
507e5e52638SMatthew G. Knepley   if (lStart) *lStart = ls;
508a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
509e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
510e5e52638SMatthew G. Knepley }
511e5e52638SMatthew G. Knepley 
5120de51b56SMatthew G. Knepley /*
5130de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
5140de51b56SMatthew G. Knepley 
5150de51b56SMatthew G. Knepley   There are several different scenarios:
5160de51b56SMatthew G. Knepley 
5170de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
5180de51b56SMatthew G. Knepley 
5190de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
5200de51b56SMatthew G. Knepley 
5210de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
5220de51b56SMatthew G. Knepley 
5230de51b56SMatthew 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.
5240de51b56SMatthew G. Knepley 
5250de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
5260de51b56SMatthew G. Knepley 
5270de51b56SMatthew 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.
5280de51b56SMatthew G. Knepley 
529636b9322SMatthew G. Knepley   4) Volumetric input mesh with boundary output mesh
530636b9322SMatthew G. Knepley 
531636b9322SMatthew G. Knepley      Here we must get a subspace for the input DS
532636b9322SMatthew G. Knepley 
5330de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
5340de51b56SMatthew G. Knepley 
5350de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
5360de51b56SMatthew G. Knepley 
5370de51b56SMatthew 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.
5380de51b56SMatthew 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
5395f790a90SMatthew 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
5400de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
5410de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
5420de51b56SMatthew G. Knepley 
5430de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
5440de51b56SMatthew G. Knepley 
5450de51b56SMatthew 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.
5460de51b56SMatthew G. Knepley */
54746fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
5481c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
5498c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
5508c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
5518c6c5593SMatthew G. Knepley {
552a6e0b375SMatthew G. Knepley   DM                 plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
553a6e0b375SMatthew G. Knepley   DMEnclosureType    encIn, encAux;
554a6e0b375SMatthew G. Knepley   PetscDS            ds = NULL, dsIn = NULL, dsAux = NULL;
555ca3d3a14SMatthew G. Knepley   Vec                localA = NULL, tv;
556aa0eca99SMatthew G. Knepley   IS                 fieldIS;
55747923291SMatthew G. Knepley   PetscSection       section;
558636b9322SMatthew G. Knepley   PetscDualSpace    *sp, *cellsp, *spIn, *cellspIn;
559ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
5608c6c5593SMatthew G. Knepley   PetscInt          *Nc;
561636b9322SMatthew G. Knepley   PetscInt           dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
5625fedec97SMatthew G. Knepley   PetscBool         *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, transform;
563c330f8ffSToby Isaac   DMField            coordField;
564c330f8ffSToby Isaac   DMLabel            depthLabel;
565c330f8ffSToby Isaac   PetscQuadrature    allPoints = NULL;
56647923291SMatthew G. Knepley   PetscErrorCode     ierr;
56747923291SMatthew G. Knepley 
56847923291SMatthew G. Knepley   PetscFunctionBegin;
569a6e0b375SMatthew G. Knepley   if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);}
570a6e0b375SMatthew G. Knepley   else        {dmIn = dm;}
57104c51a94SMatthew G. Knepley   ierr = DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, &localA);CHKERRQ(ierr);
5729a2a23afSMatthew G. Knepley   if (localA) {ierr = VecGetDM(localA, &dmAux);CHKERRQ(ierr);} else {dmAux = NULL;}
573a6e0b375SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
574a6e0b375SMatthew G. Knepley   ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr);
575a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr);
576a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
5778c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
578a6e0b375SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr);
579ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
580ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
581ca3d3a14SMatthew G. Knepley   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
5820de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
583a6e0b375SMatthew G. Knepley   if (dmAux) {
584a6e0b375SMatthew G. Knepley     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
585a6e0b375SMatthew G. Knepley     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
586a6e0b375SMatthew G. Knepley       if (!localA) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector");
587636b9322SMatthew G. Knepley     }
588636b9322SMatthew G. Knepley   }
589636b9322SMatthew G. Knepley   /* Determine height for iteration of all meshes */
590636b9322SMatthew G. Knepley   {
591636b9322SMatthew G. Knepley     DMPolytopeType ct, ctIn, ctAux;
592636b9322SMatthew G. Knepley     PetscInt       minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux;
593636b9322SMatthew G. Knepley     PetscInt       dim = -1, dimIn, dimAux;
59488aca1feSMatthew G. Knepley 
595636b9322SMatthew G. Knepley     ierr = DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd);CHKERRQ(ierr);
596636b9322SMatthew G. Knepley     if (pEnd > pStart) {
597a2c9b50fSJeremy L Thompson       ierr = DMGetFirstLabelEntry_Internal(dm, dm, label, numIds, ids, minHeight, &lStart, NULL);CHKERRQ(ierr);
598636b9322SMatthew G. Knepley       p    = lStart < 0 ? pStart : lStart;
599636b9322SMatthew G. Knepley       ierr = DMPlexGetCellType(plex, p, &ct);CHKERRQ(ierr);
600636b9322SMatthew G. Knepley       dim  = DMPolytopeTypeGetDim(ct);
601636b9322SMatthew G. Knepley       ierr = DMPlexGetVTKCellHeight(plexIn, &minHeightIn);CHKERRQ(ierr);
602636b9322SMatthew G. Knepley       ierr = DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL);CHKERRQ(ierr);
603636b9322SMatthew G. Knepley       ierr = DMPlexGetCellType(plexIn, pStartIn, &ctIn);CHKERRQ(ierr);
604636b9322SMatthew G. Knepley       dimIn = DMPolytopeTypeGetDim(ctIn);
605636b9322SMatthew G. Knepley       if (dmAux) {
606636b9322SMatthew G. Knepley         ierr = DMPlexGetVTKCellHeight(plexAux, &minHeightAux);CHKERRQ(ierr);
607636b9322SMatthew G. Knepley         ierr = DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, NULL);CHKERRQ(ierr);
608636b9322SMatthew G. Knepley         ierr = DMPlexGetCellType(plexAux, pStartAux, &ctAux);CHKERRQ(ierr);
609636b9322SMatthew G. Knepley         dimAux = DMPolytopeTypeGetDim(ctAux);
610636b9322SMatthew G. Knepley       } else dimAux = dim;
611e39fcb4eSStefano Zampini     }
612636b9322SMatthew G. Knepley     if (dim < 0) {
613636b9322SMatthew G. Knepley       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
614636b9322SMatthew G. Knepley 
615636b9322SMatthew G. Knepley       /* Fall back to determination based on being a submesh */
616636b9322SMatthew G. Knepley       ierr = DMPlexGetSubpointMap(plex,   &spmap);CHKERRQ(ierr);
617636b9322SMatthew G. Knepley       ierr = DMPlexGetSubpointMap(plexIn, &spmapIn);CHKERRQ(ierr);
618636b9322SMatthew G. Knepley       if (plexAux) {ierr = DMPlexGetSubpointMap(plexAux, &spmapAux);CHKERRQ(ierr);}
619636b9322SMatthew G. Knepley       dim    = spmap    ? 1 : 0;
620636b9322SMatthew G. Knepley       dimIn  = spmapIn  ? 1 : 0;
621636b9322SMatthew G. Knepley       dimAux = spmapAux ? 1 : 0;
62288aca1feSMatthew G. Knepley     }
623636b9322SMatthew G. Knepley     {
624636b9322SMatthew G. Knepley       PetscInt dimProj = PetscMin(PetscMin(dim, dimIn), dimAux);
625636b9322SMatthew G. Knepley 
626636b9322SMatthew G. Knepley       if (PetscAbsInt(dimProj - dim) > 1 || PetscAbsInt(dimProj - dimIn) > 1 || PetscAbsInt(dimProj - dimAux) > 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not currently support differences of more than 1 in dimension");
627636b9322SMatthew G. Knepley       if (dimProj < dim) minHeight = 1;
628636b9322SMatthew G. Knepley       htInc    =  dim    - dimProj;
629636b9322SMatthew G. Knepley       htIncIn  =  dimIn  - dimProj;
630636b9322SMatthew G. Knepley       htIncAux =  dimAux - dimProj;
6310de51b56SMatthew G. Knepley     }
632a6e0b375SMatthew G. Knepley   }
633a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
634a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
635a6e0b375SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr);
6362af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
637e39fcb4eSStefano Zampini   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
638a2c9b50fSJeremy L Thompson   ierr = DMGetFirstLabelEntry_Internal(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr);
639a6e0b375SMatthew G. Knepley   if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);}
640a2c9b50fSJeremy L Thompson   ierr = DMGetFirstLabelEntry_Internal(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr);
641a6e0b375SMatthew G. Knepley   if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);}
642a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
643a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
644aa0eca99SMatthew G. Knepley   ierr = DMGetNumFields(dm, &NfTot);CHKERRQ(ierr);
645aa0eca99SMatthew G. Knepley   ierr = DMFindRegionNum(dm, ds, &regionNum);CHKERRQ(ierr);
646aa0eca99SMatthew G. Knepley   ierr = DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL);CHKERRQ(ierr);
6475fedec97SMatthew G. Knepley   ierr = PetscDSIsCohesive(ds, &isCohesive);CHKERRQ(ierr);
6488c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
64992fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
6500c898c18SMatthew G. Knepley   if (dmAux) {
651a6e0b375SMatthew G. Knepley     ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr);
652a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
6530c898c18SMatthew G. Knepley   }
654a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
655636b9322SMatthew G. Knepley   ierr = PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn);CHKERRQ(ierr);
656636b9322SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn);CHKERRQ(ierr);}
657636b9322SMatthew G. Knepley   else               {cellsp = sp; cellspIn = spIn;}
658a6e0b375SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
6598c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
66047923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
661665f567fSMatthew G. Knepley     PetscDiscType disctype;
66247923291SMatthew G. Knepley 
663665f567fSMatthew G. Knepley     ierr = PetscDSGetDiscType_Internal(ds, f, &disctype);CHKERRQ(ierr);
664665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
665665f567fSMatthew G. Knepley       PetscFE fe;
66647923291SMatthew G. Knepley 
66747923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
668665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
669665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
67047923291SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
671665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
672665f567fSMatthew G. Knepley       PetscFV fv;
67347923291SMatthew G. Knepley 
67447923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
675665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
676665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fv);CHKERRQ(ierr);
67747923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
678665f567fSMatthew G. Knepley     } else {
679665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
680665f567fSMatthew G. Knepley       cellsp[f] = NULL;
681665f567fSMatthew G. Knepley     }
68247923291SMatthew G. Knepley   }
683636b9322SMatthew G. Knepley   for (f = 0; f < NfIn; ++f) {
684636b9322SMatthew G. Knepley     PetscDiscType disctype;
685636b9322SMatthew G. Knepley 
686636b9322SMatthew G. Knepley     ierr = PetscDSGetDiscType_Internal(dsIn, f, &disctype);CHKERRQ(ierr);
687636b9322SMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
688636b9322SMatthew G. Knepley       PetscFE fe;
689636b9322SMatthew G. Knepley 
690636b9322SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fe);CHKERRQ(ierr);
691636b9322SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &cellspIn[f]);CHKERRQ(ierr);
692636b9322SMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
693636b9322SMatthew G. Knepley       PetscFV fv;
694636b9322SMatthew G. Knepley 
695636b9322SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fv);CHKERRQ(ierr);
696636b9322SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellspIn[f]);CHKERRQ(ierr);
697636b9322SMatthew G. Knepley     } else {
698636b9322SMatthew G. Knepley       cellspIn[f] = NULL;
699636b9322SMatthew G. Knepley     }
700636b9322SMatthew G. Knepley   }
701c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
702636b9322SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
703636b9322SMatthew G. Knepley     if (!htInc) {sp[f] = cellsp[f];}
704636b9322SMatthew G. Knepley     else        {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]);CHKERRQ(ierr);}
705636b9322SMatthew G. Knepley   }
706ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
70774fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
708665f567fSMatthew G. Knepley     PetscDiscType    disctype;
7094a3e9fdbSToby Isaac     const PetscReal *points;
7104a3e9fdbSToby Isaac     PetscInt         numPoints;
7110c898c18SMatthew G. Knepley 
7122af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
713636b9322SMatthew G. Knepley     ierr = PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &allPoints);CHKERRQ(ierr);
714c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL);CHKERRQ(ierr);
715ef0bb6c7SMatthew G. Knepley     ierr = PetscMalloc2(NfIn, &T, NfAux, &TAux);CHKERRQ(ierr);
716a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
717636b9322SMatthew G. Knepley       if (!htIncIn) {spIn[f] = cellspIn[f];}
718636b9322SMatthew G. Knepley       else          {ierr = PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]);CHKERRQ(ierr);}
719636b9322SMatthew G. Knepley 
720665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscType_Internal(dsIn, f, &disctype);CHKERRQ(ierr);
721665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
722a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
723636b9322SMatthew G. Knepley       if (!htIncIn) {subfem = fem;}
724636b9322SMatthew G. Knepley       else          {ierr = PetscFEGetHeightSubspace(fem, htIncIn, &subfem);CHKERRQ(ierr);}
725ef0bb6c7SMatthew G. Knepley       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]);CHKERRQ(ierr);
7260c898c18SMatthew G. Knepley     }
7270c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
728665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscType_Internal(dsAux, f, &disctype);CHKERRQ(ierr);
729665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
730a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
731636b9322SMatthew G. Knepley       if (!htIncAux) {subfem = fem;}
732636b9322SMatthew G. Knepley       else           {ierr = PetscFEGetHeightSubspace(fem, htIncAux, &subfem);CHKERRQ(ierr);}
733ef0bb6c7SMatthew G. Knepley       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]);CHKERRQ(ierr);
7340c898c18SMatthew G. Knepley     }
7350c898c18SMatthew G. Knepley   }
73647923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
7372af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
738636b9322SMatthew G. Knepley     PetscInt     hEff     = h - minHeight + htInc;
739636b9322SMatthew G. Knepley     PetscInt     hEffIn   = h - minHeight + htIncIn;
740636b9322SMatthew G. Knepley     PetscInt     hEffAux  = h - minHeight + htIncAux;
741a6e0b375SMatthew G. Knepley     PetscDS      dsEff    = ds;
742636b9322SMatthew G. Knepley     PetscDS      dsEffIn  = dsIn;
743636b9322SMatthew G. Knepley     PetscDS      dsEffAux = dsAux;
7448c6c5593SMatthew G. Knepley     PetscScalar *values;
745b7260050SToby Isaac     PetscBool   *fieldActive;
746b7260050SToby Isaac     PetscInt     maxDegree;
747e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
748c330f8ffSToby Isaac     IS           heightIS;
7498c6c5593SMatthew G. Knepley 
750636b9322SMatthew G. Knepley     if (h > minHeight) {
751636b9322SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]);CHKERRQ(ierr);}
752636b9322SMatthew G. Knepley     }
753412e9a14SMatthew G. Knepley     ierr = DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd);CHKERRQ(ierr);
754a2c9b50fSJeremy L Thompson     ierr = DMGetFirstLabelEntry_Internal(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
755c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
756c330f8ffSToby Isaac     if (pEnd <= pStart) {
757c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
758c330f8ffSToby Isaac       continue;
759c330f8ffSToby Isaac     }
76047923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
76147923291SMatthew G. Knepley     totDim = 0;
76247923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
7635fedec97SMatthew G. Knepley       PetscBool cohesive;
7645fedec97SMatthew G. Knepley 
765665f567fSMatthew G. Knepley       if (!sp[f]) continue;
7665fedec97SMatthew G. Knepley       ierr = PetscDSGetCohesive(ds, f, &cohesive);CHKERRQ(ierr);
76747923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
7689c3cf19fSMatthew G. Knepley       totDim += spDim;
7695fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) totDim += spDim;
77047923291SMatthew G. Knepley     }
771636b9322SMatthew G. Knepley     p    = lStart < 0 ? pStart : lStart;
772636b9322SMatthew G. Knepley     ierr = DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL);CHKERRQ(ierr);
773636b9322SMatthew G. Knepley     if (numValues != totDim) SETERRQ6(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);
774c330f8ffSToby Isaac     if (!totDim) {
775c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
776c330f8ffSToby Isaac       continue;
777c330f8ffSToby Isaac     }
778636b9322SMatthew G. Knepley     if (htInc) {ierr = PetscDSGetHeightSubspace(ds, hEff, &dsEff);CHKERRQ(ierr);}
779636b9322SMatthew G. Knepley     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
780636b9322SMatthew G. Knepley     if (localU) {
781636b9322SMatthew G. Knepley       PetscInt totDimIn, pIn, numValuesIn;
782636b9322SMatthew G. Knepley 
783636b9322SMatthew G. Knepley       totDimIn = 0;
784636b9322SMatthew G. Knepley       for (f = 0; f < NfIn; ++f) {
7855fedec97SMatthew G. Knepley         PetscBool cohesive;
7865fedec97SMatthew G. Knepley 
787636b9322SMatthew G. Knepley         if (!spIn[f]) continue;
7885fedec97SMatthew G. Knepley         ierr = PetscDSGetCohesive(dsIn, f, &cohesive);CHKERRQ(ierr);
789636b9322SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(spIn[f], &spDim);CHKERRQ(ierr);
790636b9322SMatthew G. Knepley         totDimIn += spDim;
7915fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) totDimIn += spDim;
792636b9322SMatthew G. Knepley       }
793636b9322SMatthew G. Knepley       ierr = DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn);CHKERRQ(ierr);
794636b9322SMatthew G. Knepley       ierr = DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL);CHKERRQ(ierr);
795636b9322SMatthew G. Knepley       if (numValuesIn != totDimIn) SETERRQ4(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);
796636b9322SMatthew G. Knepley       if (htIncIn) {ierr = PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn);CHKERRQ(ierr);}
797636b9322SMatthew G. Knepley     }
798636b9322SMatthew G. Knepley     if (htIncAux) {ierr = PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux);CHKERRQ(ierr);}
79947923291SMatthew G. Knepley     /* Loop over points at this height */
80069291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
801aa0eca99SMatthew G. Knepley     ierr = DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive);CHKERRQ(ierr);
802aa0eca99SMatthew G. Knepley     {
803aa0eca99SMatthew G. Knepley       const PetscInt *fields;
804aa0eca99SMatthew G. Knepley 
805aa0eca99SMatthew G. Knepley       ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);
806aa0eca99SMatthew G. Knepley       for (f = 0; f < NfTot; ++f) {fieldActive[f] = PETSC_FALSE;}
807aa0eca99SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;}
808aa0eca99SMatthew G. Knepley       ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);
809aa0eca99SMatthew G. Knepley     }
8108c6c5593SMatthew G. Knepley     if (label) {
8118c6c5593SMatthew G. Knepley       PetscInt i;
81247923291SMatthew G. Knepley 
81347923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
814c330f8ffSToby Isaac         IS              pointIS, isectIS;
81547923291SMatthew G. Knepley         const PetscInt *points;
8168c6c5593SMatthew G. Knepley         PetscInt        n;
817c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
818c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
81947923291SMatthew G. Knepley 
82047923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
82147923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
822c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
823c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
824c330f8ffSToby Isaac         if (!isectIS) continue;
825c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
826c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
827b7260050SToby Isaac         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
828b7260050SToby Isaac         if (maxDegree <= 1) {
829c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
830c330f8ffSToby Isaac         }
831c330f8ffSToby Isaac         if (!quad) {
832c330f8ffSToby Isaac           if (!h && allPoints) {
833c330f8ffSToby Isaac             quad = allPoints;
834c330f8ffSToby Isaac             allPoints = NULL;
835c330f8ffSToby Isaac           } else {
8365fedec97SMatthew G. Knepley             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,isCohesive ? dim-htInc-1 : dim-htInc,funcs,&quad);CHKERRQ(ierr);
837c330f8ffSToby Isaac           }
838c330f8ffSToby Isaac         }
839636b9322SMatthew G. Knepley         ierr = DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom);CHKERRQ(ierr);
84047923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
84147923291SMatthew G. Knepley           const PetscInt  point = points[p];
84247923291SMatthew G. Knepley 
843580bdb30SBarry Smith           ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
844c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
8451b32699bSMatthew G. Knepley           ierr = DMPlexSetActivePoint(dm, point);CHKERRQ(ierr);
846636b9322SMatthew 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);
84747923291SMatthew G. Knepley           if (ierr) {
84847923291SMatthew G. Knepley             PetscErrorCode ierr2;
84969291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
85069291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
85147923291SMatthew G. Knepley             CHKERRQ(ierr);
85247923291SMatthew G. Knepley           }
853a6e0b375SMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
8545f790a90SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode);CHKERRQ(ierr);
85547923291SMatthew G. Knepley         }
856c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
857c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
858c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
859c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
860c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
86147923291SMatthew G. Knepley       }
8628c6c5593SMatthew G. Knepley     } else {
863c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
864c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
865c330f8ffSToby Isaac       IS              pointIS;
866c330f8ffSToby Isaac 
867c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
868b7260050SToby Isaac       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
869b7260050SToby Isaac       if (maxDegree <= 1) {
870c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
871c330f8ffSToby Isaac       }
872c330f8ffSToby Isaac       if (!quad) {
873c330f8ffSToby Isaac         if (!h && allPoints) {
874c330f8ffSToby Isaac           quad = allPoints;
875c330f8ffSToby Isaac           allPoints = NULL;
876c330f8ffSToby Isaac         } else {
877636b9322SMatthew G. Knepley           ierr = PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &quad);CHKERRQ(ierr);
878c330f8ffSToby Isaac         }
879c330f8ffSToby Isaac       }
880636b9322SMatthew G. Knepley       ierr = DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom);CHKERRQ(ierr);
8818c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
882580bdb30SBarry Smith         ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
883c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
8841b32699bSMatthew G. Knepley         ierr = DMPlexSetActivePoint(dm, p);CHKERRQ(ierr);
885636b9322SMatthew 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);
8868c6c5593SMatthew G. Knepley         if (ierr) {
8878c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
88869291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
88969291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
8908c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
8918c6c5593SMatthew G. Knepley         }
892a6e0b375SMatthew G. Knepley         if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
8935f790a90SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode);CHKERRQ(ierr);
8948c6c5593SMatthew G. Knepley       }
895c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
896c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
897c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
898c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
8998c6c5593SMatthew G. Knepley     }
900c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
90169291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
90269291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
90347923291SMatthew G. Knepley   }
9048c6c5593SMatthew G. Knepley   /* Cleanup */
905ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
906636b9322SMatthew G. Knepley     for (f = 0; f < NfIn;  ++f) {ierr = PetscTabulationDestroy(&T[f]);CHKERRQ(ierr);}
907636b9322SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {ierr = PetscTabulationDestroy(&TAux[f]);CHKERRQ(ierr);}
908ef0bb6c7SMatthew G. Knepley     ierr = PetscFree2(T, TAux);CHKERRQ(ierr);
9090c898c18SMatthew G. Knepley   }
910c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
911636b9322SMatthew G. Knepley   ierr = PetscFree3(isFE, sp, spIn);CHKERRQ(ierr);
912636b9322SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree2(cellsp, cellspIn);CHKERRQ(ierr);}
913a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
914a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plexIn);CHKERRQ(ierr);
915a6e0b375SMatthew G. Knepley   if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);}
9168c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
91747923291SMatthew G. Knepley }
9188c6c5593SMatthew G. Knepley 
9198c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
9208c6c5593SMatthew G. Knepley {
9218c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
9228c6c5593SMatthew G. Knepley 
9238c6c5593SMatthew G. Knepley   PetscFunctionBegin;
924636b9322SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
9258c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
9268c6c5593SMatthew G. Knepley }
9278c6c5593SMatthew G. Knepley 
9281c531cf8SMatthew 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)
9298c6c5593SMatthew G. Knepley {
9308c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
9318c6c5593SMatthew G. Knepley 
9328c6c5593SMatthew G. Knepley   PetscFunctionBegin;
933636b9322SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
93447923291SMatthew G. Knepley   PetscFunctionReturn(0);
93547923291SMatthew G. Knepley }
93647923291SMatthew G. Knepley 
9378c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
93847923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
93947923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
94047923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
941191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
94247923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
94347923291SMatthew G. Knepley {
94447923291SMatthew G. Knepley   PetscErrorCode ierr;
94547923291SMatthew G. Knepley 
94647923291SMatthew G. Knepley   PetscFunctionBegin;
9471c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
9488c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
94947923291SMatthew G. Knepley }
95047923291SMatthew G. Knepley 
9511c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
9528c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
9538c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
9548c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
955191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
9568c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
9578c6c5593SMatthew G. Knepley {
9588c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
95947923291SMatthew G. Knepley 
9608c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9611c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
96247923291SMatthew G. Knepley   PetscFunctionReturn(0);
96347923291SMatthew G. Knepley }
964ece3a9fcSMatthew G. Knepley 
965ece3a9fcSMatthew G. Knepley PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
966ece3a9fcSMatthew G. Knepley                                                void (**funcs)(PetscInt, PetscInt, PetscInt,
967ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
968ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
969ece3a9fcSMatthew G. Knepley                                                               PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
970ece3a9fcSMatthew G. Knepley                                                InsertMode mode, Vec localX)
971ece3a9fcSMatthew G. Knepley {
972ece3a9fcSMatthew G. Knepley   PetscErrorCode ierr;
973ece3a9fcSMatthew G. Knepley 
974ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
975ece3a9fcSMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
976ece3a9fcSMatthew G. Knepley   PetscFunctionReturn(0);
977ece3a9fcSMatthew G. Knepley }
978