xref: /petsc/src/dm/impls/plex/plexproject.c (revision 27f02ce888cd6d83c79d307cd15e2865fe9f0ab4)
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 
5a6e0b375SMatthew G. Knepley /*
6a6e0b375SMatthew G. Knepley   DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point
7a6e0b375SMatthew G. Knepley 
8a6e0b375SMatthew G. Knepley   Input Parameters:
9a6e0b375SMatthew G. Knepley + dm     - The output DM
10a6e0b375SMatthew G. Knepley . ds     - The output DS
11a6e0b375SMatthew G. Knepley . dmIn   - The input DM
12a6e0b375SMatthew G. Knepley . dsIn   - The input DS
13a6e0b375SMatthew G. Knepley . time   - The time for this evaluation
14a6e0b375SMatthew G. Knepley . fegeom - The FE geometry for this point
15a6e0b375SMatthew G. Knepley . fvgeom - The FV geometry for this point
16a6e0b375SMatthew G. Knepley . isFE   - Flag indicating whether each output field has an FE discretization
17a6e0b375SMatthew G. Knepley . sp     - The output PetscDualSpace for each field
18a6e0b375SMatthew G. Knepley . funcs  - The evaluation function for each field
19a6e0b375SMatthew G. Knepley - ctxs   - The user context for each field
20a6e0b375SMatthew G. Knepley 
21a6e0b375SMatthew G. Knepley   Output Parameter:
22a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space
23a6e0b375SMatthew G. Knepley 
24a6e0b375SMatthew G. Knepley   Level: developer
25a6e0b375SMatthew G. Knepley 
26a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private()
27a6e0b375SMatthew G. Knepley */
28a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[],
298c6c5593SMatthew G. Knepley                                                   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs,
308c6c5593SMatthew G. Knepley                                                   PetscScalar values[])
3147923291SMatthew G. Knepley {
32a6e0b375SMatthew G. Knepley   PetscInt       coordDim, Nf, *Nc, f, spDim, d, v, tp;
33*27f02ce8SMatthew G. Knepley   PetscBool      isAffine, isHybrid, transform;
348c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
358c6c5593SMatthew G. Knepley 
368c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
37a6e0b375SMatthew G. Knepley   ierr = DMGetCoordinateDim(dmIn, &coordDim);CHKERRQ(ierr);
38a6e0b375SMatthew G. Knepley   ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr);
39a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
40a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
41*27f02ce8SMatthew G. Knepley   ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr);
428c6c5593SMatthew G. Knepley   /* Get values for closure */
43c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
44c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
458c6c5593SMatthew G. Knepley     void * const ctx = ctxs ? ctxs[f] : NULL;
468c6c5593SMatthew G. Knepley 
478c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
488c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
498c6c5593SMatthew G. Knepley     if (funcs[f]) {
50c330f8ffSToby Isaac       if (isFE[f]) {
51c330f8ffSToby Isaac         PetscQuadrature   allPoints;
52c330f8ffSToby Isaac         PetscInt          q, dim, numPoints;
53c330f8ffSToby Isaac         const PetscReal   *points;
54c330f8ffSToby Isaac         PetscScalar       *pointEval;
55c330f8ffSToby Isaac         PetscReal         *x;
56ca3d3a14SMatthew G. Knepley         DM                rdm;
57c330f8ffSToby Isaac 
58ca3d3a14SMatthew G. Knepley         ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr);
59b4457527SToby Isaac         ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
60c330f8ffSToby Isaac         ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
61ca3d3a14SMatthew G. Knepley         ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
62ca3d3a14SMatthew G. Knepley         ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
63c330f8ffSToby Isaac         for (q = 0; q < numPoints; q++, tp++) {
64c330f8ffSToby Isaac           const PetscReal *v0;
65c330f8ffSToby Isaac 
66c330f8ffSToby Isaac           if (isAffine) {
67c330f8ffSToby Isaac             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x);
68c330f8ffSToby Isaac             v0 = x;
698c6c5593SMatthew G. Knepley           } else {
70c330f8ffSToby Isaac             v0 = &fegeom->v[tp*coordDim];
718c6c5593SMatthew G. Knepley           }
72a6e0b375SMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;}
73c330f8ffSToby Isaac           ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
74c330f8ffSToby Isaac         }
754bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
762edcad52SToby Isaac         ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr);
77c330f8ffSToby Isaac         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
78ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
79ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
80c330f8ffSToby Isaac         v += spDim;
81*27f02ce8SMatthew G. Knepley         if (isHybrid && (f < Nf-1)) {
82*27f02ce8SMatthew G. Knepley           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
83*27f02ce8SMatthew G. Knepley         }
84c330f8ffSToby Isaac       } else {
85c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
86c330f8ffSToby Isaac           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
87c330f8ffSToby Isaac         }
88c330f8ffSToby Isaac       }
89c330f8ffSToby Isaac     } else {
90*27f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
91*27f02ce8SMatthew G. Knepley       if (isHybrid && (f < Nf-1)) {
92*27f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
93*27f02ce8SMatthew G. Knepley       }
948c6c5593SMatthew G. Knepley     }
959c3cf19fSMatthew G. Knepley   }
968c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
978c6c5593SMatthew G. Knepley }
988c6c5593SMatthew G. Knepley 
99a6e0b375SMatthew G. Knepley /*
100a6e0b375SMatthew G. Knepley   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
101a6e0b375SMatthew G. Knepley 
102a6e0b375SMatthew G. Knepley   Input Parameters:
103a6e0b375SMatthew G. Knepley + dm             - The output DM
104a6e0b375SMatthew G. Knepley . ds             - The output DS
105a6e0b375SMatthew G. Knepley . dmIn           - The input DM
106a6e0b375SMatthew G. Knepley . dsIn           - The input DS
107a6e0b375SMatthew G. Knepley . dmAux          - The auxiliary DM, which is always for the input space
108a6e0b375SMatthew G. Knepley . dsAux          - The auxiliary DS, which is always for the input space
109a6e0b375SMatthew G. Knepley . time           - The time for this evaluation
110a6e0b375SMatthew G. Knepley . localU         - The local solution
111a6e0b375SMatthew G. Knepley . localA         - The local auziliary fields
112a6e0b375SMatthew G. Knepley . cgeom          - The FE geometry for this point
113a6e0b375SMatthew G. Knepley . sp             - The output PetscDualSpace for each field
114a6e0b375SMatthew G. Knepley . p              - The point in the output DM
115ef0bb6c7SMatthew G. Knepley . T              - Input basis and derviatives for each field tabulated on the quadrature points
116ef0bb6c7SMatthew G. Knepley . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
117a6e0b375SMatthew G. Knepley . funcs          - The evaluation function for each field
118a6e0b375SMatthew G. Knepley - ctxs           - The user context for each field
119a6e0b375SMatthew G. Knepley 
120a6e0b375SMatthew G. Knepley   Output Parameter:
121a6e0b375SMatthew G. Knepley . values         - The value for each dual basis vector in the output dual space
122a6e0b375SMatthew G. Knepley 
123a6e0b375SMatthew G. Knepley   Note: Not supported for FV
124a6e0b375SMatthew G. Knepley 
125a6e0b375SMatthew G. Knepley   Level: developer
126a6e0b375SMatthew G. Knepley 
127a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private()
128a6e0b375SMatthew G. Knepley */
129a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p,
130ef0bb6c7SMatthew G. Knepley                                                    PetscTabulation *T, PetscTabulation *TAux,
1318c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
1328c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1338c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
134191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
1358c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
1368c6c5593SMatthew G. Knepley {
1378c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
1384bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
1398c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
1408c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
141191494d9SMatthew G. Knepley   const PetscScalar *constants;
1428c6c5593SMatthew G. Knepley   PetscReal         *x;
143ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
1444bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
1454bee2e38SMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed;
146ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
147*27f02ce8SMatthew G. Knepley   PetscBool          isAffine, isHybrid, transform;
1488c6c5593SMatthew G. Knepley   PetscErrorCode     ierr;
1498c6c5593SMatthew G. Knepley 
1508c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
151a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
152a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
153*27f02ce8SMatthew G. Knepley   ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr);
154a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
155a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(dsIn, &uOff);CHKERRQ(ierr);
156a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x);CHKERRQ(ierr);
157a6e0b375SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
158a6e0b375SMatthew G. Knepley   ierr = PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
159a6e0b375SMatthew G. Knepley   ierr = PetscDSGetConstants(dsIn, &numConstants, &constants);CHKERRQ(ierr);
160a6e0b375SMatthew G. Knepley   ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr);
161a6e0b375SMatthew G. Knepley   ierr = DMGetLocalSection(dmIn, &section);CHKERRQ(ierr);
162a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp);CHKERRQ(ierr);
163a6e0b375SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
1648c6c5593SMatthew G. Knepley   if (dmAux) {
16544171101SMatthew G. Knepley     PetscInt subp;
1661ba34526SMatthew G. Knepley 
167a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
168a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
16992fd8e1eSJed Brown     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
170a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
171a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
172a6e0b375SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
1731ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
1748c6c5593SMatthew G. Knepley   }
1758c6c5593SMatthew G. Knepley   /* Get values for closure */
1764bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
177*27f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
178*27f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
1794bee2e38SMatthew G. Knepley   if (isAffine) {
1804bee2e38SMatthew G. Knepley     fegeom.v    = x;
1814bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
1824bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
1834bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
1844bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
1854bee2e38SMatthew G. Knepley   }
1868c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
187c330f8ffSToby Isaac     PetscQuadrature   allPoints;
188c330f8ffSToby Isaac     PetscInt          q, dim, numPoints;
189c330f8ffSToby Isaac     const PetscReal   *points;
190c330f8ffSToby Isaac     PetscScalar       *pointEval;
191c330f8ffSToby Isaac     DM                dm;
192c330f8ffSToby Isaac 
1938c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
1948c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
195c330f8ffSToby Isaac     if (!funcs[f]) {
196be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
197*27f02ce8SMatthew G. Knepley       if (isHybrid && (f < Nf-1)) {
198*27f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
199*27f02ce8SMatthew G. Knepley       }
200c330f8ffSToby Isaac       continue;
201c330f8ffSToby Isaac     }
202c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
203b4457527SToby Isaac     ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
204c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
205c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
2060c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
207c330f8ffSToby Isaac       if (isAffine) {
208ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
2091c531cf8SMatthew G. Knepley       } else {
2104bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp*dE];
2114bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp*dE*dE];
2124bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
2134bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
2148c6c5593SMatthew G. Knepley       }
215ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
216ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
217a6e0b375SMatthew G. Knepley       if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);}
218a6e0b375SMatthew G. Knepley       (*funcs[f])(dE, NfIn, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, time, fegeom.v, numConstants, constants, &pointEval[Nc[f]*q]);
2191c531cf8SMatthew G. Knepley     }
220c330f8ffSToby Isaac     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
221c330f8ffSToby Isaac     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
222c330f8ffSToby Isaac     v += spDim;
223*27f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
224*27f02ce8SMatthew G. Knepley     if (isHybrid && (f < Nf-1)) {
225*27f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
226*27f02ce8SMatthew G. Knepley     }
2278c6c5593SMatthew G. Knepley   }
228a6e0b375SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
2298c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
2308c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2318c6c5593SMatthew G. Knepley }
2328c6c5593SMatthew G. Knepley 
233a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p,
234ef0bb6c7SMatthew G. Knepley                                                      PetscTabulation *T, PetscTabulation *TAux,
235b18799e0SMatthew G. Knepley                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
236b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
237b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
238b18799e0SMatthew G. Knepley                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
239b18799e0SMatthew G. Knepley                                                      PetscScalar values[])
240b18799e0SMatthew G. Knepley {
241b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
2424bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
243b18799e0SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
244b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
245b18799e0SMatthew G. Knepley   const PetscScalar *constants;
246b18799e0SMatthew G. Knepley   PetscReal         *x;
247ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
2489f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
2494bee2e38SMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed;
250ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
251b18799e0SMatthew G. Knepley   PetscBool          isAffine;
252b18799e0SMatthew G. Knepley   PetscErrorCode     ierr;
253b18799e0SMatthew G. Knepley 
254b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
255a6e0b375SMatthew G. Knepley   if (dm != dmIn) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
256a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
257a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
258a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
259a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
260a6e0b375SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
261a6e0b375SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
262a6e0b375SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
26392fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
264a6e0b375SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
265b18799e0SMatthew G. Knepley   if (dmAux) {
266a6e0b375SMatthew G. Knepley     PetscInt subp;
267b18799e0SMatthew G. Knepley 
268a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
269a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
27092fd8e1eSJed Brown     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
271a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
272a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
273a6e0b375SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
274b18799e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
275b18799e0SMatthew G. Knepley   }
276b18799e0SMatthew G. Knepley   /* Get values for closure */
2774bee2e38SMatthew G. Knepley   isAffine = fgeom->isAffine;
2784dcf62a8SSatish Balay   fegeom.n  = 0;
2794dcf62a8SSatish Balay   fegeom.J  = 0;
2804dcf62a8SSatish Balay   fegeom.v  = 0;
2814dcf62a8SSatish Balay   fegeom.xi = 0;
282a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
283a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
2844bee2e38SMatthew G. Knepley   if (isAffine) {
2854bee2e38SMatthew G. Knepley     fegeom.v    = x;
2864bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
2874bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
2884bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
2894bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
2904bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
2919f209ee4SMatthew G. Knepley 
2929f209ee4SMatthew G. Knepley     cgeom.J     = fgeom->suppJ[0];
2939f209ee4SMatthew G. Knepley     cgeom.invJ  = fgeom->suppInvJ[0];
2949f209ee4SMatthew G. Knepley     cgeom.detJ  = fgeom->suppDetJ[0];
2954bee2e38SMatthew G. Knepley   }
296b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
297b18799e0SMatthew G. Knepley     PetscQuadrature   allPoints;
298b18799e0SMatthew G. Knepley     PetscInt          q, dim, numPoints;
299b18799e0SMatthew G. Knepley     const PetscReal   *points;
300b18799e0SMatthew G. Knepley     PetscScalar       *pointEval;
301b18799e0SMatthew G. Knepley     DM                dm;
302b18799e0SMatthew G. Knepley 
303b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
304b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
305b18799e0SMatthew G. Knepley     if (!funcs[f]) {
306b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
307b18799e0SMatthew G. Knepley       continue;
308b18799e0SMatthew G. Knepley     }
309b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
310b4457527SToby Isaac     ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
311b18799e0SMatthew G. Knepley     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
312b18799e0SMatthew G. Knepley     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
313b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
314b18799e0SMatthew G. Knepley       if (isAffine) {
315ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
316b18799e0SMatthew G. Knepley       } else {
3173fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp*dE];
3189f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp*dE*dE];
3199f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp*dE*dE];
3204bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
3214bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp*dE];
3229f209ee4SMatthew G. Knepley 
3239f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][tp*dE*dE];
3249f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][tp*dE*dE];
3259f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][tp];
326b18799e0SMatthew G. Knepley       }
327a6e0b375SMatthew G. Knepley       /* TODO We should use cgeom here, instead of fegeom, however the geometry coming in through fgeom does not have the support cell geometry */
328ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
329ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
3304bee2e38SMatthew 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]);
331b18799e0SMatthew G. Knepley     }
332b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
333b18799e0SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
334b18799e0SMatthew G. Knepley     v += spDim;
335b18799e0SMatthew G. Knepley   }
336a6e0b375SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
337b18799e0SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
338b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
339b18799e0SMatthew G. Knepley }
340b18799e0SMatthew G. Knepley 
341a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscFEGeom *fegeom, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[],
342ef0bb6c7SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux,
3438c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
3448c6c5593SMatthew G. Knepley {
3458c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
3468c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
3478c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
3488c6c5593SMatthew G. Knepley 
3498c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
3508c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
3518c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
3528c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
3538c6c5593SMatthew G. Knepley   switch (type) {
3548c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
3558c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
356a6e0b375SMatthew G. Knepley     ierr = DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break;
3578c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
3588c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
359ef0bb6c7SMatthew G. Knepley     ierr = DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
3600c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
3610c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3620c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
363191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
364b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
365ef0bb6c7SMatthew G. Knepley     ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
366b18799e0SMatthew G. Knepley                                           (void (**)(PetscInt, PetscInt, PetscInt,
367b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
368b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
369b18799e0SMatthew G. Knepley                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
3708c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
3718c6c5593SMatthew G. Knepley   }
3728c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
3738c6c5593SMatthew G. Knepley }
3748c6c5593SMatthew G. Knepley 
375df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
376df5c1128SToby Isaac {
377df5c1128SToby Isaac   PetscReal      *points;
378df5c1128SToby Isaac   PetscInt       f, numPoints;
379df5c1128SToby Isaac   PetscErrorCode ierr;
380df5c1128SToby Isaac 
381df5c1128SToby Isaac   PetscFunctionBegin;
382df5c1128SToby Isaac   numPoints = 0;
383df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
384df5c1128SToby Isaac     if (funcs[f]) {
385df5c1128SToby Isaac       PetscQuadrature fAllPoints;
386df5c1128SToby Isaac       PetscInt        fNumPoints;
387df5c1128SToby Isaac 
388b4457527SToby Isaac       ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr);
389df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
390df5c1128SToby Isaac       numPoints += fNumPoints;
391df5c1128SToby Isaac     }
392df5c1128SToby Isaac   }
393df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
394df5c1128SToby Isaac   numPoints = 0;
395df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
396df5c1128SToby Isaac     PetscInt spDim;
397df5c1128SToby Isaac 
398df5c1128SToby Isaac     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
399df5c1128SToby Isaac     if (funcs[f]) {
400df5c1128SToby Isaac       PetscQuadrature fAllPoints;
40154ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
402df5c1128SToby Isaac       const PetscReal *fPoints;
403df5c1128SToby Isaac 
404b4457527SToby Isaac       ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr);
40554ee0cceSMatthew G. Knepley       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
40654ee0cceSMatthew 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);
407df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
408df5c1128SToby Isaac       numPoints += fNumPoints;
409df5c1128SToby Isaac     }
410df5c1128SToby Isaac   }
411df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
412df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
413df5c1128SToby Isaac   PetscFunctionReturn(0);
414df5c1128SToby Isaac }
415df5c1128SToby Isaac 
416a6e0b375SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *ds)
417e5e52638SMatthew G. Knepley {
418a6e0b375SMatthew G. Knepley   DM              plex;
419a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
420e5e52638SMatthew G. Knepley   DMLabel         depthLabel;
421e5e52638SMatthew G. Knepley   PetscInt        dim, cdepth, ls = -1, i;
422e5e52638SMatthew G. Knepley   PetscErrorCode  ierr;
423e5e52638SMatthew G. Knepley 
424e5e52638SMatthew G. Knepley   PetscFunctionBegin;
425e5e52638SMatthew G. Knepley   if (lStart) *lStart = -1;
426e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
427a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dm, odm, &enc);CHKERRQ(ierr);
428e5e52638SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
429a6e0b375SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
430a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
431e5e52638SMatthew G. Knepley   cdepth = dim - height;
432e5e52638SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
433e5e52638SMatthew G. Knepley     IS              pointIS;
434e5e52638SMatthew G. Knepley     const PetscInt *points;
435a6e0b375SMatthew G. Knepley     PetscInt        pdepth, point;
436e5e52638SMatthew G. Knepley 
437e5e52638SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
438e5e52638SMatthew G. Knepley     if (!pointIS) continue; /* No points with that id on this process */
439e5e52638SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
440a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dm, odm, enc, points[0], &point);CHKERRQ(ierr);
441a6e0b375SMatthew G. Knepley     ierr = DMLabelGetValue(depthLabel, point, &pdepth);CHKERRQ(ierr);
442e5e52638SMatthew G. Knepley     if (pdepth == cdepth) {
443a6e0b375SMatthew G. Knepley       ls = point;
444a6e0b375SMatthew G. Knepley       if (ds) {ierr = DMGetCellDS(dm, ls, ds);CHKERRQ(ierr);}
445e5e52638SMatthew G. Knepley     }
446e5e52638SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
447e5e52638SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
448e5e52638SMatthew G. Knepley     if (ls >= 0) break;
449e5e52638SMatthew G. Knepley   }
450e5e52638SMatthew G. Knepley   if (lStart) *lStart = ls;
451a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
452e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
453e5e52638SMatthew G. Knepley }
454e5e52638SMatthew G. Knepley 
4550de51b56SMatthew G. Knepley /*
4560de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
4570de51b56SMatthew G. Knepley 
4580de51b56SMatthew G. Knepley   There are several different scenarios:
4590de51b56SMatthew G. Knepley 
4600de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
4610de51b56SMatthew G. Knepley 
4620de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
4630de51b56SMatthew G. Knepley 
4640de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
4650de51b56SMatthew G. Knepley 
4660de51b56SMatthew 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.
4670de51b56SMatthew G. Knepley 
4680de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
4690de51b56SMatthew G. Knepley 
4700de51b56SMatthew 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.
4710de51b56SMatthew G. Knepley 
4720de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
4730de51b56SMatthew G. Knepley 
4740de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
4750de51b56SMatthew G. Knepley 
4760de51b56SMatthew 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.
4770de51b56SMatthew 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
4785f790a90SMatthew G. Knepley       is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When the effective height is nonzero, we need to extract
4790de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
4800de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
4810de51b56SMatthew G. Knepley 
4820de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
4830de51b56SMatthew G. Knepley 
4840de51b56SMatthew 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.
4850de51b56SMatthew G. Knepley */
48646fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
4871c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
4888c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
4898c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
4908c6c5593SMatthew G. Knepley {
491a6e0b375SMatthew G. Knepley   DM                 plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
492a6e0b375SMatthew G. Knepley   DMEnclosureType    encIn, encAux;
493a6e0b375SMatthew G. Knepley   PetscDS            ds = NULL, dsIn = NULL, dsAux = NULL;
494ca3d3a14SMatthew G. Knepley   Vec                localA = NULL, tv;
49547923291SMatthew G. Knepley   PetscSection       section;
4968c6c5593SMatthew G. Knepley   PetscDualSpace    *sp, *cellsp;
497ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
4988c6c5593SMatthew G. Knepley   PetscInt          *Nc;
499a6e0b375SMatthew G. Knepley   PetscInt           dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfIn, NfAux = 0, f;
500*27f02ce8SMatthew G. Knepley   PetscBool         *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, isHybrid = PETSC_FALSE, transform;
501c330f8ffSToby Isaac   DMField            coordField;
502c330f8ffSToby Isaac   DMLabel            depthLabel;
503c330f8ffSToby Isaac   PetscQuadrature    allPoints = NULL;
50447923291SMatthew G. Knepley   PetscErrorCode     ierr;
50547923291SMatthew G. Knepley 
50647923291SMatthew G. Knepley   PetscFunctionBegin;
507a6e0b375SMatthew G. Knepley   if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);}
508a6e0b375SMatthew G. Knepley   else        {dmIn = dm;}
50988aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
51088aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
511a6e0b375SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
512a6e0b375SMatthew G. Knepley   ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr);
513a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr);
514a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
5158c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
516a6e0b375SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr);
517ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
518ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
519ca3d3a14SMatthew G. Knepley   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
5200de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
521a6e0b375SMatthew G. Knepley   if (dmAux) {
522a6e0b375SMatthew G. Knepley     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
523a6e0b375SMatthew G. Knepley     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
524a6e0b375SMatthew G. Knepley       if (!localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
525a6e0b375SMatthew G. Knepley       if (!minHeight) {
52688aca1feSMatthew G. Knepley         DMLabel spmap;
52788aca1feSMatthew G. Knepley 
52888aca1feSMatthew G. Knepley         /* If dmAux is a surface, then force the projection to take place over a surface */
529a6e0b375SMatthew G. Knepley         ierr = DMPlexGetSubpointMap(plexAux, &spmap);CHKERRQ(ierr);
530e39fcb4eSStefano Zampini         if (spmap) {
531a6e0b375SMatthew G. Knepley           ierr = DMPlexGetVTKCellHeight(plexAux, &minHeight);CHKERRQ(ierr);
532e39fcb4eSStefano Zampini           auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
533e39fcb4eSStefano Zampini         }
53488aca1feSMatthew G. Knepley       }
5350de51b56SMatthew G. Knepley     }
536a6e0b375SMatthew G. Knepley   }
537a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
538a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
539a6e0b375SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr);
5402af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
541e39fcb4eSStefano Zampini   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
542a6e0b375SMatthew G. Knepley   ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr);
543a6e0b375SMatthew G. Knepley   if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);}
544a6e0b375SMatthew G. Knepley   ierr = DMGetFirstLabelEntry_Private(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr);
545a6e0b375SMatthew G. Knepley   if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);}
546a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
547a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
548*27f02ce8SMatthew G. Knepley   ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr);
5498c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
55092fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
5510c898c18SMatthew G. Knepley   if (dmAux) {
552a6e0b375SMatthew G. Knepley     ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr);
553a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
5540c898c18SMatthew G. Knepley   }
555a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
556496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
5578c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
5588c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
559a6e0b375SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
5608c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
56147923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
56247923291SMatthew G. Knepley     PetscObject  obj;
56347923291SMatthew G. Knepley     PetscClassId id;
56447923291SMatthew G. Knepley 
565a6e0b375SMatthew G. Knepley     ierr = PetscDSGetDiscretization(ds, f, &obj);CHKERRQ(ierr);
56647923291SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
56747923291SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
56847923291SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
56947923291SMatthew G. Knepley 
57047923291SMatthew G. Knepley       hasFE   = PETSC_TRUE;
57147923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
57247923291SMatthew G. Knepley       ierr  = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
57347923291SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
57447923291SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
57547923291SMatthew G. Knepley 
57647923291SMatthew G. Knepley       hasFV   = PETSC_TRUE;
57747923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
57847923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
57947923291SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
58047923291SMatthew G. Knepley   }
581c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
582ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
58374fc349aSMatthew G. Knepley     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
58474fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
585a6e0b375SMatthew G. Knepley     PetscBool        isfe;
5864a3e9fdbSToby Isaac     const PetscReal *points;
5874a3e9fdbSToby Isaac     PetscInt         numPoints;
5880c898c18SMatthew G. Knepley 
5892af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
59054ee0cceSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
59154ee0cceSMatthew G. Knepley       if (!effectiveHeight) {sp[f] = cellsp[f];}
59254ee0cceSMatthew G. Knepley       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
59354ee0cceSMatthew G. Knepley     }
59454ee0cceSMatthew G. Knepley     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
595c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
596ef0bb6c7SMatthew G. Knepley     ierr = PetscMalloc2(NfIn, &T, NfAux, &TAux);CHKERRQ(ierr);
597a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
598a6e0b375SMatthew G. Knepley       ierr = PetscDSIsFE_Internal(dsIn, f, &isfe);CHKERRQ(ierr);
599a6e0b375SMatthew G. Knepley       if (!isfe) continue;
600a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
60174fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
60274fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
603ef0bb6c7SMatthew G. Knepley       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]);CHKERRQ(ierr);
6040c898c18SMatthew G. Knepley     }
6050c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
606a6e0b375SMatthew G. Knepley       ierr = PetscDSIsFE_Internal(dsAux, f, &isfe);CHKERRQ(ierr);
607a6e0b375SMatthew G. Knepley       if (!isfe) continue;
608a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
60974fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
61074fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
611ef0bb6c7SMatthew G. Knepley       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]);CHKERRQ(ierr);
6120c898c18SMatthew G. Knepley     }
6130c898c18SMatthew G. Knepley   }
61447923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
6152af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
61688aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
617a6e0b375SMatthew G. Knepley     PetscDS      dsEff         = ds;
6188c6c5593SMatthew G. Knepley     PetscScalar *values;
619b7260050SToby Isaac     PetscBool   *fieldActive;
620b7260050SToby Isaac     PetscInt     maxDegree;
621e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
622c330f8ffSToby Isaac     IS           heightIS;
6238c6c5593SMatthew G. Knepley 
624a6e0b375SMatthew G. Knepley     /* Note we assume that dm and dmIn share the same topology */
625412e9a14SMatthew G. Knepley     ierr = DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd);CHKERRQ(ierr);
626a6e0b375SMatthew G. Knepley     ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
627c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
628c330f8ffSToby Isaac     if (pEnd <= pStart) {
629c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
630c330f8ffSToby Isaac       continue;
631c330f8ffSToby Isaac     }
63247923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
63347923291SMatthew G. Knepley     totDim = 0;
63447923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
635bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
63647923291SMatthew G. Knepley         sp[f] = cellsp[f];
63747923291SMatthew G. Knepley       } else {
638bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
63947923291SMatthew G. Knepley         if (!sp[f]) continue;
64047923291SMatthew G. Knepley       }
64147923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
6429c3cf19fSMatthew G. Knepley       totDim += spDim;
643*27f02ce8SMatthew G. Knepley       if (isHybrid && (f < Nf-1)) totDim += spDim;
64447923291SMatthew G. Knepley     }
645a6e0b375SMatthew G. Knepley     ierr = DMPlexVecGetClosure(plex, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr);
6468c6c5593SMatthew 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);
647c330f8ffSToby Isaac     if (!totDim) {
648c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
649c330f8ffSToby Isaac       continue;
650c330f8ffSToby Isaac     }
651a6e0b375SMatthew G. Knepley     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(ds, effectiveHeight, &dsEff);CHKERRQ(ierr);}
65247923291SMatthew G. Knepley     /* Loop over points at this height */
65369291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
65469291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
6558c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
6568c6c5593SMatthew G. Knepley     if (label) {
6578c6c5593SMatthew G. Knepley       PetscInt i;
65847923291SMatthew G. Knepley 
65947923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
660c330f8ffSToby Isaac         IS              pointIS, isectIS;
66147923291SMatthew G. Knepley         const PetscInt *points;
6628c6c5593SMatthew G. Knepley         PetscInt        n;
663c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
664c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
66547923291SMatthew G. Knepley 
66647923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
66747923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
668c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
669c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
670c330f8ffSToby Isaac         if (!isectIS) continue;
671c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
672c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
673b7260050SToby Isaac         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
674b7260050SToby Isaac         if (maxDegree <= 1) {
675c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
676c330f8ffSToby Isaac         }
677c330f8ffSToby Isaac         if (!quad) {
678c330f8ffSToby Isaac           if (!h && allPoints) {
679c330f8ffSToby Isaac             quad = allPoints;
680c330f8ffSToby Isaac             allPoints = NULL;
681c330f8ffSToby Isaac           } else {
682*27f02ce8SMatthew G. Knepley             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,isHybrid ? dim-h-1 : dim-h,funcs,&quad);CHKERRQ(ierr);
683c330f8ffSToby Isaac           }
684c330f8ffSToby Isaac         }
685a6e0b375SMatthew G. Knepley         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
68647923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
68747923291SMatthew G. Knepley           const PetscInt  point = points[p];
68847923291SMatthew G. Knepley 
689580bdb30SBarry Smith           ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
690c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
691ef0bb6c7SMatthew G. Knepley           ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values);
69247923291SMatthew G. Knepley           if (ierr) {
69347923291SMatthew G. Knepley             PetscErrorCode ierr2;
69469291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
69569291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
69647923291SMatthew G. Knepley             CHKERRQ(ierr);
69747923291SMatthew G. Knepley           }
698a6e0b375SMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
6995f790a90SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode);CHKERRQ(ierr);
70047923291SMatthew G. Knepley         }
701c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
702c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
703c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
704c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
705c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
70647923291SMatthew G. Knepley       }
7078c6c5593SMatthew G. Knepley     } else {
708c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
709c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
710c330f8ffSToby Isaac       IS              pointIS;
711c330f8ffSToby Isaac 
712c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
713b7260050SToby Isaac       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
714b7260050SToby Isaac       if (maxDegree <= 1) {
715c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
716c330f8ffSToby Isaac       }
717c330f8ffSToby Isaac       if (!quad) {
718c330f8ffSToby Isaac         if (!h && allPoints) {
719c330f8ffSToby Isaac           quad = allPoints;
720c330f8ffSToby Isaac           allPoints = NULL;
721c330f8ffSToby Isaac         } else {
722c330f8ffSToby Isaac           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
723c330f8ffSToby Isaac         }
724c330f8ffSToby Isaac       }
725a6e0b375SMatthew G. Knepley       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
7268c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
727580bdb30SBarry Smith         ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
728c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
729ef0bb6c7SMatthew G. Knepley         ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values);
7308c6c5593SMatthew G. Knepley         if (ierr) {
7318c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
73269291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
73369291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
7348c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
7358c6c5593SMatthew G. Knepley         }
736a6e0b375SMatthew G. Knepley         if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
7375f790a90SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode);CHKERRQ(ierr);
7388c6c5593SMatthew G. Knepley       }
739c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
740c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
741c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
742c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
7438c6c5593SMatthew G. Knepley     }
744c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
74569291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
74669291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
74747923291SMatthew G. Knepley   }
7488c6c5593SMatthew G. Knepley   /* Cleanup */
749ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
75074fc349aSMatthew G. Knepley     PetscInt  effectiveHeight = auxBd ? minHeight : 0;
75174fc349aSMatthew G. Knepley     PetscFE   fem, subfem;
752a6e0b375SMatthew G. Knepley     PetscBool isfe;
7530c898c18SMatthew G. Knepley 
754a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
755a6e0b375SMatthew G. Knepley       ierr = PetscDSIsFE_Internal(dsIn, f, &isfe);CHKERRQ(ierr);
756a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
75774fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
75874fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
759ef0bb6c7SMatthew G. Knepley       ierr = PetscTabulationDestroy(&T[f]);CHKERRQ(ierr);
7600c898c18SMatthew G. Knepley     }
7610c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
762a6e0b375SMatthew G. Knepley       ierr = PetscDSIsFE_Internal(dsAux, f, &isfe);CHKERRQ(ierr);
763a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
76474fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
76574fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
766ef0bb6c7SMatthew G. Knepley       ierr = PetscTabulationDestroy(&TAux[f]);CHKERRQ(ierr);
7670c898c18SMatthew G. Knepley     }
768ef0bb6c7SMatthew G. Knepley     ierr = PetscFree2(T, TAux);CHKERRQ(ierr);
7690c898c18SMatthew G. Knepley   }
770c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
771496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
7728c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
773a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
774a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plexIn);CHKERRQ(ierr);
775a6e0b375SMatthew G. Knepley   if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);}
7768c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
77747923291SMatthew G. Knepley }
7788c6c5593SMatthew G. Knepley 
7798c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
7808c6c5593SMatthew G. Knepley {
7818c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
7828c6c5593SMatthew G. Knepley 
7838c6c5593SMatthew G. Knepley   PetscFunctionBegin;
7841c531cf8SMatthew 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);
7858c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
7868c6c5593SMatthew G. Knepley }
7878c6c5593SMatthew G. Knepley 
7881c531cf8SMatthew 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)
7898c6c5593SMatthew G. Knepley {
7908c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
7918c6c5593SMatthew G. Knepley 
7928c6c5593SMatthew G. Knepley   PetscFunctionBegin;
7931c531cf8SMatthew 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);
79447923291SMatthew G. Knepley   PetscFunctionReturn(0);
79547923291SMatthew G. Knepley }
79647923291SMatthew G. Knepley 
7978c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
79847923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
79947923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
80047923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
801191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
80247923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
80347923291SMatthew G. Knepley {
80447923291SMatthew G. Knepley   PetscErrorCode ierr;
80547923291SMatthew G. Knepley 
80647923291SMatthew G. Knepley   PetscFunctionBegin;
8071c531cf8SMatthew 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);
8088c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
80947923291SMatthew G. Knepley }
81047923291SMatthew G. Knepley 
8111c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
8128c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
8138c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
8148c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
815191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
8168c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
8178c6c5593SMatthew G. Knepley {
8188c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
81947923291SMatthew G. Knepley 
8208c6c5593SMatthew G. Knepley   PetscFunctionBegin;
8211c531cf8SMatthew 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);
82247923291SMatthew G. Knepley   PetscFunctionReturn(0);
82347923291SMatthew G. Knepley }
824ece3a9fcSMatthew G. Knepley 
825ece3a9fcSMatthew G. Knepley PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
826ece3a9fcSMatthew G. Knepley                                                void (**funcs)(PetscInt, PetscInt, PetscInt,
827ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
828ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
829ece3a9fcSMatthew G. Knepley                                                               PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
830ece3a9fcSMatthew G. Knepley                                                InsertMode mode, Vec localX)
831ece3a9fcSMatthew G. Knepley {
832ece3a9fcSMatthew G. Knepley   PetscErrorCode ierr;
833ece3a9fcSMatthew G. Knepley 
834ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
835ece3a9fcSMatthew G. Knepley   ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr);
836ece3a9fcSMatthew G. Knepley   PetscFunctionReturn(0);
837ece3a9fcSMatthew G. Knepley }
838