xref: /petsc/src/dm/impls/plex/plexproject.c (revision 083401c691dfef2393fb5fbd2434cf41374a5bcb)
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;
3327f02ce8SMatthew 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);
4127f02ce8SMatthew 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) {
67665f567fSMatthew G. Knepley             const PetscReal *refpoint = &points[q*dim];
68665f567fSMatthew G. Knepley             PetscReal        injpoint[3] = {0., 0., 0.};
69665f567fSMatthew G. Knepley 
70665f567fSMatthew G. Knepley             if (dim != fegeom->dim) {
71665f567fSMatthew G. Knepley               if (isHybrid) {
72665f567fSMatthew G. Knepley                 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
73665f567fSMatthew G. Knepley                 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
74665f567fSMatthew G. Knepley                 refpoint = injpoint;
75665f567fSMatthew G. Knepley               } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %D != %D dual basis spatial dimension", fegeom->dim, dim);
76665f567fSMatthew G. Knepley             }
77665f567fSMatthew G. Knepley             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x);
78c330f8ffSToby Isaac             v0 = x;
798c6c5593SMatthew G. Knepley           } else {
80c330f8ffSToby Isaac             v0 = &fegeom->v[tp*coordDim];
818c6c5593SMatthew G. Knepley           }
82a6e0b375SMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;}
83c330f8ffSToby Isaac           ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr);
84c330f8ffSToby Isaac         }
854bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
862edcad52SToby Isaac         ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr);
87c330f8ffSToby Isaac         ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
88ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr);
89ca3d3a14SMatthew G. Knepley         ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
90c330f8ffSToby Isaac         v += spDim;
9127f02ce8SMatthew G. Knepley         if (isHybrid && (f < Nf-1)) {
9227f02ce8SMatthew G. Knepley           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
9327f02ce8SMatthew G. Knepley         }
94c330f8ffSToby Isaac       } else {
95c330f8ffSToby Isaac         for (d = 0; d < spDim; ++d, ++v) {
96c330f8ffSToby Isaac           ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
97c330f8ffSToby Isaac         }
98c330f8ffSToby Isaac       }
99c330f8ffSToby Isaac     } else {
10027f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
10127f02ce8SMatthew G. Knepley       if (isHybrid && (f < Nf-1)) {
10227f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
10327f02ce8SMatthew G. Knepley       }
1048c6c5593SMatthew G. Knepley     }
1059c3cf19fSMatthew G. Knepley   }
1068c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1078c6c5593SMatthew G. Knepley }
1088c6c5593SMatthew G. Knepley 
109a6e0b375SMatthew G. Knepley /*
110a6e0b375SMatthew G. Knepley   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
111a6e0b375SMatthew G. Knepley 
112a6e0b375SMatthew G. Knepley   Input Parameters:
113a6e0b375SMatthew G. Knepley + dm             - The output DM
114a6e0b375SMatthew G. Knepley . ds             - The output DS
115a6e0b375SMatthew G. Knepley . dmIn           - The input DM
116a6e0b375SMatthew G. Knepley . dsIn           - The input DS
117a6e0b375SMatthew G. Knepley . dmAux          - The auxiliary DM, which is always for the input space
118a6e0b375SMatthew G. Knepley . dsAux          - The auxiliary DS, which is always for the input space
119a6e0b375SMatthew G. Knepley . time           - The time for this evaluation
120a6e0b375SMatthew G. Knepley . localU         - The local solution
121a6e0b375SMatthew G. Knepley . localA         - The local auziliary fields
122a6e0b375SMatthew G. Knepley . cgeom          - The FE geometry for this point
123a6e0b375SMatthew G. Knepley . sp             - The output PetscDualSpace for each field
124a6e0b375SMatthew G. Knepley . p              - The point in the output DM
125ef0bb6c7SMatthew G. Knepley . T              - Input basis and derviatives for each field tabulated on the quadrature points
126ef0bb6c7SMatthew G. Knepley . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
127a6e0b375SMatthew G. Knepley . funcs          - The evaluation function for each field
128a6e0b375SMatthew G. Knepley - ctxs           - The user context for each field
129a6e0b375SMatthew G. Knepley 
130a6e0b375SMatthew G. Knepley   Output Parameter:
131a6e0b375SMatthew G. Knepley . values         - The value for each dual basis vector in the output dual space
132a6e0b375SMatthew G. Knepley 
133a6e0b375SMatthew G. Knepley   Note: Not supported for FV
134a6e0b375SMatthew G. Knepley 
135a6e0b375SMatthew G. Knepley   Level: developer
136a6e0b375SMatthew G. Knepley 
137a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private()
138a6e0b375SMatthew G. Knepley */
139a6e0b375SMatthew 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,
140ef0bb6c7SMatthew G. Knepley                                                    PetscTabulation *T, PetscTabulation *TAux,
1418c6c5593SMatthew G. Knepley                                                    void (**funcs)(PetscInt, PetscInt, PetscInt,
1428c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1438c6c5593SMatthew G. Knepley                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
144191494d9SMatthew G. Knepley                                                                   PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
1458c6c5593SMatthew G. Knepley                                                    PetscScalar values[])
1468c6c5593SMatthew G. Knepley {
1478c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
1484bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
1498c6c5593SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
1508c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
151191494d9SMatthew G. Knepley   const PetscScalar *constants;
1528c6c5593SMatthew G. Knepley   PetscReal         *x;
153ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
1544bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
1554bee2e38SMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed;
156ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
15727f02ce8SMatthew G. Knepley   PetscBool          isAffine, isHybrid, transform;
1588c6c5593SMatthew G. Knepley   PetscErrorCode     ierr;
1598c6c5593SMatthew G. Knepley 
1608c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
161a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
162a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
16327f02ce8SMatthew G. Knepley   ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr);
164a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
165a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(dsIn, &uOff);CHKERRQ(ierr);
166a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x);CHKERRQ(ierr);
167a6e0b375SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
168a6e0b375SMatthew G. Knepley   ierr = PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
169a6e0b375SMatthew G. Knepley   ierr = PetscDSGetConstants(dsIn, &numConstants, &constants);CHKERRQ(ierr);
170a6e0b375SMatthew G. Knepley   ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr);
171a6e0b375SMatthew G. Knepley   ierr = DMGetLocalSection(dmIn, &section);CHKERRQ(ierr);
172a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp);CHKERRQ(ierr);
173a6e0b375SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
1748c6c5593SMatthew G. Knepley   if (dmAux) {
17544171101SMatthew G. Knepley     PetscInt subp;
1761ba34526SMatthew G. Knepley 
177a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
178a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
17992fd8e1eSJed Brown     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
180a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
181a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
182a6e0b375SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
1831ba34526SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
1848c6c5593SMatthew G. Knepley   }
1858c6c5593SMatthew G. Knepley   /* Get values for closure */
1864bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
18727f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
18827f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
1894bee2e38SMatthew G. Knepley   if (isAffine) {
1904bee2e38SMatthew G. Knepley     fegeom.v    = x;
1914bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
1924bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
1934bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
1944bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
1954bee2e38SMatthew G. Knepley   }
1968c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
197c330f8ffSToby Isaac     PetscQuadrature   allPoints;
198c330f8ffSToby Isaac     PetscInt          q, dim, numPoints;
199c330f8ffSToby Isaac     const PetscReal   *points;
200c330f8ffSToby Isaac     PetscScalar       *pointEval;
201c330f8ffSToby Isaac     DM                dm;
202c330f8ffSToby Isaac 
2038c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2048c6c5593SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
205c330f8ffSToby Isaac     if (!funcs[f]) {
206be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
20727f02ce8SMatthew G. Knepley       if (isHybrid && (f < Nf-1)) {
20827f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
20927f02ce8SMatthew G. Knepley       }
210c330f8ffSToby Isaac       continue;
211c330f8ffSToby Isaac     }
212c330f8ffSToby Isaac     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
213b4457527SToby Isaac     ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
214c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
215c330f8ffSToby Isaac     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
2160c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
217c330f8ffSToby Isaac       if (isAffine) {
218ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x);
2191c531cf8SMatthew G. Knepley       } else {
2204bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp*dE];
2214bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp*dE*dE];
2224bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp*dE*dE];
2234bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
2248c6c5593SMatthew G. Knepley       }
225ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
226ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
227a6e0b375SMatthew G. Knepley       if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);}
228a6e0b375SMatthew 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]);
2291c531cf8SMatthew G. Knepley     }
230c330f8ffSToby Isaac     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
231c330f8ffSToby Isaac     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
232c330f8ffSToby Isaac     v += spDim;
23327f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
23427f02ce8SMatthew G. Knepley     if (isHybrid && (f < Nf-1)) {
23527f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
23627f02ce8SMatthew G. Knepley     }
2378c6c5593SMatthew G. Knepley   }
238a6e0b375SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr);
2398c6c5593SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
2408c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2418c6c5593SMatthew G. Knepley }
2428c6c5593SMatthew G. Knepley 
243a6e0b375SMatthew 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,
244ef0bb6c7SMatthew G. Knepley                                                      PetscTabulation *T, PetscTabulation *TAux,
245b18799e0SMatthew G. Knepley                                                      void (**funcs)(PetscInt, PetscInt, PetscInt,
246b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
247b18799e0SMatthew G. Knepley                                                                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
248b18799e0SMatthew G. Knepley                                                                     PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs,
249b18799e0SMatthew G. Knepley                                                      PetscScalar values[])
250b18799e0SMatthew G. Knepley {
251b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
2524bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
253b18799e0SMatthew G. Knepley   PetscScalar       *coefficients   = NULL, *coefficientsAux   = NULL;
254b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
255b18799e0SMatthew G. Knepley   const PetscScalar *constants;
256b18799e0SMatthew G. Knepley   PetscReal         *x;
257ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
2589f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
2594bee2e38SMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed;
260ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
261b18799e0SMatthew G. Knepley   PetscBool          isAffine;
262b18799e0SMatthew G. Knepley   PetscErrorCode     ierr;
263b18799e0SMatthew G. Knepley 
264b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
265a6e0b375SMatthew G. Knepley   if (dm != dmIn) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
266a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
267a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
268a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
269a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
270a6e0b375SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr);
271a6e0b375SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
272a6e0b375SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
27392fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
274a6e0b375SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
275b18799e0SMatthew G. Knepley   if (dmAux) {
276a6e0b375SMatthew G. Knepley     PetscInt subp;
277b18799e0SMatthew G. Knepley 
278a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr);
279a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
28092fd8e1eSJed Brown     ierr = DMGetLocalSection(dmAux, &sectionAux);CHKERRQ(ierr);
281a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
282a6e0b375SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
283a6e0b375SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr);
284b18799e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr);
285b18799e0SMatthew G. Knepley   }
286b18799e0SMatthew G. Knepley   /* Get values for closure */
2874bee2e38SMatthew G. Knepley   isAffine = fgeom->isAffine;
2884dcf62a8SSatish Balay   fegeom.n  = 0;
2894dcf62a8SSatish Balay   fegeom.J  = 0;
2904dcf62a8SSatish Balay   fegeom.v  = 0;
2914dcf62a8SSatish Balay   fegeom.xi = 0;
292a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
293a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
2944bee2e38SMatthew G. Knepley   if (isAffine) {
2954bee2e38SMatthew G. Knepley     fegeom.v    = x;
2964bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
2974bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
2984bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
2994bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
3004bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
3019f209ee4SMatthew G. Knepley 
3029f209ee4SMatthew G. Knepley     cgeom.J     = fgeom->suppJ[0];
3039f209ee4SMatthew G. Knepley     cgeom.invJ  = fgeom->suppInvJ[0];
3049f209ee4SMatthew G. Knepley     cgeom.detJ  = fgeom->suppDetJ[0];
3054bee2e38SMatthew G. Knepley   }
306b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
307b18799e0SMatthew G. Knepley     PetscQuadrature   allPoints;
308b18799e0SMatthew G. Knepley     PetscInt          q, dim, numPoints;
309b18799e0SMatthew G. Knepley     const PetscReal   *points;
310b18799e0SMatthew G. Knepley     PetscScalar       *pointEval;
311b18799e0SMatthew G. Knepley     DM                dm;
312b18799e0SMatthew G. Knepley 
313b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
314b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
315b18799e0SMatthew G. Knepley     if (!funcs[f]) {
316b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
317b18799e0SMatthew G. Knepley       continue;
318b18799e0SMatthew G. Knepley     }
319b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr);
320b4457527SToby Isaac     ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr);
321b18799e0SMatthew G. Knepley     ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
322b18799e0SMatthew G. Knepley     ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
323b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
324b18799e0SMatthew G. Knepley       if (isAffine) {
325ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x);
326b18799e0SMatthew G. Knepley       } else {
3273fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp*dE];
3289f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp*dE*dE];
3299f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp*dE*dE];
3304bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
3314bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp*dE];
3329f209ee4SMatthew G. Knepley 
3339f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][tp*dE*dE];
3349f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][tp*dE*dE];
3359f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][tp];
336b18799e0SMatthew G. Knepley       }
337a6e0b375SMatthew 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 */
338ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr);
339ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);}
3404bee2e38SMatthew 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]);
341b18799e0SMatthew G. Knepley     }
342b18799e0SMatthew G. Knepley     ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr);
343b18799e0SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr);
344b18799e0SMatthew G. Knepley     v += spDim;
345b18799e0SMatthew G. Knepley   }
346a6e0b375SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr);
347b18799e0SMatthew G. Knepley   if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);}
348b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
349b18799e0SMatthew G. Knepley }
350b18799e0SMatthew G. Knepley 
351a6e0b375SMatthew 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[],
352ef0bb6c7SMatthew G. Knepley                                              PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux,
3538c6c5593SMatthew G. Knepley                                              DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
3548c6c5593SMatthew G. Knepley {
3558c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
3568c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
3578c6c5593SMatthew G. Knepley   PetscErrorCode  ierr;
3588c6c5593SMatthew G. Knepley 
3598c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
3608c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
3618c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
3628c6c5593SMatthew G. Knepley   if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);}
3638c6c5593SMatthew G. Knepley   switch (type) {
3648c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
3658c6c5593SMatthew G. Knepley   case DM_BC_NATURAL:
366a6e0b375SMatthew 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;
3678c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
3688c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
369ef0bb6c7SMatthew G. Knepley     ierr = DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
3700c898c18SMatthew G. Knepley                                         (void (**)(PetscInt, PetscInt, PetscInt,
3710c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
3720c898c18SMatthew G. Knepley                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
373191494d9SMatthew G. Knepley                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
374b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
375ef0bb6c7SMatthew G. Knepley     ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux,
376b18799e0SMatthew G. Knepley                                           (void (**)(PetscInt, PetscInt, PetscInt,
377b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
378b18799e0SMatthew G. Knepley                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
379b18799e0SMatthew G. Knepley                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break;
3808c6c5593SMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type);
3818c6c5593SMatthew G. Knepley   }
3828c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
3838c6c5593SMatthew G. Knepley }
3848c6c5593SMatthew G. Knepley 
385df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
386df5c1128SToby Isaac {
387df5c1128SToby Isaac   PetscReal      *points;
388df5c1128SToby Isaac   PetscInt       f, numPoints;
389df5c1128SToby Isaac   PetscErrorCode ierr;
390df5c1128SToby Isaac 
391df5c1128SToby Isaac   PetscFunctionBegin;
392df5c1128SToby Isaac   numPoints = 0;
393df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
394df5c1128SToby Isaac     if (funcs[f]) {
395df5c1128SToby Isaac       PetscQuadrature fAllPoints;
396df5c1128SToby Isaac       PetscInt        fNumPoints;
397df5c1128SToby Isaac 
398b4457527SToby Isaac       ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr);
399df5c1128SToby Isaac       ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr);
400df5c1128SToby Isaac       numPoints += fNumPoints;
401df5c1128SToby Isaac     }
402df5c1128SToby Isaac   }
403df5c1128SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
404df5c1128SToby Isaac   numPoints = 0;
405df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
406df5c1128SToby Isaac     if (funcs[f]) {
407df5c1128SToby Isaac       PetscQuadrature fAllPoints;
40854ee0cceSMatthew G. Knepley       PetscInt        qdim, fNumPoints, q;
409df5c1128SToby Isaac       const PetscReal *fPoints;
410df5c1128SToby Isaac 
411b4457527SToby Isaac       ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr);
41254ee0cceSMatthew G. Knepley       ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr);
41354ee0cceSMatthew 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);
414df5c1128SToby Isaac       for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q];
415df5c1128SToby Isaac       numPoints += fNumPoints;
416df5c1128SToby Isaac     }
417df5c1128SToby Isaac   }
418df5c1128SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
419df5c1128SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
420df5c1128SToby Isaac   PetscFunctionReturn(0);
421df5c1128SToby Isaac }
422df5c1128SToby Isaac 
423a6e0b375SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *ds)
424e5e52638SMatthew G. Knepley {
425a6e0b375SMatthew G. Knepley   DM              plex;
426a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
427e5e52638SMatthew G. Knepley   DMLabel         depthLabel;
428e5e52638SMatthew G. Knepley   PetscInt        dim, cdepth, ls = -1, i;
429e5e52638SMatthew G. Knepley   PetscErrorCode  ierr;
430e5e52638SMatthew G. Knepley 
431e5e52638SMatthew G. Knepley   PetscFunctionBegin;
432e5e52638SMatthew G. Knepley   if (lStart) *lStart = -1;
433e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
434a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dm, odm, &enc);CHKERRQ(ierr);
435e5e52638SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
436a6e0b375SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
437a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
438e5e52638SMatthew G. Knepley   cdepth = dim - height;
439e5e52638SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
440e5e52638SMatthew G. Knepley     IS              pointIS;
441e5e52638SMatthew G. Knepley     const PetscInt *points;
442a6e0b375SMatthew G. Knepley     PetscInt        pdepth, point;
443e5e52638SMatthew G. Knepley 
444e5e52638SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
445e5e52638SMatthew G. Knepley     if (!pointIS) continue; /* No points with that id on this process */
446e5e52638SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
447a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosurePoint(dm, odm, enc, points[0], &point);CHKERRQ(ierr);
448a6e0b375SMatthew G. Knepley     ierr = DMLabelGetValue(depthLabel, point, &pdepth);CHKERRQ(ierr);
449e5e52638SMatthew G. Knepley     if (pdepth == cdepth) {
450a6e0b375SMatthew G. Knepley       ls = point;
451a6e0b375SMatthew G. Knepley       if (ds) {ierr = DMGetCellDS(dm, ls, ds);CHKERRQ(ierr);}
452e5e52638SMatthew G. Knepley     }
453e5e52638SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
454e5e52638SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
455e5e52638SMatthew G. Knepley     if (ls >= 0) break;
456e5e52638SMatthew G. Knepley   }
457e5e52638SMatthew G. Knepley   if (lStart) *lStart = ls;
458a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
459e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
460e5e52638SMatthew G. Knepley }
461e5e52638SMatthew G. Knepley 
4620de51b56SMatthew G. Knepley /*
4630de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
4640de51b56SMatthew G. Knepley 
4650de51b56SMatthew G. Knepley   There are several different scenarios:
4660de51b56SMatthew G. Knepley 
4670de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
4680de51b56SMatthew G. Knepley 
4690de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
4700de51b56SMatthew G. Knepley 
4710de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
4720de51b56SMatthew G. Knepley 
4730de51b56SMatthew 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.
4740de51b56SMatthew G. Knepley 
4750de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
4760de51b56SMatthew G. Knepley 
4770de51b56SMatthew 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.
4780de51b56SMatthew G. Knepley 
4790de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
4800de51b56SMatthew G. Knepley 
4810de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
4820de51b56SMatthew G. Knepley 
4830de51b56SMatthew 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.
4840de51b56SMatthew 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
4855f790a90SMatthew 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
4860de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
4870de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
4880de51b56SMatthew G. Knepley 
4890de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
4900de51b56SMatthew G. Knepley 
4910de51b56SMatthew 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.
4920de51b56SMatthew G. Knepley */
49346fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU,
4941c531cf8SMatthew G. Knepley                                                   PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[],
4958c6c5593SMatthew G. Knepley                                                   DMBoundaryConditionType type, void (**funcs)(void), void **ctxs,
4968c6c5593SMatthew G. Knepley                                                   InsertMode mode, Vec localX)
4978c6c5593SMatthew G. Knepley {
498a6e0b375SMatthew G. Knepley   DM                 plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
499a6e0b375SMatthew G. Knepley   DMEnclosureType    encIn, encAux;
500a6e0b375SMatthew G. Knepley   PetscDS            ds = NULL, dsIn = NULL, dsAux = NULL;
501ca3d3a14SMatthew G. Knepley   Vec                localA = NULL, tv;
50247923291SMatthew G. Knepley   PetscSection       section;
5038c6c5593SMatthew G. Knepley   PetscDualSpace    *sp, *cellsp;
504ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
5058c6c5593SMatthew G. Knepley   PetscInt          *Nc;
506a6e0b375SMatthew G. Knepley   PetscInt           dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfIn, NfAux = 0, f;
50727f02ce8SMatthew G. Knepley   PetscBool         *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, isHybrid = PETSC_FALSE, transform;
508c330f8ffSToby Isaac   DMField            coordField;
509c330f8ffSToby Isaac   DMLabel            depthLabel;
510c330f8ffSToby Isaac   PetscQuadrature    allPoints = NULL;
51147923291SMatthew G. Knepley   PetscErrorCode     ierr;
51247923291SMatthew G. Knepley 
51347923291SMatthew G. Knepley   PetscFunctionBegin;
514a6e0b375SMatthew G. Knepley   if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);}
515a6e0b375SMatthew G. Knepley   else        {dmIn = dm;}
51688aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
51788aca1feSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr);
518a6e0b375SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
519a6e0b375SMatthew G. Knepley   ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr);
520a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr);
521a6e0b375SMatthew G. Knepley   ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
5228c6c5593SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
523a6e0b375SMatthew G. Knepley   ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr);
524ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr);
525ca3d3a14SMatthew G. Knepley   ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr);
526ca3d3a14SMatthew G. Knepley   ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr);
5270de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
528a6e0b375SMatthew G. Knepley   if (dmAux) {
529a6e0b375SMatthew G. Knepley     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
530a6e0b375SMatthew G. Knepley     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
531a6e0b375SMatthew G. Knepley       if (!localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector");
532a6e0b375SMatthew G. Knepley       if (!minHeight) {
53388aca1feSMatthew G. Knepley         DMLabel spmap;
53488aca1feSMatthew G. Knepley 
53588aca1feSMatthew G. Knepley         /* If dmAux is a surface, then force the projection to take place over a surface */
536a6e0b375SMatthew G. Knepley         ierr = DMPlexGetSubpointMap(plexAux, &spmap);CHKERRQ(ierr);
537e39fcb4eSStefano Zampini         if (spmap) {
538a6e0b375SMatthew G. Knepley           ierr = DMPlexGetVTKCellHeight(plexAux, &minHeight);CHKERRQ(ierr);
539e39fcb4eSStefano Zampini           auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE;
540e39fcb4eSStefano Zampini         }
54188aca1feSMatthew G. Knepley       }
5420de51b56SMatthew G. Knepley     }
543a6e0b375SMatthew G. Knepley   }
544a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
545a6e0b375SMatthew G. Knepley   ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr);
546a6e0b375SMatthew G. Knepley   ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr);
5472af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
548e39fcb4eSStefano Zampini   if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim);
549a6e0b375SMatthew G. Knepley   ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr);
550a6e0b375SMatthew G. Knepley   if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);}
551a6e0b375SMatthew G. Knepley   ierr = DMGetFirstLabelEntry_Private(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr);
552a6e0b375SMatthew G. Knepley   if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);}
553a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
554a6e0b375SMatthew G. Knepley   ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr);
55527f02ce8SMatthew G. Knepley   ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr);
5568c6c5593SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
55792fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
5580c898c18SMatthew G. Knepley   if (dmAux) {
559a6e0b375SMatthew G. Knepley     ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr);
560a6e0b375SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
5610c898c18SMatthew G. Knepley   }
562a6e0b375SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
563496733e2SMatthew G. Knepley   ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr);
5648c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);}
5658c6c5593SMatthew G. Knepley   else               {cellsp = sp;}
566a6e0b375SMatthew G. Knepley   if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);}
5678c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
56847923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
569665f567fSMatthew G. Knepley     PetscDiscType disctype;
57047923291SMatthew G. Knepley 
571665f567fSMatthew G. Knepley     ierr = PetscDSGetDiscType_Internal(ds, f, &disctype);CHKERRQ(ierr);
572665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
573665f567fSMatthew G. Knepley       PetscFE fe;
57447923291SMatthew G. Knepley 
57547923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
576665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
577665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
57847923291SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr);
579665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
580665f567fSMatthew G. Knepley       PetscFV fv;
58147923291SMatthew G. Knepley 
58247923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
583665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
584665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fv);CHKERRQ(ierr);
58547923291SMatthew G. Knepley       ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr);
586665f567fSMatthew G. Knepley     } else {
587665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
588665f567fSMatthew G. Knepley       cellsp[f] = NULL;
589665f567fSMatthew G. Knepley     }
59047923291SMatthew G. Knepley   }
591c330f8ffSToby Isaac   ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr);
592ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
59374fc349aSMatthew G. Knepley     PetscInt         effectiveHeight = auxBd ? minHeight : 0;
59474fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
595665f567fSMatthew G. Knepley     PetscDiscType    disctype;
5964a3e9fdbSToby Isaac     const PetscReal *points;
5974a3e9fdbSToby Isaac     PetscInt         numPoints;
5980c898c18SMatthew G. Knepley 
5992af58022SMatthew G. Knepley     if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
60054ee0cceSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
60154ee0cceSMatthew G. Knepley       if (!effectiveHeight) {sp[f] = cellsp[f];}
60254ee0cceSMatthew G. Knepley       else                  {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);}
60354ee0cceSMatthew G. Knepley     }
60454ee0cceSMatthew G. Knepley     ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr);
605c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr);
606ef0bb6c7SMatthew G. Knepley     ierr = PetscMalloc2(NfIn, &T, NfAux, &TAux);CHKERRQ(ierr);
607a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
608665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscType_Internal(dsIn, f, &disctype);CHKERRQ(ierr);
609665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
610a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
61174fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
61274fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
613ef0bb6c7SMatthew G. Knepley       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]);CHKERRQ(ierr);
6140c898c18SMatthew G. Knepley     }
6150c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
616665f567fSMatthew G. Knepley       ierr = PetscDSGetDiscType_Internal(dsAux, f, &disctype);CHKERRQ(ierr);
617665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
618a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
61974fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
62074fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
621ef0bb6c7SMatthew G. Knepley       ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]);CHKERRQ(ierr);
6220c898c18SMatthew G. Knepley     }
6230c898c18SMatthew G. Knepley   }
62447923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
6252af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
62688aca1feSMatthew G. Knepley     PetscInt     effectiveHeight = h - (auxBd ? 0 : minHeight);
627a6e0b375SMatthew G. Knepley     PetscDS      dsEff         = ds;
6288c6c5593SMatthew G. Knepley     PetscScalar *values;
629b7260050SToby Isaac     PetscBool   *fieldActive;
630b7260050SToby Isaac     PetscInt     maxDegree;
631e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
632c330f8ffSToby Isaac     IS           heightIS;
6338c6c5593SMatthew G. Knepley 
634a6e0b375SMatthew G. Knepley     /* Note we assume that dm and dmIn share the same topology */
635412e9a14SMatthew G. Knepley     ierr = DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd);CHKERRQ(ierr);
636a6e0b375SMatthew G. Knepley     ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr);
637c330f8ffSToby Isaac     ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr);
638c330f8ffSToby Isaac     if (pEnd <= pStart) {
639c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
640c330f8ffSToby Isaac       continue;
641c330f8ffSToby Isaac     }
64247923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
64347923291SMatthew G. Knepley     totDim = 0;
64447923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
645bec1f4aaSMatthew G. Knepley       if (!effectiveHeight) {
64647923291SMatthew G. Knepley         sp[f] = cellsp[f];
64747923291SMatthew G. Knepley       } else {
648bec1f4aaSMatthew G. Knepley         ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);
64947923291SMatthew G. Knepley       }
650665f567fSMatthew G. Knepley       if (!sp[f]) continue;
65147923291SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
6529c3cf19fSMatthew G. Knepley       totDim += spDim;
65327f02ce8SMatthew G. Knepley       if (isHybrid && (f < Nf-1)) totDim += spDim;
65447923291SMatthew G. Knepley     }
655a6e0b375SMatthew G. Knepley     ierr = DMPlexVecGetClosure(plex, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr);
656*083401c6SMatthew G. Knepley     if (numValues != totDim) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point (%D) closure size %D != dual space dimension %D", lStart < 0 ? pStart : lStart, numValues, totDim);
657c330f8ffSToby Isaac     if (!totDim) {
658c330f8ffSToby Isaac       ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
659c330f8ffSToby Isaac       continue;
660c330f8ffSToby Isaac     }
661a6e0b375SMatthew G. Knepley     if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(ds, effectiveHeight, &dsEff);CHKERRQ(ierr);}
66247923291SMatthew G. Knepley     /* Loop over points at this height */
66369291d52SBarry Smith     ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
66469291d52SBarry Smith     ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
6658c6c5593SMatthew G. Knepley     for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
6668c6c5593SMatthew G. Knepley     if (label) {
6678c6c5593SMatthew G. Knepley       PetscInt i;
66847923291SMatthew G. Knepley 
66947923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
670c330f8ffSToby Isaac         IS              pointIS, isectIS;
67147923291SMatthew G. Knepley         const PetscInt *points;
6728c6c5593SMatthew G. Knepley         PetscInt        n;
673c330f8ffSToby Isaac         PetscFEGeom  *fegeom = NULL, *chunkgeom = NULL;
674c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
67547923291SMatthew G. Knepley 
67647923291SMatthew G. Knepley         ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
67747923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
678c330f8ffSToby Isaac         ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr);
679c330f8ffSToby Isaac         ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
680c330f8ffSToby Isaac         if (!isectIS) continue;
681c330f8ffSToby Isaac         ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr);
682c330f8ffSToby Isaac         ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr);
683b7260050SToby Isaac         ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr);
684b7260050SToby Isaac         if (maxDegree <= 1) {
685c330f8ffSToby Isaac           ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr);
686c330f8ffSToby Isaac         }
687c330f8ffSToby Isaac         if (!quad) {
688c330f8ffSToby Isaac           if (!h && allPoints) {
689c330f8ffSToby Isaac             quad = allPoints;
690c330f8ffSToby Isaac             allPoints = NULL;
691c330f8ffSToby Isaac           } else {
69227f02ce8SMatthew G. Knepley             ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,isHybrid ? dim-h-1 : dim-h,funcs,&quad);CHKERRQ(ierr);
693c330f8ffSToby Isaac           }
694c330f8ffSToby Isaac         }
695a6e0b375SMatthew G. Knepley         ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
69647923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
69747923291SMatthew G. Knepley           const PetscInt  point = points[p];
69847923291SMatthew G. Knepley 
699580bdb30SBarry Smith           ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
700c330f8ffSToby Isaac           ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
701ef0bb6c7SMatthew 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);
70247923291SMatthew G. Knepley           if (ierr) {
70347923291SMatthew G. Knepley             PetscErrorCode ierr2;
70469291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
70569291d52SBarry Smith             ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
70647923291SMatthew G. Knepley             CHKERRQ(ierr);
70747923291SMatthew G. Knepley           }
708a6e0b375SMatthew G. Knepley           if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
7095f790a90SMatthew G. Knepley           ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode);CHKERRQ(ierr);
71047923291SMatthew G. Knepley         }
711c330f8ffSToby Isaac         ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr);
712c330f8ffSToby Isaac         ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
713c330f8ffSToby Isaac         ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
714c330f8ffSToby Isaac         ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr);
715c330f8ffSToby Isaac         ierr = ISDestroy(&isectIS);CHKERRQ(ierr);
71647923291SMatthew G. Knepley       }
7178c6c5593SMatthew G. Knepley     } else {
718c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
719c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
720c330f8ffSToby Isaac       IS              pointIS;
721c330f8ffSToby Isaac 
722c330f8ffSToby Isaac       ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr);
723b7260050SToby Isaac       ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr);
724b7260050SToby Isaac       if (maxDegree <= 1) {
725c330f8ffSToby Isaac         ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr);
726c330f8ffSToby Isaac       }
727c330f8ffSToby Isaac       if (!quad) {
728c330f8ffSToby Isaac         if (!h && allPoints) {
729c330f8ffSToby Isaac           quad = allPoints;
730c330f8ffSToby Isaac           allPoints = NULL;
731c330f8ffSToby Isaac         } else {
732c330f8ffSToby Isaac           ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr);
733c330f8ffSToby Isaac         }
734c330f8ffSToby Isaac       }
735a6e0b375SMatthew G. Knepley       ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr);
7368c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
737580bdb30SBarry Smith         ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr);
738c330f8ffSToby Isaac         ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr);
739ef0bb6c7SMatthew 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);
7408c6c5593SMatthew G. Knepley         if (ierr) {
7418c6c5593SMatthew G. Knepley           PetscErrorCode ierr2;
74269291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2);
74369291d52SBarry Smith           ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2);
7448c6c5593SMatthew G. Knepley           CHKERRQ(ierr);
7458c6c5593SMatthew G. Knepley         }
746a6e0b375SMatthew G. Knepley         if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);}
7475f790a90SMatthew G. Knepley         ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode);CHKERRQ(ierr);
7488c6c5593SMatthew G. Knepley       }
749c330f8ffSToby Isaac       ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr);
750c330f8ffSToby Isaac       ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr);
751c330f8ffSToby Isaac       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
752c330f8ffSToby Isaac       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
7538c6c5593SMatthew G. Knepley     }
754c330f8ffSToby Isaac     ierr = ISDestroy(&heightIS);CHKERRQ(ierr);
75569291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
75669291d52SBarry Smith     ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr);
75747923291SMatthew G. Knepley   }
7588c6c5593SMatthew G. Knepley   /* Cleanup */
759ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
76074fc349aSMatthew G. Knepley     PetscInt effectiveHeight = auxBd ? minHeight : 0;
76174fc349aSMatthew G. Knepley     PetscFE  fem, subfem;
7620c898c18SMatthew G. Knepley 
763a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
764a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr);
76574fc349aSMatthew G. Knepley       if (!effectiveHeight) {subfem = fem;}
76674fc349aSMatthew G. Knepley       else                  {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
767ef0bb6c7SMatthew G. Knepley       ierr = PetscTabulationDestroy(&T[f]);CHKERRQ(ierr);
7680c898c18SMatthew G. Knepley     }
7690c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
770a6e0b375SMatthew G. Knepley       ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr);
77174fc349aSMatthew G. Knepley       if (!effectiveHeight || auxBd) {subfem = fem;}
77274fc349aSMatthew G. Knepley       else                           {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);}
773ef0bb6c7SMatthew G. Knepley       ierr = PetscTabulationDestroy(&TAux[f]);CHKERRQ(ierr);
7740c898c18SMatthew G. Knepley     }
775ef0bb6c7SMatthew G. Knepley     ierr = PetscFree2(T, TAux);CHKERRQ(ierr);
7760c898c18SMatthew G. Knepley   }
777c330f8ffSToby Isaac   ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr);
778496733e2SMatthew G. Knepley   ierr = PetscFree2(isFE, sp);CHKERRQ(ierr);
7798c6c5593SMatthew G. Knepley   if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);}
780a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
781a6e0b375SMatthew G. Knepley   ierr = DMDestroy(&plexIn);CHKERRQ(ierr);
782a6e0b375SMatthew G. Knepley   if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);}
7838c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
78447923291SMatthew G. Knepley }
7858c6c5593SMatthew G. Knepley 
7868c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
7878c6c5593SMatthew G. Knepley {
7888c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
7898c6c5593SMatthew G. Knepley 
7908c6c5593SMatthew G. Knepley   PetscFunctionBegin;
7911c531cf8SMatthew 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);
7928c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
7938c6c5593SMatthew G. Knepley }
7948c6c5593SMatthew G. Knepley 
7951c531cf8SMatthew 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)
7968c6c5593SMatthew G. Knepley {
7978c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
7988c6c5593SMatthew G. Knepley 
7998c6c5593SMatthew G. Knepley   PetscFunctionBegin;
8001c531cf8SMatthew 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);
80147923291SMatthew G. Knepley   PetscFunctionReturn(0);
80247923291SMatthew G. Knepley }
80347923291SMatthew G. Knepley 
8048c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU,
80547923291SMatthew G. Knepley                                         void (**funcs)(PetscInt, PetscInt, PetscInt,
80647923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
80747923291SMatthew G. Knepley                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
808191494d9SMatthew G. Knepley                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
80947923291SMatthew G. Knepley                                         InsertMode mode, Vec localX)
81047923291SMatthew G. Knepley {
81147923291SMatthew G. Knepley   PetscErrorCode ierr;
81247923291SMatthew G. Knepley 
81347923291SMatthew G. Knepley   PetscFunctionBegin;
8141c531cf8SMatthew 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);
8158c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
81647923291SMatthew G. Knepley }
81747923291SMatthew G. Knepley 
8181c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
8198c6c5593SMatthew G. Knepley                                              void (**funcs)(PetscInt, PetscInt, PetscInt,
8208c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
8218c6c5593SMatthew G. Knepley                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
822191494d9SMatthew G. Knepley                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
8238c6c5593SMatthew G. Knepley                                              InsertMode mode, Vec localX)
8248c6c5593SMatthew G. Knepley {
8258c6c5593SMatthew G. Knepley   PetscErrorCode ierr;
82647923291SMatthew G. Knepley 
8278c6c5593SMatthew G. Knepley   PetscFunctionBegin;
8281c531cf8SMatthew 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);
82947923291SMatthew G. Knepley   PetscFunctionReturn(0);
83047923291SMatthew G. Knepley }
831ece3a9fcSMatthew G. Knepley 
832ece3a9fcSMatthew G. Knepley PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU,
833ece3a9fcSMatthew G. Knepley                                                void (**funcs)(PetscInt, PetscInt, PetscInt,
834ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
835ece3a9fcSMatthew G. Knepley                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
836ece3a9fcSMatthew G. Knepley                                                               PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
837ece3a9fcSMatthew G. Knepley                                                InsertMode mode, Vec localX)
838ece3a9fcSMatthew G. Knepley {
839ece3a9fcSMatthew G. Knepley   PetscErrorCode ierr;
840ece3a9fcSMatthew G. Knepley 
841ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
842ece3a9fcSMatthew 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);
843ece3a9fcSMatthew G. Knepley   PetscFunctionReturn(0);
844ece3a9fcSMatthew G. Knepley }
845