xref: /petsc/src/dm/impls/plex/plexproject.c (revision 4bee2e389ac4efdf19d1420f70098a911b40ccd1)
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 
5c330f8ffSToby Isaac static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS prob, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
68c6c5593SMatthew G. Knepley                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
78c6c5593SMatthew G. Knepley                                                   PetscScalar values[])
847923291SMatthew G. Knepley {
9c330f8ffSToby Isaac   PetscInt       coordDim, Nf, *Nc, f, totDim, spDim, d, v, tp;
10ca3d3a14SMatthew G. Knepley   PetscBool      isAffine, transform;
118c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
128c6c5593SMatthew G. Knepley 
138c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
14c330f8ffSToby Isaac   ierr = DMGetCoordinateDim(dm,&coordDim);CHKERRQ(ierr);
15ca3d3a14SMatthew G. Knepley   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
168c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
17496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
189c3cf19fSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
198c6c5593SMatthew G. Knepley   /* Get values for closure */
20c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
21c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
228c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
238c6c5593SMatthew G. Knepley 
248c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
258c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
268c6c5593SMatthew G. Knepley     if (funcs[f]) {
27c330f8ffSToby Isaac       if (isFE[f]) {
28c330f8ffSToby Isaac         PetscQuadrature   allPoints;
29c330f8ffSToby Isaac         PetscInt          q, dim, numPoints;
30c330f8ffSToby Isaac         const PetscReal   *points;
31c330f8ffSToby Isaac         PetscScalar       *pointEval;
32c330f8ffSToby Isaac         PetscReal         *x;
33ca3d3a14SMatthew G. Knepley         DM                rdm;
34c330f8ffSToby Isaac 
35ca3d3a14SMatthew G. Knepley         ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr);
36c330f8ffSToby Isaac         ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
37c330f8ffSToby Isaac         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
38ca3d3a14SMatthew G. Knepley         ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
39ca3d3a14SMatthew G. Knepley         ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
40c330f8ffSToby Isaac         for (q = 0; q < numPoints; q++, tp++) {
41c330f8ffSToby Isaac           const PetscReal *v0;
42c330f8ffSToby Isaac 
43c330f8ffSToby Isaac           if (isAffine) {
44c330f8ffSToby Isaac             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
45c330f8ffSToby Isaac             v0 = x;
468c6c5593SMatthew G. Knepley           } else {
47c330f8ffSToby Isaac             v0 = &fegeom->v[tp*coordDim];
488c6c5593SMatthew G. Knepley           }
49ec277c0fSMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;}
50c330f8ffSToby Isaac           ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
51c330f8ffSToby Isaac         }
52*4bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
53*4bee2e38SMatthew G. Knepley         ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr);
54c330f8ffSToby Isaac         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
55ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
56ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
57c330f8ffSToby Isaac         v += spDim;
58c330f8ffSToby Isaac       } else {
59c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
60c330f8ffSToby Isaac           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
61c330f8ffSToby Isaac         }
62c330f8ffSToby Isaac       }
63c330f8ffSToby Isaac     } else {
64c330f8ffSToby Isaac       for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;}
658c6c5593SMatthew G. Knepley     }
669c3cf19fSMatthew G. Knepley   }
678c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
688c6c5593SMatthew G. Knepley }
698c6c5593SMatthew G. Knepley 
70*4bee2e38SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
710c898c18SMatthew G. Knepley                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
728c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
738c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
748c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
75191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
768c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
778c6c5593SMatthew G. Knepley {
788c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
79*4bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
808c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
818c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
82191494d9SMatthew G. Knepley   const PetscScalar *constants;
838c6c5593SMatthew G. Knepley   PetscReal         *x;
84496733e2SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
85*4bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
86*4bee2e38SMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed;
87c330f8ffSToby Isaac   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
88ca3d3a14SMatthew G. Knepley   PetscBool          isAffine, transform;
898c6c5593SMatthew G. Knepley   PetscErrorCode     ierr;
908c6c5593SMatthew G. Knepley 
918c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
928c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
93496733e2SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
94496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
958c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
968c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
979c3cf19fSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
98*4bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
99191494d9SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
100ca3d3a14SMatthew G. Knepley   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
101e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
102496733e2SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1038c6c5593SMatthew G. Knepley   if (dmAux) {
10444171101SMatthew G. Knepley     PetscInt subp;
1051ba34526SMatthew G. Knepley 
10644171101SMatthew G. Knepley     ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, p, &subp);CHKERRQ(ierr);
1071ba34526SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
1088c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
109496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
110496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
111e87a4003SBarry Smith     ierr = DMGetSection(dmAux, &sectionAux);CHKERRQ(ierr);
1128c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
1138c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
11411d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
1151ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
1168c6c5593SMatthew G. Knepley   }
1178c6c5593SMatthew G. Knepley   /* Get values for closure */
118*4bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
119*4bee2e38SMatthew G. Knepley   if (isAffine) {
120*4bee2e38SMatthew G. Knepley     fegeom.v    = x;
121*4bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
122*4bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
123*4bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
124*4bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
125*4bee2e38SMatthew G. Knepley   }
1268c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
127c330f8ffSToby Isaac     PetscQuadrature   allPoints;
128c330f8ffSToby Isaac     PetscInt          q, dim, numPoints;
129c330f8ffSToby Isaac     const PetscReal   *points;
130c330f8ffSToby Isaac     PetscScalar       *pointEval;
131c330f8ffSToby Isaac     DM                dm;
132c330f8ffSToby Isaac 
1338c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
1348c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
135c330f8ffSToby Isaac     if (!funcs[f]) {
136be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
137c330f8ffSToby Isaac       continue;
138c330f8ffSToby Isaac     }
139c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
140c330f8ffSToby Isaac     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
141c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
142c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
1430c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
144c330f8ffSToby Isaac       if (isAffine) {
145*4bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
1461c531cf8SMatthew G. Knepley       } else {
147*4bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp*dE];
148*4bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp*dE*dE];
149*4bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
150*4bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
1518c6c5593SMatthew G. Knepley       }
152*4bee2e38SMatthew G. Knepley       ierr = EvaluateFieldJets(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
153*4bee2e38SMatthew G. Knepley       if (probAux) {ierr = EvaluateFieldJets(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
154*4bee2e38SMatthew G. Knepley       if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dm, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);}
155*4bee2e38SMatthew 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, numConstants, constants, &pointEval[Nc[f]*q]);
1561c531cf8SMatthew G. Knepley     }
157c330f8ffSToby Isaac     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
158c330f8ffSToby Isaac     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
159c330f8ffSToby Isaac     v += spDim;
1608c6c5593SMatthew G. Knepley   }
1618c6c5593SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1628c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
1638c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1648c6c5593SMatthew G. Knepley }
1658c6c5593SMatthew G. Knepley 
166*4bee2e38SMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
167b18799e0SMatthew G. Knepley                                                      PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
168b18799e0SMatthew G. Knepley                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
169b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
170b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
171b18799e0SMatthew G. Knepley                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
172b18799e0SMatthew G. Knepley                                                      PetscScalar values[])
173b18799e0SMatthew G. Knepley {
174b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
175*4bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
176b18799e0SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
177b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
178b18799e0SMatthew G. Knepley   const PetscScalar *constants;
179b18799e0SMatthew G. Knepley   PetscReal         *x;
180b18799e0SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
181*4bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
182*4bee2e38SMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed;
183b18799e0SMatthew G. Knepley   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
184b18799e0SMatthew G. Knepley   PetscBool          isAffine;
185b18799e0SMatthew G. Knepley   PetscErrorCode     ierr;
186b18799e0SMatthew G. Knepley 
187b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
188b18799e0SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
189b18799e0SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
190b18799e0SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
191b18799e0SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
192b18799e0SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
193b18799e0SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
194*4bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(prob, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
195b18799e0SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
196b18799e0SMatthew G. Knepley   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
197b18799e0SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
198b18799e0SMatthew G. Knepley   if (dmAux) {
199b18799e0SMatthew G. Knepley     DMLabel  spmap;
200b18799e0SMatthew G. Knepley     PetscInt subp = p;
201b18799e0SMatthew G. Knepley 
202b18799e0SMatthew G. Knepley     /* If dm is a submesh, do not get subpoint */
203b18799e0SMatthew G. Knepley     ierr = DMPlexGetSubpointMap(dm, &spmap);CHKERRQ(ierr);
204b18799e0SMatthew G. Knepley     if (!spmap) {ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr);}
205b18799e0SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
206b18799e0SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
207b18799e0SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
208b18799e0SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
209b18799e0SMatthew G. Knepley     ierr = DMGetSection(dmAux, &sectionAux);CHKERRQ(ierr);
210b18799e0SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
211b18799e0SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
212b18799e0SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
213b18799e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
214b18799e0SMatthew G. Knepley   }
215b18799e0SMatthew G. Knepley   /* Get values for closure */
216*4bee2e38SMatthew G. Knepley   isAffine = fgeom->isAffine;
217*4bee2e38SMatthew G. Knepley   if (isAffine) {
218*4bee2e38SMatthew G. Knepley     fegeom.v    = x;
219*4bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
220*4bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
221*4bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
222*4bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
223*4bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
224*4bee2e38SMatthew G. Knepley   }
225b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
226b18799e0SMatthew G. Knepley     PetscQuadrature   allPoints;
227b18799e0SMatthew G. Knepley     PetscInt          q, dim, numPoints;
228b18799e0SMatthew G. Knepley     const PetscReal   *points;
229b18799e0SMatthew G. Knepley     PetscScalar       *pointEval;
230b18799e0SMatthew G. Knepley     DM                dm;
231b18799e0SMatthew G. Knepley 
232b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
233b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
234b18799e0SMatthew G. Knepley     if (!funcs[f]) {
235b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
236b18799e0SMatthew G. Knepley       continue;
237b18799e0SMatthew G. Knepley     }
238b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
239b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
240b18799e0SMatthew G. Knepley     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
241b18799e0SMatthew G. Knepley     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
242b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
243b18799e0SMatthew G. Knepley       if (isAffine) {
244*4bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
245b18799e0SMatthew G. Knepley       } else {
246*4bee2e38SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp*dE];
247*4bee2e38SMatthew G. Knepley         fegeom.invJ = &fgeom->suppInvJ[0][tp*dE*dE];
248*4bee2e38SMatthew G. Knepley         //fegeom.invJ = &fgeom->invJ[tp*dE*dE];
249*4bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
250*4bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp*dE];
251b18799e0SMatthew G. Knepley       }
252*4bee2e38SMatthew G. Knepley       ierr = EvaluateFieldJets(prob, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
253*4bee2e38SMatthew G. Knepley       if (probAux) {ierr = EvaluateFieldJets(probAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
254*4bee2e38SMatthew 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]);
255b18799e0SMatthew G. Knepley     }
256b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
257b18799e0SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
258b18799e0SMatthew G. Knepley     v += spDim;
259b18799e0SMatthew G. Knepley   }
260b18799e0SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
261b18799e0SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
262b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
263b18799e0SMatthew G. Knepley }
264b18799e0SMatthew G. Knepley 
265c330f8ffSToby Isaac static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS prob, PetscFEGeom *fegeom, DM dmAux, PetscDS probAux, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
2661c531cf8SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
2671c531cf8SMatthew G. Knepley                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
2688c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
2698c6c5593SMatthew G. Knepley {
2708c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
2718c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
2728c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
2738c6c5593SMatthew G. Knepley 
2748c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
2758c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2768c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
2778c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
2788c6c5593SMatthew G. Knepley   switch (type) {
2798c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
2808c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
281c330f8ffSToby Isaac     ierr = DMProjectPoint_Func_Private(dm, prob, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
2828c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
2838c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
284c330f8ffSToby Isaac     ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
2850c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
2860c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
2870c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
2880c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
289191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
290b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
291b18799e0SMatthew G. Knepley     ierr = DMProjectPoint_BdField_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
292b18799e0SMatthew G. Knepley                                           basisTab, basisDerTab, basisTabAux, basisDerTabAux,
293b18799e0SMatthew G. Knepley                                           (void (**)(PetscInt, PetscInt, PetscInt,
294b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
295b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
296b18799e0SMatthew G. Knepley                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
2978c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
2988c6c5593SMatthew G. Knepley   }
2998c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
3008c6c5593SMatthew G. Knepley }
3018c6c5593SMatthew G. Knepley 
302df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
303df5c1128SToby Isaac {
304df5c1128SToby Isaac   PetscReal      *points;
305df5c1128SToby Isaac   PetscInt       f, numPoints;
306df5c1128SToby Isaac   PetscErrorCode ierr;
307df5c1128SToby Isaac 
308df5c1128SToby Isaac   PetscFunctionBegin;
309df5c1128SToby Isaac   numPoints = 0;
310df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
311df5c1128SToby Isaac     if (funcs[f]) {
312df5c1128SToby Isaac       PetscQuadrature fAllPoints;
313df5c1128SToby Isaac       PetscInt        fNumPoints;
314df5c1128SToby Isaac 
315df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
316df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
317df5c1128SToby Isaac       numPoints += fNumPoints;
318df5c1128SToby Isaac     }
319df5c1128SToby Isaac   }
320df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
321df5c1128SToby Isaac   numPoints = 0;
322df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
323df5c1128SToby Isaac     PetscInt spDim;
324df5c1128SToby Isaac 
325df5c1128SToby Isaac     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
326df5c1128SToby Isaac     if (funcs[f]) {
327df5c1128SToby Isaac       PetscQuadrature fAllPoints;
32854ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
329df5c1128SToby Isaac       const PetscReal *fPoints;
330df5c1128SToby Isaac 
331df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
33254ee0cceSMatthew G. Knepley       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
33354ee0cceSMatthew 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);
334df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
335df5c1128SToby Isaac       numPoints += fNumPoints;
336df5c1128SToby Isaac     }
337df5c1128SToby Isaac   }
338df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
339df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
340df5c1128SToby Isaac   PetscFunctionReturn(0);
341df5c1128SToby Isaac }
342df5c1128SToby Isaac 
343e5e52638SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *prob)
344e5e52638SMatthew G. Knepley {
345e5e52638SMatthew G. Knepley   DMLabel        depthLabel;
346e5e52638SMatthew G. Knepley   PetscInt       dim, cdepth, ls = -1, i;
347e5e52638SMatthew G. Knepley   PetscErrorCode ierr;
348e5e52638SMatthew G. Knepley 
349e5e52638SMatthew G. Knepley   PetscFunctionBegin;
350e5e52638SMatthew G. Knepley   if (lStart) *lStart = -1;
351e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
352e5e52638SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
353e5e52638SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
354e5e52638SMatthew G. Knepley   cdepth = dim - height;
355e5e52638SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
356e5e52638SMatthew G. Knepley     IS              pointIS;
357e5e52638SMatthew G. Knepley     const PetscInt *points;
358e5e52638SMatthew G. Knepley     PetscInt        pdepth;
359e5e52638SMatthew G. Knepley 
360e5e52638SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
361e5e52638SMatthew G. Knepley     if (!pointIS) continue; /* No points with that id on this process */
362e5e52638SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
363e5e52638SMatthew G. Knepley     ierr = DMLabelGetValue(depthLabel, points[0], &pdepth);CHKERRQ(ierr);
364e5e52638SMatthew G. Knepley     if (pdepth == cdepth) {
365e5e52638SMatthew G. Knepley       ls = points[0];
366e5e52638SMatthew G. Knepley       if (prob) {ierr = DMGetCellDS(dm, ls, prob);CHKERRQ(ierr);}
367e5e52638SMatthew G. Knepley     }
368e5e52638SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
369e5e52638SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
370e5e52638SMatthew G. Knepley     if (ls >= 0) break;
371e5e52638SMatthew G. Knepley   }
372e5e52638SMatthew G. Knepley   if (lStart) *lStart = ls;
373e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
374e5e52638SMatthew G. Knepley }
375e5e52638SMatthew G. Knepley 
3760de51b56SMatthew G. Knepley /*
3770de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
3780de51b56SMatthew G. Knepley 
3790de51b56SMatthew G. Knepley   There are several different scenarios:
3800de51b56SMatthew G. Knepley 
3810de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
3820de51b56SMatthew G. Knepley 
3830de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
3840de51b56SMatthew G. Knepley 
3850de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
3860de51b56SMatthew G. Knepley 
3870de51b56SMatthew 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.
3880de51b56SMatthew G. Knepley 
3890de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
3900de51b56SMatthew G. Knepley 
3910de51b56SMatthew 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.
3920de51b56SMatthew G. Knepley 
3930de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
3940de51b56SMatthew G. Knepley 
3950de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
3960de51b56SMatthew G. Knepley 
3970de51b56SMatthew 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.
3980de51b56SMatthew 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
3990de51b56SMatthew G. Knepley       is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When th effective height is nonzero, we need to extract
4000de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
4010de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
4020de51b56SMatthew G. Knepley 
4030de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
4040de51b56SMatthew G. Knepley 
4050de51b56SMatthew 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.
4060de51b56SMatthew G. Knepley */
40746fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
4081c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
4098c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
4108c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
4118c6c5593SMatthew G. Knepley {
412ca3d3a14SMatthew G. Knepley   DM              dmAux = NULL, tdm;
413e5e52638SMatthew G. Knepley   PetscDS         prob = NULL, probAux = NULL;
414ca3d3a14SMatthew G. Knepley   Vec             localA = NULL, tv;
41547923291SMatthew G. Knepley   PetscSection    section;
4168c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
4170c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
4188c6c5593SMatthew G. Knepley   PetscInt       *Nc;
419c330f8ffSToby Isaac   PetscInt        dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f;
420ca3d3a14SMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, transform;
421c330f8ffSToby Isaac   DMField         coordField;
422c330f8ffSToby Isaac   DMLabel         depthLabel;
423c330f8ffSToby Isaac   PetscQuadrature allPoints = NULL;
42447923291SMatthew G. Knepley   PetscErrorCode  ierr;
42547923291SMatthew G. Knepley 
42647923291SMatthew G. Knepley   PetscFunctionBegin;
42788aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
42888aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
4298c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4302af58022SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
431ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
432ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
433ca3d3a14SMatthew G. Knepley   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
4340de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
4350de51b56SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
436e39fcb4eSStefano Zampini     if (dmAux && !localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
43788aca1feSMatthew G. Knepley     if (!minHeight && dmAux) {
43888aca1feSMatthew G. Knepley       DMLabel spmap;
43988aca1feSMatthew G. Knepley 
44088aca1feSMatthew G. Knepley       /* If dmAux is a surface, then force the projection to take place over a surface */
44188aca1feSMatthew G. Knepley       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
442e39fcb4eSStefano Zampini       if (spmap) {
443e39fcb4eSStefano Zampini         ierr = DMPlexGetVTKCellHeight(dmAux, &minHeight);CHKERRQ(ierr);
444e39fcb4eSStefano Zampini         auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
445e39fcb4eSStefano Zampini       }
44688aca1feSMatthew G. Knepley     }
4470de51b56SMatthew G. Knepley   }
448e5e52638SMatthew G. Knepley   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
449e5e52638SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
4508c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
4512af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
452e39fcb4eSStefano Zampini   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
453e5e52638SMatthew G. Knepley   ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, 0, NULL, &prob);CHKERRQ(ierr);
454e5e52638SMatthew G. Knepley   if (!prob) {ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);}
455e5e52638SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
4568c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
457e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
4580c898c18SMatthew G. Knepley   if (dmAux) {
4590c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
4600c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
4610c898c18SMatthew G. Knepley   }
462496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
463496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
4648c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
4658c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
4668c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
4678c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
46847923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
46947923291SMatthew G. Knepley     PetscObject  obj;
47047923291SMatthew G. Knepley     PetscClassId id;
47147923291SMatthew G. Knepley 
472e5e52638SMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
47347923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
47447923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
47547923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
47647923291SMatthew G. Knepley 
47747923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
47847923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
47947923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
48047923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
48147923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
48247923291SMatthew G. Knepley 
48347923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
48447923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
48547923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
48647923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
48747923291SMatthew G. Knepley   }
488c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
4890c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
49074fc349aSMatthew G. Knepley     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
49174fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
4924a3e9fdbSToby Isaac     const PetscReal *points;
4934a3e9fdbSToby Isaac     PetscInt         numPoints;
4940c898c18SMatthew G. Knepley 
4952af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
49654ee0cceSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
49754ee0cceSMatthew G. Knepley       if (!effectiveHeight) {sp[f] = cellsp[f];}
49854ee0cceSMatthew G. Knepley       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
49954ee0cceSMatthew G. Knepley     }
50054ee0cceSMatthew G. Knepley     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
501c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
5020c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
5030c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
5040c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
5050c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
50674fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
50774fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
50874fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
5090c898c18SMatthew G. Knepley     }
5100c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
5110c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
51274fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
51374fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
51474fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
5150c898c18SMatthew G. Knepley     }
5160c898c18SMatthew G. Knepley   }
51747923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
5182af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
51988aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
5200de51b56SMatthew G. Knepley     PetscDS      probEff         = prob;
5218c6c5593SMatthew G. Knepley     PetscScalar *values;
522b7260050SToby Isaac     PetscBool   *fieldActive;
523b7260050SToby Isaac     PetscInt     maxDegree;
524e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
525c330f8ffSToby Isaac     IS           heightIS;
5268c6c5593SMatthew G. Knepley 
52747923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
528e5e52638SMatthew G. Knepley     ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
529c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
5308c6c5593SMatthew G. Knepley     if (!h) {
5318c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
5328c6c5593SMatthew G. Knepley 
5338c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
5348c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
5358c6c5593SMatthew G. Knepley     }
536c330f8ffSToby Isaac     if (pEnd <= pStart) {
537c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
538c330f8ffSToby Isaac       continue;
539c330f8ffSToby Isaac     }
54047923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
54147923291SMatthew G. Knepley     totDim = 0;
54247923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
543bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
54447923291SMatthew G. Knepley         sp[f] = cellsp[f];
54547923291SMatthew G. Knepley       } else {
546bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
54747923291SMatthew G. Knepley         if (!sp[f]) continue;
54847923291SMatthew G. Knepley       }
54947923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
5509c3cf19fSMatthew G. Knepley       totDim += spDim;
55147923291SMatthew G. Knepley     }
552e5e52638SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr);
5538c6c5593SMatthew G. Knepley     if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim);
554c330f8ffSToby Isaac     if (!totDim) {
555c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
556c330f8ffSToby Isaac       continue;
557c330f8ffSToby Isaac     }
558c9731667SStefano Zampini     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);}
55947923291SMatthew G. Knepley     /* Loop over points at this height */
56069291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
56169291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
5628c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
5638c6c5593SMatthew G. Knepley     if (label) {
5648c6c5593SMatthew G. Knepley       PetscInt i;
56547923291SMatthew G. Knepley 
56647923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
567c330f8ffSToby Isaac         IS              pointIS, isectIS;
56847923291SMatthew G. Knepley         const PetscInt *points;
5698c6c5593SMatthew G. Knepley         PetscInt        n;
570c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
571c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
57247923291SMatthew G. Knepley 
57347923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
57447923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
575c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
576c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
577c330f8ffSToby Isaac         if (!isectIS) continue;
578c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
579c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
580b7260050SToby Isaac         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
581b7260050SToby Isaac         if (maxDegree <= 1) {
582c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
583c330f8ffSToby Isaac         }
584c330f8ffSToby Isaac         if (!quad) {
585c330f8ffSToby Isaac           if (!h && allPoints) {
586c330f8ffSToby Isaac             quad = allPoints;
587c330f8ffSToby Isaac             allPoints = NULL;
588c330f8ffSToby Isaac           } else {
589c330f8ffSToby Isaac             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
590c330f8ffSToby Isaac           }
591c330f8ffSToby Isaac         }
592c330f8ffSToby Isaac         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
59347923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
59447923291SMatthew G. Knepley           const PetscInt  point = points[p];
59547923291SMatthew G. Knepley 
596b420b6ccSMatthew G. Knepley           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
597c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
598c330f8ffSToby Isaac           ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
59947923291SMatthew G. Knepley           if (ierr) {
60047923291SMatthew G. Knepley             PetscErrorCode ierr2;
60169291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
60269291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
60347923291SMatthew G. Knepley             CHKERRQ(ierr);
60447923291SMatthew G. Knepley           }
605ca3d3a14SMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
606ba322698SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
60747923291SMatthew G. Knepley         }
608c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
609c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
610c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
611c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
612c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
61347923291SMatthew G. Knepley       }
6148c6c5593SMatthew G. Knepley     } else {
615c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
616c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
617c330f8ffSToby Isaac       IS              pointIS;
618c330f8ffSToby Isaac 
619c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
620b7260050SToby Isaac       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
621b7260050SToby Isaac       if (maxDegree <= 1) {
622c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
623c330f8ffSToby Isaac       }
624c330f8ffSToby Isaac       if (!quad) {
625c330f8ffSToby Isaac         if (!h && allPoints) {
626c330f8ffSToby Isaac           quad = allPoints;
627c330f8ffSToby Isaac           allPoints = NULL;
628c330f8ffSToby Isaac         } else {
629c330f8ffSToby Isaac           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
630c330f8ffSToby Isaac         }
631c330f8ffSToby Isaac       }
632c330f8ffSToby Isaac       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
6338c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
634b420b6ccSMatthew G. Knepley         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
635c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
636c330f8ffSToby Isaac         ierr = DMProjectPoint_Private(dm, probEff, chunkgeom, dmAux, probAux, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, Ncc, comps, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values);
6378c6c5593SMatthew G. Knepley         if (ierr) {
6388c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
63969291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
64069291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
6418c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
6428c6c5593SMatthew G. Knepley         }
643ca3d3a14SMatthew G. Knepley         if (transform) {ierr = DMPlexBasisTransformPoint_Internal(dm, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
644ba322698SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
6458c6c5593SMatthew G. Knepley       }
646c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
647c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
648c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
649c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
6508c6c5593SMatthew G. Knepley     }
651c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
65269291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
65369291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
65447923291SMatthew G. Knepley   }
6558c6c5593SMatthew G. Knepley   /* Cleanup */
6560c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
65774fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
65874fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
6590c898c18SMatthew G. Knepley 
6600c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
6610c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
6620c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
66374fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
66474fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
66574fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
6660c898c18SMatthew G. Knepley     }
6670c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
6680c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
66974fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
67074fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
67174fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
6720c898c18SMatthew G. Knepley     }
6730c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
6740c898c18SMatthew G. Knepley   }
675c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
676496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
6778c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
6788c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
67947923291SMatthew G. Knepley }
6808c6c5593SMatthew G. Knepley 
6818c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
6828c6c5593SMatthew G. Knepley {
6838c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
6848c6c5593SMatthew G. Knepley 
6858c6c5593SMatthew G. Knepley   PetscFunctionBegin;
6861c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
6878c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
6888c6c5593SMatthew G. Knepley }
6898c6c5593SMatthew G. Knepley 
6901c531cf8SMatthew 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)
6918c6c5593SMatthew G. Knepley {
6928c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
6938c6c5593SMatthew G. Knepley 
6948c6c5593SMatthew G. Knepley   PetscFunctionBegin;
6951c531cf8SMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr);
69647923291SMatthew G. Knepley   PetscFunctionReturn(0);
69747923291SMatthew G. Knepley }
69847923291SMatthew G. Knepley 
6998c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
70047923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
70147923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
70247923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
703191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
70447923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
70547923291SMatthew G. Knepley {
70647923291SMatthew G. Knepley   PetscErrorCode ierr;
70747923291SMatthew G. Knepley 
70847923291SMatthew G. Knepley   PetscFunctionBegin;
7091c531cf8SMatthew 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);
7108c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
71147923291SMatthew G. Knepley }
71247923291SMatthew G. Knepley 
7131c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
7148c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
7158c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
7168c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
717191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
7188c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
7198c6c5593SMatthew G. Knepley {
7208c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
72147923291SMatthew G. Knepley 
7228c6c5593SMatthew G. Knepley   PetscFunctionBegin;
7231c531cf8SMatthew 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);
72447923291SMatthew G. Knepley   PetscFunctionReturn(0);
72547923291SMatthew G. Knepley }
726