xref: /petsc/src/dm/impls/plex/plexproject.c (revision b18799e03d1720dd1c748f04b7f4e1e607858c24)
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;
10c330f8ffSToby Isaac   PetscBool      isAffine;
118c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
128c6c5593SMatthew G. Knepley 
138c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
14c330f8ffSToby Isaac   ierr = DMGetCoordinateDim(dm,&coordDim);CHKERRQ(ierr);
158c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
16496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
179c3cf19fSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
188c6c5593SMatthew G. Knepley   /* Get values for closure */
19c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
20c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
218c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
228c6c5593SMatthew G. Knepley 
238c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
248c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
258c6c5593SMatthew G. Knepley     if (funcs[f]) {
26c330f8ffSToby Isaac       if (isFE[f]) {
27c330f8ffSToby Isaac         PetscQuadrature   allPoints;
28c330f8ffSToby Isaac         PetscInt          q, dim, numPoints;
29c330f8ffSToby Isaac         const PetscReal   *points;
30c330f8ffSToby Isaac         PetscScalar       *pointEval;
31c330f8ffSToby Isaac         PetscReal         *x;
32c330f8ffSToby Isaac         DM                dm;
33c330f8ffSToby Isaac 
34c330f8ffSToby Isaac         ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
35c330f8ffSToby Isaac         ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
36c330f8ffSToby Isaac         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
37c330f8ffSToby Isaac         ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
38c330f8ffSToby Isaac         ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
39c330f8ffSToby Isaac         ierr = DMGetWorkArray(dm,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           }
49c330f8ffSToby Isaac           ierr = (*funcs[f])(coordDim,time,v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
50c330f8ffSToby Isaac         }
51c330f8ffSToby Isaac         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
52c330f8ffSToby Isaac         ierr = DMRestoreWorkArray(dm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
53c330f8ffSToby Isaac         ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
54c330f8ffSToby Isaac         v += spDim;
55c330f8ffSToby Isaac       } else {
56c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
57c330f8ffSToby Isaac           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
58c330f8ffSToby Isaac         }
59c330f8ffSToby Isaac       }
60c330f8ffSToby Isaac     } else {
61c330f8ffSToby Isaac       for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;}
628c6c5593SMatthew G. Knepley     }
639c3cf19fSMatthew G. Knepley   }
648c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
658c6c5593SMatthew G. Knepley }
668c6c5593SMatthew G. Knepley 
67c330f8ffSToby Isaac static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fegeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
680c898c18SMatthew G. Knepley                                                    PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
698c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
708c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
718c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
72191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
738c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
748c6c5593SMatthew G. Knepley {
758c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
769c3cf19fSMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL, *bc;
778c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
788c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
79191494d9SMatthew G. Knepley   const PetscScalar *constants;
808c6c5593SMatthew G. Knepley   PetscReal         *x;
81496733e2SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
82c330f8ffSToby Isaac   const PetscInt     dE = fegeom->dimEmbed;
83c330f8ffSToby Isaac   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
84c330f8ffSToby Isaac   PetscBool          isAffine;
858c6c5593SMatthew G. Knepley   PetscErrorCode     ierr;
868c6c5593SMatthew G. Knepley 
878c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
888c6c5593SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
89496733e2SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
90496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
918c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
928c6c5593SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
939c3cf19fSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
94496733e2SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
95191494d9SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
96e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
97496733e2SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
988c6c5593SMatthew G. Knepley   if (dmAux) {
9944171101SMatthew G. Knepley     PetscInt subp;
1001ba34526SMatthew G. Knepley 
10144171101SMatthew G. Knepley     ierr = DMPlexGetAuxiliaryPoint(dm, dmAux, p, &subp);CHKERRQ(ierr);
1021ba34526SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
1038c6c5593SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
104496733e2SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
105496733e2SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
106e87a4003SBarry Smith     ierr = DMGetSection(dmAux, &sectionAux);CHKERRQ(ierr);
1078c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
1088c6c5593SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
10911d189daSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
110496733e2SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
1111ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
1128c6c5593SMatthew G. Knepley   }
1138c6c5593SMatthew G. Knepley   /* Get values for closure */
114c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
1158c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
116c330f8ffSToby Isaac     PetscQuadrature   allPoints;
117c330f8ffSToby Isaac     PetscInt          q, dim, numPoints;
118c330f8ffSToby Isaac     const PetscReal   *points;
119c330f8ffSToby Isaac     PetscScalar       *pointEval;
120c330f8ffSToby Isaac     DM                dm;
121c330f8ffSToby Isaac 
1228c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
1238c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
124c330f8ffSToby Isaac     if (!funcs[f]) {
125be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
126c330f8ffSToby Isaac       continue;
127c330f8ffSToby Isaac     }
128c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
129c330f8ffSToby Isaac     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
130c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
131c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
1320c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
133c330f8ffSToby Isaac       const PetscReal *v0;
134c330f8ffSToby Isaac       const PetscReal *invJ;
135c330f8ffSToby Isaac 
136c330f8ffSToby Isaac       if (isAffine) {
137be1504a2SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
138c330f8ffSToby Isaac         v0 = x;
139c330f8ffSToby Isaac         invJ = fegeom->invJ;
1401c531cf8SMatthew G. Knepley       } else {
141c330f8ffSToby Isaac         v0 = &fegeom->v[tp*dE];
142c330f8ffSToby Isaac         invJ = &fegeom->invJ[tp*dE*dE];
1438c6c5593SMatthew G. Knepley       }
144c330f8ffSToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, invJ, coefficients, coefficients_t, u, u_x, u_t);
145e39fcb4eSStefano Zampini       if (probAux) EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);
146be1504a2SMatthew G. Knepley       (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, v0, numConstants, constants, &pointEval[Nc[f]*q]);
1471c531cf8SMatthew G. Knepley     }
148c330f8ffSToby Isaac     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
149c330f8ffSToby Isaac     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
150c330f8ffSToby Isaac     v += spDim;
1518c6c5593SMatthew G. Knepley   }
1528c6c5593SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
1538c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
1548c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1558c6c5593SMatthew G. Knepley }
1568c6c5593SMatthew G. Knepley 
157*b18799e0SMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS prob, DM dmAux, PetscDS probAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fegeom, PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
158*b18799e0SMatthew G. Knepley                                                      PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
159*b18799e0SMatthew G. Knepley                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
160*b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
161*b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
162*b18799e0SMatthew G. Knepley                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
163*b18799e0SMatthew G. Knepley                                                      PetscScalar values[])
164*b18799e0SMatthew G. Knepley {
165*b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
166*b18799e0SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *refSpaceDer, *refSpaceDerAux = NULL, *bc;
167*b18799e0SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
168*b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
169*b18799e0SMatthew G. Knepley   const PetscScalar *constants;
170*b18799e0SMatthew G. Knepley   PetscReal         *x;
171*b18799e0SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
172*b18799e0SMatthew G. Knepley   const PetscInt     dE = fegeom->dimEmbed;
173*b18799e0SMatthew G. Knepley   PetscInt           dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
174*b18799e0SMatthew G. Knepley   PetscBool          isAffine;
175*b18799e0SMatthew G. Knepley   PetscErrorCode     ierr;
176*b18799e0SMatthew G. Knepley 
177*b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
178*b18799e0SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
179*b18799e0SMatthew G. Knepley   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
180*b18799e0SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
181*b18799e0SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
182*b18799e0SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
183*b18799e0SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
184*b18799e0SMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
185*b18799e0SMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
186*b18799e0SMatthew G. Knepley   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
187*b18799e0SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
188*b18799e0SMatthew G. Knepley   if (dmAux) {
189*b18799e0SMatthew G. Knepley     DMLabel  spmap;
190*b18799e0SMatthew G. Knepley     PetscInt subp = p;
191*b18799e0SMatthew G. Knepley 
192*b18799e0SMatthew G. Knepley     /* If dm is a submesh, do not get subpoint */
193*b18799e0SMatthew G. Knepley     ierr = DMPlexGetSubpointMap(dm, &spmap);CHKERRQ(ierr);
194*b18799e0SMatthew G. Knepley     if (!spmap) {ierr = DMPlexGetSubpoint(dmAux, p, &subp);CHKERRQ(ierr);}
195*b18799e0SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
196*b18799e0SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
197*b18799e0SMatthew G. Knepley     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
198*b18799e0SMatthew G. Knepley     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
199*b18799e0SMatthew G. Knepley     ierr = DMGetSection(dmAux, &sectionAux);CHKERRQ(ierr);
200*b18799e0SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
201*b18799e0SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
202*b18799e0SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
203*b18799e0SMatthew G. Knepley     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
204*b18799e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
205*b18799e0SMatthew G. Knepley   }
206*b18799e0SMatthew G. Knepley   /* Get values for closure */
207*b18799e0SMatthew G. Knepley   isAffine = fegeom->isAffine;
208*b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
209*b18799e0SMatthew G. Knepley     PetscQuadrature   allPoints;
210*b18799e0SMatthew G. Knepley     PetscInt          q, dim, numPoints;
211*b18799e0SMatthew G. Knepley     const PetscReal   *points;
212*b18799e0SMatthew G. Knepley     PetscScalar       *pointEval;
213*b18799e0SMatthew G. Knepley     DM                dm;
214*b18799e0SMatthew G. Knepley 
215*b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
216*b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
217*b18799e0SMatthew G. Knepley     if (!funcs[f]) {
218*b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
219*b18799e0SMatthew G. Knepley       continue;
220*b18799e0SMatthew G. Knepley     }
221*b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
222*b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr);
223*b18799e0SMatthew G. Knepley     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
224*b18799e0SMatthew G. Knepley     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
225*b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
226*b18799e0SMatthew G. Knepley       const PetscReal *v0;
227*b18799e0SMatthew G. Knepley       const PetscReal *invJ;
228*b18799e0SMatthew G. Knepley       const PetscReal *n;
229*b18799e0SMatthew G. Knepley 
230*b18799e0SMatthew G. Knepley       if (isAffine) {
231*b18799e0SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
232*b18799e0SMatthew G. Knepley         v0 = x;
233*b18799e0SMatthew G. Knepley         invJ = fegeom->invJ;
234*b18799e0SMatthew G. Knepley         n  = fegeom->n;
235*b18799e0SMatthew G. Knepley       } else {
236*b18799e0SMatthew G. Knepley         v0 = &fegeom->v[tp*dE];
237*b18799e0SMatthew G. Knepley         invJ = &fegeom->invJ[tp*dE*dE];
238*b18799e0SMatthew G. Knepley         n  = &fegeom->n[tp*dE];
239*b18799e0SMatthew G. Knepley       }
240*b18799e0SMatthew G. Knepley       EvaluateFieldJets(dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, refSpaceDer, invJ, coefficients, coefficients_t, u, u_x, u_t);
241*b18799e0SMatthew G. Knepley       if (probAux) {EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, refSpaceDerAux, invJ, coefficientsAux, coefficientsAux_t, a, a_x, a_t);}
242*b18799e0SMatthew G. Knepley       (*funcs[f])(dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, v0, n, numConstants, constants, &pointEval[Nc[f]*q]);
243*b18799e0SMatthew G. Knepley     }
244*b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
245*b18799e0SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
246*b18799e0SMatthew G. Knepley     v += spDim;
247*b18799e0SMatthew G. Knepley   }
248*b18799e0SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dm, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
249*b18799e0SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
250*b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
251*b18799e0SMatthew G. Knepley }
252*b18799e0SMatthew G. Knepley 
253c330f8ffSToby 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[],
2541c531cf8SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscInt Ncc, const PetscInt comps[],
2551c531cf8SMatthew G. Knepley                                              PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux,
2568c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
2578c6c5593SMatthew G. Knepley {
2588c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
2598c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
2608c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
2618c6c5593SMatthew G. Knepley 
2628c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
2638c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2648c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
2658c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
2668c6c5593SMatthew G. Knepley   switch (type) {
2678c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
2688c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
269c330f8ffSToby 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;
2708c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
2718c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
272c330f8ffSToby Isaac     ierr = DMProjectPoint_Field_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
2730c898c18SMatthew G. Knepley                                         basisTab, basisDerTab, basisTabAux, basisDerTabAux,
2740c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
2750c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
2760c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
277191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
278*b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
279*b18799e0SMatthew G. Knepley     ierr = DMProjectPoint_BdField_Private(dm, prob, dmAux, probAux, time, localU, localA, fegeom, sp, p, Ncc, comps,
280*b18799e0SMatthew G. Knepley                                           basisTab, basisDerTab, basisTabAux, basisDerTabAux,
281*b18799e0SMatthew G. Knepley                                           (void (**)(PetscInt, PetscInt, PetscInt,
282*b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
283*b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
284*b18799e0SMatthew G. Knepley                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
2858c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
2868c6c5593SMatthew G. Knepley   }
2878c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2888c6c5593SMatthew G. Knepley }
2898c6c5593SMatthew G. Knepley 
290df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
291df5c1128SToby Isaac {
292df5c1128SToby Isaac   PetscReal      *points;
293df5c1128SToby Isaac   PetscInt       f, numPoints;
294df5c1128SToby Isaac   PetscErrorCode ierr;
295df5c1128SToby Isaac 
296df5c1128SToby Isaac   PetscFunctionBegin;
297df5c1128SToby Isaac   numPoints = 0;
298df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
299df5c1128SToby Isaac     if (funcs[f]) {
300df5c1128SToby Isaac       PetscQuadrature fAllPoints;
301df5c1128SToby Isaac       PetscInt        fNumPoints;
302df5c1128SToby Isaac 
303df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
304df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
305df5c1128SToby Isaac       numPoints += fNumPoints;
306df5c1128SToby Isaac     }
307df5c1128SToby Isaac   }
308df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
309df5c1128SToby Isaac   numPoints = 0;
310df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
311df5c1128SToby Isaac     PetscInt spDim;
312df5c1128SToby Isaac 
313df5c1128SToby Isaac     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
314df5c1128SToby Isaac     if (funcs[f]) {
315df5c1128SToby Isaac       PetscQuadrature fAllPoints;
31654ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
317df5c1128SToby Isaac       const PetscReal *fPoints;
318df5c1128SToby Isaac 
319df5c1128SToby Isaac       ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr);
32054ee0cceSMatthew G. Knepley       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
32154ee0cceSMatthew 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);
322df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
323df5c1128SToby Isaac       numPoints += fNumPoints;
324df5c1128SToby Isaac     }
325df5c1128SToby Isaac   }
326df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
327df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
328df5c1128SToby Isaac   PetscFunctionReturn(0);
329df5c1128SToby Isaac }
330df5c1128SToby Isaac 
331e5e52638SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *prob)
332e5e52638SMatthew G. Knepley {
333e5e52638SMatthew G. Knepley   DMLabel        depthLabel;
334e5e52638SMatthew G. Knepley   PetscInt       dim, cdepth, ls = -1, i;
335e5e52638SMatthew G. Knepley   PetscErrorCode ierr;
336e5e52638SMatthew G. Knepley 
337e5e52638SMatthew G. Knepley   PetscFunctionBegin;
338e5e52638SMatthew G. Knepley   if (lStart) *lStart = -1;
339e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
340e5e52638SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
341e5e52638SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
342e5e52638SMatthew G. Knepley   cdepth = dim - height;
343e5e52638SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
344e5e52638SMatthew G. Knepley     IS              pointIS;
345e5e52638SMatthew G. Knepley     const PetscInt *points;
346e5e52638SMatthew G. Knepley     PetscInt        pdepth;
347e5e52638SMatthew G. Knepley 
348e5e52638SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
349e5e52638SMatthew G. Knepley     if (!pointIS) continue; /* No points with that id on this process */
350e5e52638SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
351e5e52638SMatthew G. Knepley     ierr = DMLabelGetValue(depthLabel, points[0], &pdepth);CHKERRQ(ierr);
352e5e52638SMatthew G. Knepley     if (pdepth == cdepth) {
353e5e52638SMatthew G. Knepley       ls = points[0];
354e5e52638SMatthew G. Knepley       if (prob) {ierr = DMGetCellDS(dm, ls, prob);CHKERRQ(ierr);}
355e5e52638SMatthew G. Knepley     }
356e5e52638SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
357e5e52638SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
358e5e52638SMatthew G. Knepley     if (ls >= 0) break;
359e5e52638SMatthew G. Knepley   }
360e5e52638SMatthew G. Knepley   if (lStart) *lStart = ls;
361e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
362e5e52638SMatthew G. Knepley }
363e5e52638SMatthew G. Knepley 
3640de51b56SMatthew G. Knepley /*
3650de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
3660de51b56SMatthew G. Knepley 
3670de51b56SMatthew G. Knepley   There are several different scenarios:
3680de51b56SMatthew G. Knepley 
3690de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
3700de51b56SMatthew G. Knepley 
3710de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
3720de51b56SMatthew G. Knepley 
3730de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
3740de51b56SMatthew G. Knepley 
3750de51b56SMatthew 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.
3760de51b56SMatthew G. Knepley 
3770de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
3780de51b56SMatthew G. Knepley 
3790de51b56SMatthew 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.
3800de51b56SMatthew G. Knepley 
3810de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
3820de51b56SMatthew G. Knepley 
3830de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
3840de51b56SMatthew G. Knepley 
3850de51b56SMatthew 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.
3860de51b56SMatthew 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
3870de51b56SMatthew 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
3880de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
3890de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
3900de51b56SMatthew G. Knepley 
3910de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
3920de51b56SMatthew G. Knepley 
3930de51b56SMatthew 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.
3940de51b56SMatthew G. Knepley */
39546fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
3961c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
3978c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
3988c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
3998c6c5593SMatthew G. Knepley {
4008c6c5593SMatthew G. Knepley   DM              dmAux = NULL;
401e5e52638SMatthew G. Knepley   PetscDS         prob = NULL, probAux = NULL;
4028c6c5593SMatthew G. Knepley   Vec             localA = NULL;
40347923291SMatthew G. Knepley   PetscSection    section;
4048c6c5593SMatthew G. Knepley   PetscDualSpace *sp, *cellsp;
4050c898c18SMatthew G. Knepley   PetscReal     **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL;
4068c6c5593SMatthew G. Knepley   PetscInt       *Nc;
407c330f8ffSToby Isaac   PetscInt        dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfAux = 0, f;
40888aca1feSMatthew G. Knepley   PetscBool      *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE;
409c330f8ffSToby Isaac   DMField         coordField;
410c330f8ffSToby Isaac   DMLabel         depthLabel;
411c330f8ffSToby Isaac   PetscQuadrature allPoints = NULL;
41247923291SMatthew G. Knepley   PetscErrorCode  ierr;
41347923291SMatthew G. Knepley 
41447923291SMatthew G. Knepley   PetscFunctionBegin;
41588aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
41688aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
4178c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4182af58022SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(dm, &minHeight);CHKERRQ(ierr);
4190de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
4200de51b56SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
421e39fcb4eSStefano Zampini     if (dmAux && !localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
42288aca1feSMatthew G. Knepley     if (!minHeight && dmAux) {
42388aca1feSMatthew G. Knepley       DMLabel spmap;
42488aca1feSMatthew G. Knepley 
42588aca1feSMatthew G. Knepley       /* If dmAux is a surface, then force the projection to take place over a surface */
42688aca1feSMatthew G. Knepley       ierr = DMPlexGetSubpointMap(dmAux, &spmap);CHKERRQ(ierr);
427e39fcb4eSStefano Zampini       if (spmap) {
428e39fcb4eSStefano Zampini         ierr = DMPlexGetVTKCellHeight(dmAux, &minHeight);CHKERRQ(ierr);
429e39fcb4eSStefano Zampini         auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
430e39fcb4eSStefano Zampini       }
43188aca1feSMatthew G. Knepley     }
4320de51b56SMatthew G. Knepley   }
433e5e52638SMatthew G. Knepley   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
434e5e52638SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(dm,&depthLabel);CHKERRQ(ierr);
4358c6c5593SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(dm, &maxHeight);CHKERRQ(ierr);
4362af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
437e39fcb4eSStefano Zampini   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
438e5e52638SMatthew G. Knepley   ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, 0, NULL, &prob);CHKERRQ(ierr);
439e5e52638SMatthew G. Knepley   if (!prob) {ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);}
440e5e52638SMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
4418c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
442e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
4430c898c18SMatthew G. Knepley   if (dmAux) {
4440c898c18SMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
4450c898c18SMatthew G. Knepley     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
4460c898c18SMatthew G. Knepley   }
447496733e2SMatthew G. Knepley   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
448496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
4498c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
4508c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
4518c6c5593SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
4528c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
45347923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
45447923291SMatthew G. Knepley     PetscObject  obj;
45547923291SMatthew G. Knepley     PetscClassId id;
45647923291SMatthew G. Knepley 
457e5e52638SMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
45847923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
45947923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
46047923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
46147923291SMatthew G. Knepley 
46247923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
46347923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
46447923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
46547923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
46647923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
46747923291SMatthew G. Knepley 
46847923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
46947923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
47047923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
47147923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
47247923291SMatthew G. Knepley   }
473c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
4740c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
47574fc349aSMatthew G. Knepley     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
47674fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
4774a3e9fdbSToby Isaac     const PetscReal *points;
4784a3e9fdbSToby Isaac     PetscInt         numPoints;
4790c898c18SMatthew G. Knepley 
4802af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
48154ee0cceSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
48254ee0cceSMatthew G. Knepley       if (!effectiveHeight) {sp[f] = cellsp[f];}
48354ee0cceSMatthew G. Knepley       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
48454ee0cceSMatthew G. Knepley     }
48554ee0cceSMatthew G. Knepley     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
486c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
4870c898c18SMatthew G. Knepley     ierr = PetscMalloc4(Nf, &basisTab, Nf, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr);
4880c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
4890c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
4900c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
49174fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
49274fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
49374fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
4940c898c18SMatthew G. Knepley     }
4950c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
4960c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
49774fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
49874fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
49974fc349aSMatthew G. Knepley       ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
5000c898c18SMatthew G. Knepley     }
5010c898c18SMatthew G. Knepley   }
50247923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
5032af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
50488aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
5050de51b56SMatthew G. Knepley     PetscDS      probEff         = prob;
5068c6c5593SMatthew G. Knepley     PetscScalar *values;
507b7260050SToby Isaac     PetscBool   *fieldActive;
508b7260050SToby Isaac     PetscInt     maxDegree;
509e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
510c330f8ffSToby Isaac     IS           heightIS;
5118c6c5593SMatthew G. Knepley 
51247923291SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);CHKERRQ(ierr);
513e5e52638SMatthew G. Knepley     ierr = DMGetFirstLabelEntry_Private(dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
514c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
5158c6c5593SMatthew G. Knepley     if (!h) {
5168c6c5593SMatthew G. Knepley       PetscInt cEndInterior;
5178c6c5593SMatthew G. Knepley 
5188c6c5593SMatthew G. Knepley       ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
5198c6c5593SMatthew G. Knepley       pEnd = cEndInterior < 0 ? pEnd : cEndInterior;
5208c6c5593SMatthew G. Knepley     }
521c330f8ffSToby Isaac     if (pEnd <= pStart) {
522c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
523c330f8ffSToby Isaac       continue;
524c330f8ffSToby Isaac     }
52547923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
52647923291SMatthew G. Knepley     totDim = 0;
52747923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
528bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
52947923291SMatthew G. Knepley         sp[f] = cellsp[f];
53047923291SMatthew G. Knepley       } else {
531bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
53247923291SMatthew G. Knepley         if (!sp[f]) continue;
53347923291SMatthew G. Knepley       }
53447923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
5359c3cf19fSMatthew G. Knepley       totDim += spDim;
53647923291SMatthew G. Knepley     }
537e5e52638SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr);
5388c6c5593SMatthew 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);
539c330f8ffSToby Isaac     if (!totDim) {
540c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
541c330f8ffSToby Isaac       continue;
542c330f8ffSToby Isaac     }
543c9731667SStefano Zampini     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(prob, effectiveHeight, &probEff);CHKERRQ(ierr);}
54447923291SMatthew G. Knepley     /* Loop over points at this height */
54569291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
54669291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
5478c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
5488c6c5593SMatthew G. Knepley     if (label) {
5498c6c5593SMatthew G. Knepley       PetscInt i;
55047923291SMatthew G. Knepley 
55147923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
552c330f8ffSToby Isaac         IS              pointIS, isectIS;
55347923291SMatthew G. Knepley         const PetscInt *points;
5548c6c5593SMatthew G. Knepley         PetscInt        n;
555c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
556c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
55747923291SMatthew G. Knepley 
55847923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
55947923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
560c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
561c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
562c330f8ffSToby Isaac         if (!isectIS) continue;
563c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
564c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
565b7260050SToby Isaac         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
566b7260050SToby Isaac         if (maxDegree <= 1) {
567c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
568c330f8ffSToby Isaac         }
569c330f8ffSToby Isaac         if (!quad) {
570c330f8ffSToby Isaac           if (!h && allPoints) {
571c330f8ffSToby Isaac             quad = allPoints;
572c330f8ffSToby Isaac             allPoints = NULL;
573c330f8ffSToby Isaac           } else {
574c330f8ffSToby Isaac             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
575c330f8ffSToby Isaac           }
576c330f8ffSToby Isaac         }
577c330f8ffSToby Isaac         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
57847923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
57947923291SMatthew G. Knepley           const PetscInt  point = points[p];
58047923291SMatthew G. Knepley 
581b420b6ccSMatthew G. Knepley           ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
582c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
583c330f8ffSToby 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);
58447923291SMatthew G. Knepley           if (ierr) {
58547923291SMatthew G. Knepley             PetscErrorCode ierr2;
58669291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
58769291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
58847923291SMatthew G. Knepley             CHKERRQ(ierr);
58947923291SMatthew G. Knepley           }
590ba322698SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr);
59147923291SMatthew G. Knepley         }
592c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
593c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
594c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
595c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
596c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
59747923291SMatthew G. Knepley       }
5988c6c5593SMatthew G. Knepley     } else {
599c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
600c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
601c330f8ffSToby Isaac       IS              pointIS;
602c330f8ffSToby Isaac 
603c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
604b7260050SToby Isaac       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
605b7260050SToby Isaac       if (maxDegree <= 1) {
606c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
607c330f8ffSToby Isaac       }
608c330f8ffSToby Isaac       if (!quad) {
609c330f8ffSToby Isaac         if (!h && allPoints) {
610c330f8ffSToby Isaac           quad = allPoints;
611c330f8ffSToby Isaac           allPoints = NULL;
612c330f8ffSToby Isaac         } else {
613c330f8ffSToby Isaac           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
614c330f8ffSToby Isaac         }
615c330f8ffSToby Isaac       }
616c330f8ffSToby Isaac       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&fegeom);CHKERRQ(ierr);
6178c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
618b420b6ccSMatthew G. Knepley         ierr = PetscMemzero(values, numValues * sizeof(PetscScalar));CHKERRQ(ierr);
619c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
620c330f8ffSToby 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);
6218c6c5593SMatthew G. Knepley         if (ierr) {
6228c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
62369291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
62469291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
6258c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
6268c6c5593SMatthew G. Knepley         }
627ba322698SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr);
6288c6c5593SMatthew G. Knepley       }
629c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
630c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
631c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
632c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
6338c6c5593SMatthew G. Knepley     }
634c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
63569291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
63669291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
63747923291SMatthew G. Knepley   }
6388c6c5593SMatthew G. Knepley   /* Cleanup */
6390c898c18SMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) {
64074fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
64174fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
6420c898c18SMatthew G. Knepley 
6430c898c18SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
6440c898c18SMatthew G. Knepley       if (!isFE[f]) continue;
6450c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fem);CHKERRQ(ierr);
64674fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
64774fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
64874fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr);
6490c898c18SMatthew G. Knepley     }
6500c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
6510c898c18SMatthew G. Knepley       ierr = PetscDSGetDiscretization(probAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
65274fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
65374fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
65474fc349aSMatthew G. Knepley       ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr);
6550c898c18SMatthew G. Knepley     }
6560c898c18SMatthew G. Knepley     ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr);
6570c898c18SMatthew G. Knepley   }
658c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
659496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
6608c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
6618c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
66247923291SMatthew G. Knepley }
6638c6c5593SMatthew G. Knepley 
6648c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
6658c6c5593SMatthew G. Knepley {
6668c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
6678c6c5593SMatthew G. Knepley 
6688c6c5593SMatthew G. Knepley   PetscFunctionBegin;
6691c531cf8SMatthew 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);
6708c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
6718c6c5593SMatthew G. Knepley }
6728c6c5593SMatthew G. Knepley 
6731c531cf8SMatthew 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)
6748c6c5593SMatthew G. Knepley {
6758c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
6768c6c5593SMatthew G. Knepley 
6778c6c5593SMatthew G. Knepley   PetscFunctionBegin;
6781c531cf8SMatthew 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);
67947923291SMatthew G. Knepley   PetscFunctionReturn(0);
68047923291SMatthew G. Knepley }
68147923291SMatthew G. Knepley 
6828c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
68347923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
68447923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
68547923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
686191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
68747923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
68847923291SMatthew G. Knepley {
68947923291SMatthew G. Knepley   PetscErrorCode ierr;
69047923291SMatthew G. Knepley 
69147923291SMatthew G. Knepley   PetscFunctionBegin;
6921c531cf8SMatthew 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);
6938c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
69447923291SMatthew G. Knepley }
69547923291SMatthew G. Knepley 
6961c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
6978c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
6988c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
6998c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
700191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
7018c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
7028c6c5593SMatthew G. Knepley {
7038c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
70447923291SMatthew G. Knepley 
7058c6c5593SMatthew G. Knepley   PetscFunctionBegin;
7061c531cf8SMatthew 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);
70747923291SMatthew G. Knepley   PetscFunctionReturn(0);
70847923291SMatthew G. Knepley }
709