xref: /petsc/src/dm/impls/plex/plexproject.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
147923291SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
247923291SMatthew G. Knepley 
38c6c5593SMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
48c6c5593SMatthew G. Knepley 
51b32699bSMatthew G. Knepley /*@
61b32699bSMatthew G. Knepley   DMPlexGetActivePoint - Get the point on which projection is currently working
71b32699bSMatthew G. Knepley 
81b32699bSMatthew G. Knepley   Not collective
91b32699bSMatthew G. Knepley 
104165533cSJose E. Roman   Input Parameter:
111b32699bSMatthew G. Knepley . dm   - the DM
121b32699bSMatthew G. Knepley 
134165533cSJose E. Roman   Output Parameter:
141b32699bSMatthew G. Knepley . point - The mesh point involved in the current projection
151b32699bSMatthew G. Knepley 
161b32699bSMatthew G. Knepley   Level: developer
171b32699bSMatthew G. Knepley 
18db781477SPatrick Sanan .seealso: `DMPlexSetActivePoint()`
191b32699bSMatthew G. Knepley @*/
20*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point)
21*d71ae5a4SJacob Faibussowitsch {
221b32699bSMatthew G. Knepley   PetscFunctionBeginHot;
231b32699bSMatthew G. Knepley   *point = ((DM_Plex *)dm->data)->activePoint;
241b32699bSMatthew G. Knepley   PetscFunctionReturn(0);
251b32699bSMatthew G. Knepley }
261b32699bSMatthew G. Knepley 
271b32699bSMatthew G. Knepley /*@
281b32699bSMatthew G. Knepley   DMPlexSetActivePoint - Set the point on which projection is currently working
291b32699bSMatthew G. Knepley 
301b32699bSMatthew G. Knepley   Not collective
311b32699bSMatthew G. Knepley 
324165533cSJose E. Roman   Input Parameters:
331b32699bSMatthew G. Knepley + dm   - the DM
341b32699bSMatthew G. Knepley - point - The mesh point involved in the current projection
351b32699bSMatthew G. Knepley 
361b32699bSMatthew G. Knepley   Level: developer
371b32699bSMatthew G. Knepley 
38db781477SPatrick Sanan .seealso: `DMPlexGetActivePoint()`
391b32699bSMatthew G. Knepley @*/
40*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point)
41*d71ae5a4SJacob Faibussowitsch {
421b32699bSMatthew G. Knepley   PetscFunctionBeginHot;
431b32699bSMatthew G. Knepley   ((DM_Plex *)dm->data)->activePoint = point;
441b32699bSMatthew G. Knepley   PetscFunctionReturn(0);
451b32699bSMatthew G. Knepley }
461b32699bSMatthew G. Knepley 
47a6e0b375SMatthew G. Knepley /*
48a6e0b375SMatthew G. Knepley   DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point
49a6e0b375SMatthew G. Knepley 
50a6e0b375SMatthew G. Knepley   Input Parameters:
51a6e0b375SMatthew G. Knepley + dm     - The output DM
52a6e0b375SMatthew G. Knepley . ds     - The output DS
53a6e0b375SMatthew G. Knepley . dmIn   - The input DM
54a6e0b375SMatthew G. Knepley . dsIn   - The input DS
55a6e0b375SMatthew G. Knepley . time   - The time for this evaluation
56a6e0b375SMatthew G. Knepley . fegeom - The FE geometry for this point
57a6e0b375SMatthew G. Knepley . fvgeom - The FV geometry for this point
58a6e0b375SMatthew G. Knepley . isFE   - Flag indicating whether each output field has an FE discretization
59a6e0b375SMatthew G. Knepley . sp     - The output PetscDualSpace for each field
60a6e0b375SMatthew G. Knepley . funcs  - The evaluation function for each field
61a6e0b375SMatthew G. Knepley - ctxs   - The user context for each field
62a6e0b375SMatthew G. Knepley 
63a6e0b375SMatthew G. Knepley   Output Parameter:
64a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space
65a6e0b375SMatthew G. Knepley 
66a6e0b375SMatthew G. Knepley   Level: developer
67a6e0b375SMatthew G. Knepley 
68db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()`
69a6e0b375SMatthew G. Knepley */
70*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, PetscScalar values[])
71*d71ae5a4SJacob Faibussowitsch {
72a6e0b375SMatthew G. Knepley   PetscInt  coordDim, Nf, *Nc, f, spDim, d, v, tp;
735fedec97SMatthew G. Knepley   PetscBool isAffine, isCohesive, transform;
748c6c5593SMatthew G. Knepley 
758c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
769566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dmIn, &coordDim));
779566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
809566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
818c6c5593SMatthew G. Knepley   /* Get values for closure */
82c330f8ffSToby Isaac   isAffine = fegeom->isAffine;
83c330f8ffSToby Isaac   for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
848c6c5593SMatthew G. Knepley     void *const ctx = ctxs ? ctxs[f] : NULL;
855fedec97SMatthew G. Knepley     PetscBool   cohesive;
868c6c5593SMatthew G. Knepley 
878c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
899566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
908c6c5593SMatthew G. Knepley     if (funcs[f]) {
91c330f8ffSToby Isaac       if (isFE[f]) {
92c330f8ffSToby Isaac         PetscQuadrature  allPoints;
93c330f8ffSToby Isaac         PetscInt         q, dim, numPoints;
94c330f8ffSToby Isaac         const PetscReal *points;
95c330f8ffSToby Isaac         PetscScalar     *pointEval;
96c330f8ffSToby Isaac         PetscReal       *x;
97ca3d3a14SMatthew G. Knepley         DM               rdm;
98c330f8ffSToby Isaac 
999566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDM(sp[f], &rdm));
1009566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
1019566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
1029566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
1039566063dSJacob Faibussowitsch         PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x));
104737e23dcSMatthew G. Knepley         PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f]));
105c330f8ffSToby Isaac         for (q = 0; q < numPoints; q++, tp++) {
106c330f8ffSToby Isaac           const PetscReal *v0;
107c330f8ffSToby Isaac 
108c330f8ffSToby Isaac           if (isAffine) {
109665f567fSMatthew G. Knepley             const PetscReal *refpoint    = &points[q * dim];
110665f567fSMatthew G. Knepley             PetscReal        injpoint[3] = {0., 0., 0.};
111665f567fSMatthew G. Knepley 
112665f567fSMatthew G. Knepley             if (dim != fegeom->dim) {
1135fedec97SMatthew G. Knepley               if (isCohesive) {
114665f567fSMatthew G. Knepley                 /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
115665f567fSMatthew G. Knepley                 for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
116665f567fSMatthew G. Knepley                 refpoint = injpoint;
11763a3b9bcSJacob Faibussowitsch               } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim);
118665f567fSMatthew G. Knepley             }
119665f567fSMatthew G. Knepley             CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x);
120c330f8ffSToby Isaac             v0 = x;
1218c6c5593SMatthew G. Knepley           } else {
122c330f8ffSToby Isaac             v0 = &fegeom->v[tp * coordDim];
1238c6c5593SMatthew G. Knepley           }
1249371c9d4SSatish Balay           if (transform) {
1259371c9d4SSatish Balay             PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx));
1269371c9d4SSatish Balay             v0 = x;
1279371c9d4SSatish Balay           }
1289566063dSJacob Faibussowitsch           PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx));
129c330f8ffSToby Isaac         }
1304bee2e38SMatthew G. Knepley         /* Transform point evaluations pointEval[q,c] */
1319566063dSJacob Faibussowitsch         PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
1329566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
1339566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x));
1349566063dSJacob Faibussowitsch         PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
135c330f8ffSToby Isaac         v += spDim;
1365fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) {
13727f02ce8SMatthew G. Knepley           for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
13827f02ce8SMatthew G. Knepley         }
139c330f8ffSToby Isaac       } else {
14048a46eb9SPierre Jolivet         for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
141c330f8ffSToby Isaac       }
142c330f8ffSToby Isaac     } else {
14327f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
1445fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
14527f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
14627f02ce8SMatthew G. Knepley       }
1478c6c5593SMatthew G. Knepley     }
1489c3cf19fSMatthew G. Knepley   }
1498c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
1508c6c5593SMatthew G. Knepley }
1518c6c5593SMatthew G. Knepley 
152a6e0b375SMatthew G. Knepley /*
153a6e0b375SMatthew G. Knepley   DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
154a6e0b375SMatthew G. Knepley 
155a6e0b375SMatthew G. Knepley   Input Parameters:
156a6e0b375SMatthew G. Knepley + dm             - The output DM
157a6e0b375SMatthew G. Knepley . ds             - The output DS
158a6e0b375SMatthew G. Knepley . dmIn           - The input DM
159a6e0b375SMatthew G. Knepley . dsIn           - The input DS
160a6e0b375SMatthew G. Knepley . dmAux          - The auxiliary DM, which is always for the input space
161a6e0b375SMatthew G. Knepley . dsAux          - The auxiliary DS, which is always for the input space
162a6e0b375SMatthew G. Knepley . time           - The time for this evaluation
163a6e0b375SMatthew G. Knepley . localU         - The local solution
164a6e0b375SMatthew G. Knepley . localA         - The local auziliary fields
165a6e0b375SMatthew G. Knepley . cgeom          - The FE geometry for this point
166a6e0b375SMatthew G. Knepley . sp             - The output PetscDualSpace for each field
167a6e0b375SMatthew G. Knepley . p              - The point in the output DM
168a5b23f4aSJose E. Roman . T              - Input basis and derivatives for each field tabulated on the quadrature points
169ef0bb6c7SMatthew G. Knepley . TAux           - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
170a6e0b375SMatthew G. Knepley . funcs          - The evaluation function for each field
171a6e0b375SMatthew G. Knepley - ctxs           - The user context for each field
172a6e0b375SMatthew G. Knepley 
173a6e0b375SMatthew G. Knepley   Output Parameter:
174a6e0b375SMatthew G. Knepley . values         - The value for each dual basis vector in the output dual space
175a6e0b375SMatthew G. Knepley 
176a6e0b375SMatthew G. Knepley   Note: Not supported for FV
177a6e0b375SMatthew G. Knepley 
178a6e0b375SMatthew G. Knepley   Level: developer
179a6e0b375SMatthew G. Knepley 
180db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()`
181a6e0b375SMatthew G. Knepley */
182*d71ae5a4SJacob Faibussowitsch 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, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[])
183*d71ae5a4SJacob Faibussowitsch {
1848c6c5593SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
1854bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
1868c6c5593SMatthew G. Knepley   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
1878c6c5593SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
188191494d9SMatthew G. Knepley   const PetscScalar *constants;
1898c6c5593SMatthew G. Knepley   PetscReal         *x;
190ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
1914bee2e38SMatthew G. Knepley   PetscFEGeom        fegeom;
1924bee2e38SMatthew G. Knepley   const PetscInt     dE = cgeom->dimEmbed;
193ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
1945fedec97SMatthew G. Knepley   PetscBool          isAffine, isCohesive, transform;
1958c6c5593SMatthew G. Knepley 
1968c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
1979566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1989566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
1999566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
2009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
2019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
2029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
2039566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
2049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
2059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
2069566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dmIn, &transform));
2079566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dmIn, &section));
2089566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
2099566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
2108c6c5593SMatthew G. Knepley   if (dmAux) {
21144171101SMatthew G. Knepley     PetscInt subp;
2121ba34526SMatthew G. Knepley 
2139566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
2149566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2159566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
2169566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2179566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2189566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
2199566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
2208c6c5593SMatthew G. Knepley   }
2218c6c5593SMatthew G. Knepley   /* Get values for closure */
2224bee2e38SMatthew G. Knepley   isAffine        = cgeom->isAffine;
22327f02ce8SMatthew G. Knepley   fegeom.dim      = cgeom->dim;
22427f02ce8SMatthew G. Knepley   fegeom.dimEmbed = cgeom->dimEmbed;
2254bee2e38SMatthew G. Knepley   if (isAffine) {
2264bee2e38SMatthew G. Knepley     fegeom.v    = x;
2274bee2e38SMatthew G. Knepley     fegeom.xi   = cgeom->xi;
2284bee2e38SMatthew G. Knepley     fegeom.J    = cgeom->J;
2294bee2e38SMatthew G. Knepley     fegeom.invJ = cgeom->invJ;
2304bee2e38SMatthew G. Knepley     fegeom.detJ = cgeom->detJ;
2314bee2e38SMatthew G. Knepley   }
2328c6c5593SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
233c330f8ffSToby Isaac     PetscQuadrature  allPoints;
234c330f8ffSToby Isaac     PetscInt         q, dim, numPoints;
235c330f8ffSToby Isaac     const PetscReal *points;
236c330f8ffSToby Isaac     PetscScalar     *pointEval;
2375fedec97SMatthew G. Knepley     PetscBool        cohesive;
238c330f8ffSToby Isaac     DM               dm;
239c330f8ffSToby Isaac 
2408c6c5593SMatthew G. Knepley     if (!sp[f]) continue;
2419566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
2429566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
243c330f8ffSToby Isaac     if (!funcs[f]) {
244be1504a2SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
2455fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) {
24627f02ce8SMatthew G. Knepley         for (d = 0; d < spDim; d++, v++) values[v] = 0.;
24727f02ce8SMatthew G. Knepley       }
248c330f8ffSToby Isaac       continue;
249c330f8ffSToby Isaac     }
2509566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
2519566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
2529566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
2539566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
2540c898c18SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
255c330f8ffSToby Isaac       if (isAffine) {
256ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x);
2571c531cf8SMatthew G. Knepley       } else {
2584bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[tp * dE];
2594bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[tp * dE * dE];
2604bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[tp * dE * dE];
2614bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[tp];
2628c6c5593SMatthew G. Knepley       }
2639566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
2649566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
2659566063dSJacob Faibussowitsch       if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
266a6e0b375SMatthew 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]);
2671c531cf8SMatthew G. Knepley     }
2689566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
2699566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
270c330f8ffSToby Isaac     v += spDim;
27127f02ce8SMatthew G. Knepley     /* TODO: For now, set both sides equal, but this should use info from other support cell */
2725fedec97SMatthew G. Knepley     if (isCohesive && !cohesive) {
27327f02ce8SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
27427f02ce8SMatthew G. Knepley     }
2758c6c5593SMatthew G. Knepley   }
2769566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
2779566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
2788c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
2798c6c5593SMatthew G. Knepley }
2808c6c5593SMatthew G. Knepley 
281*d71ae5a4SJacob Faibussowitsch 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, PetscTabulation *T, PetscTabulation *TAux, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, PetscScalar values[])
282*d71ae5a4SJacob Faibussowitsch {
283b18799e0SMatthew G. Knepley   PetscSection       section, sectionAux = NULL;
2844bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
285b18799e0SMatthew G. Knepley   PetscScalar       *coefficients = NULL, *coefficientsAux = NULL;
286b18799e0SMatthew G. Knepley   PetscScalar       *coefficients_t = NULL, *coefficientsAux_t = NULL;
287b18799e0SMatthew G. Knepley   const PetscScalar *constants;
288b18799e0SMatthew G. Knepley   PetscReal         *x;
289ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc;
2909f209ee4SMatthew G. Knepley   PetscFEGeom        fegeom, cgeom;
2914bee2e38SMatthew G. Knepley   const PetscInt     dE = fgeom->dimEmbed;
292ef0bb6c7SMatthew G. Knepley   PetscInt           numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0;
293b18799e0SMatthew G. Knepley   PetscBool          isAffine;
294b18799e0SMatthew G. Knepley 
295b18799e0SMatthew G. Knepley   PetscFunctionBeginHot;
29608401ef6SPierre Jolivet   PetscCheck(dm == dmIn, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM");
2979566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2989566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
2999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
3009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
3019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x));
3029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
3039566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
3049566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
3059566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients));
306b18799e0SMatthew G. Knepley   if (dmAux) {
307a6e0b375SMatthew G. Knepley     PetscInt subp;
308b18799e0SMatthew G. Knepley 
3099566063dSJacob Faibussowitsch     PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
3109566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3119566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmAux, &sectionAux));
3129566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3139566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3149566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
3159566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
316b18799e0SMatthew G. Knepley   }
317b18799e0SMatthew G. Knepley   /* Get values for closure */
3184bee2e38SMatthew G. Knepley   isAffine       = fgeom->isAffine;
319ea78f98cSLisandro Dalcin   fegeom.n       = NULL;
320ea78f98cSLisandro Dalcin   fegeom.J       = NULL;
321ea78f98cSLisandro Dalcin   fegeom.v       = NULL;
322ea78f98cSLisandro Dalcin   fegeom.xi      = NULL;
323a6e0b375SMatthew G. Knepley   cgeom.dim      = fgeom->dim;
324a6e0b375SMatthew G. Knepley   cgeom.dimEmbed = fgeom->dimEmbed;
3254bee2e38SMatthew G. Knepley   if (isAffine) {
3264bee2e38SMatthew G. Knepley     fegeom.v    = x;
3274bee2e38SMatthew G. Knepley     fegeom.xi   = fgeom->xi;
3284bee2e38SMatthew G. Knepley     fegeom.J    = fgeom->J;
3294bee2e38SMatthew G. Knepley     fegeom.invJ = fgeom->invJ;
3304bee2e38SMatthew G. Knepley     fegeom.detJ = fgeom->detJ;
3314bee2e38SMatthew G. Knepley     fegeom.n    = fgeom->n;
3329f209ee4SMatthew G. Knepley 
3339f209ee4SMatthew G. Knepley     cgeom.J    = fgeom->suppJ[0];
3349f209ee4SMatthew G. Knepley     cgeom.invJ = fgeom->suppInvJ[0];
3359f209ee4SMatthew G. Knepley     cgeom.detJ = fgeom->suppDetJ[0];
3364bee2e38SMatthew G. Knepley   }
337b18799e0SMatthew G. Knepley   for (f = 0, v = 0; f < Nf; ++f) {
338b18799e0SMatthew G. Knepley     PetscQuadrature  allPoints;
339b18799e0SMatthew G. Knepley     PetscInt         q, dim, numPoints;
340b18799e0SMatthew G. Knepley     const PetscReal *points;
341b18799e0SMatthew G. Knepley     PetscScalar     *pointEval;
342b18799e0SMatthew G. Knepley     DM               dm;
343b18799e0SMatthew G. Knepley 
344b18799e0SMatthew G. Knepley     if (!sp[f]) continue;
3459566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
346b18799e0SMatthew G. Knepley     if (!funcs[f]) {
347b18799e0SMatthew G. Knepley       for (d = 0; d < spDim; d++, v++) values[v] = 0.;
348b18799e0SMatthew G. Knepley       continue;
349b18799e0SMatthew G. Knepley     }
3509566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
3519566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
3529566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
3539566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
354b18799e0SMatthew G. Knepley     for (q = 0; q < numPoints; ++q, ++tp) {
355b18799e0SMatthew G. Knepley       if (isAffine) {
356ef0bb6c7SMatthew G. Knepley         CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x);
357b18799e0SMatthew G. Knepley       } else {
3583fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[tp * dE];
3599f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[tp * dE * dE];
3609f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[tp * dE * dE];
3614bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[tp];
3624bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[tp * dE];
3639f209ee4SMatthew G. Knepley 
3649f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][tp * dE * dE];
3659f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE];
3669f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][tp];
367b18799e0SMatthew G. Knepley       }
368a6e0b375SMatthew 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 */
3699566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
3709566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
3714bee2e38SMatthew 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]);
372b18799e0SMatthew G. Knepley     }
3739566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
3749566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
375b18799e0SMatthew G. Knepley     v += spDim;
376b18799e0SMatthew G. Knepley   }
3779566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients));
3789566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
379b18799e0SMatthew G. Knepley   PetscFunctionReturn(0);
380b18799e0SMatthew G. Knepley }
381b18799e0SMatthew G. Knepley 
382*d71ae5a4SJacob Faibussowitsch 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[], PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[])
383*d71ae5a4SJacob Faibussowitsch {
3848c6c5593SMatthew G. Knepley   PetscFVCellGeom fvgeom;
3858c6c5593SMatthew G. Knepley   PetscInt        dim, dimEmbed;
3868c6c5593SMatthew G. Knepley 
3878c6c5593SMatthew G. Knepley   PetscFunctionBeginHot;
3889566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
3899566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
3909566063dSJacob Faibussowitsch   if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
3918c6c5593SMatthew G. Knepley   switch (type) {
3928c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL:
393*d71ae5a4SJacob Faibussowitsch   case DM_BC_NATURAL:
394*d71ae5a4SJacob Faibussowitsch     PetscCall(DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode(**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))funcs, ctxs, values));
395*d71ae5a4SJacob Faibussowitsch     break;
3968c6c5593SMatthew G. Knepley   case DM_BC_ESSENTIAL_FIELD:
3978c6c5593SMatthew G. Knepley   case DM_BC_NATURAL_FIELD:
3989371c9d4SSatish Balay     PetscCall(DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
399d0609cedSBarry Smith     break;
400b18799e0SMatthew G. Knepley   case DM_BC_ESSENTIAL_BD_FIELD:
4019371c9d4SSatish Balay     PetscCall(DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, (void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
402d0609cedSBarry Smith     break;
403*d71ae5a4SJacob Faibussowitsch   default:
404*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type);
4058c6c5593SMatthew G. Knepley   }
4068c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
4078c6c5593SMatthew G. Knepley }
4088c6c5593SMatthew G. Knepley 
409*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints)
410*d71ae5a4SJacob Faibussowitsch {
411df5c1128SToby Isaac   PetscReal *points;
412df5c1128SToby Isaac   PetscInt   f, numPoints;
413df5c1128SToby Isaac 
414df5c1128SToby Isaac   PetscFunctionBegin;
41519552267SMatthew G. Knepley   if (!dim) {
41619552267SMatthew G. Knepley     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
41719552267SMatthew G. Knepley     PetscFunctionReturn(0);
41819552267SMatthew G. Knepley   }
419df5c1128SToby Isaac   numPoints = 0;
420df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
421df5c1128SToby Isaac     if (funcs[f]) {
422df5c1128SToby Isaac       PetscQuadrature fAllPoints;
423df5c1128SToby Isaac       PetscInt        fNumPoints;
424df5c1128SToby Isaac 
4259566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
4269566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
427df5c1128SToby Isaac       numPoints += fNumPoints;
428df5c1128SToby Isaac     }
429df5c1128SToby Isaac   }
4309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * numPoints, &points));
431df5c1128SToby Isaac   numPoints = 0;
432df5c1128SToby Isaac   for (f = 0; f < Nf; ++f) {
433df5c1128SToby Isaac     if (funcs[f]) {
434df5c1128SToby Isaac       PetscQuadrature  fAllPoints;
43554ee0cceSMatthew G. Knepley       PetscInt         qdim, fNumPoints, q;
436df5c1128SToby Isaac       const PetscReal *fPoints;
437df5c1128SToby Isaac 
4389566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
4399566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
44063a3b9bcSJacob Faibussowitsch       PetscCheck(qdim == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %" PetscInt_FMT " for dual basis does not match input dimension %" PetscInt_FMT, qdim, dim);
441df5c1128SToby Isaac       for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q];
442df5c1128SToby Isaac       numPoints += fNumPoints;
443df5c1128SToby Isaac     }
444df5c1128SToby Isaac   }
4459566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
4469566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL));
447df5c1128SToby Isaac   PetscFunctionReturn(0);
448df5c1128SToby Isaac }
449df5c1128SToby Isaac 
4505f15299fSJeremy L Thompson /*@C
4515f15299fSJeremy L Thompson   DMGetFirstLabeledPoint - Find first labeled point p_o in odm such that the corresponding point p in dm has the specified height. Return p and the corresponding ds.
4525f15299fSJeremy L Thompson 
4535f15299fSJeremy L Thompson   Input Parameters:
4545f15299fSJeremy L Thompson   dm - the DM
4555f15299fSJeremy L Thompson   odm - the enclosing DM
4565f15299fSJeremy L Thompson   label - label for DM domain, or NULL for whole domain
4575f15299fSJeremy L Thompson   numIds - the number of ids
4585f15299fSJeremy L Thompson   ids - An array of the label ids in sequence for the domain
4595f15299fSJeremy L Thompson   height - Height of target cells in DMPlex topology
4605f15299fSJeremy L Thompson 
4615f15299fSJeremy L Thompson   Output Parameters:
4625f15299fSJeremy L Thompson   point - the first labeled point
4635f15299fSJeremy L Thompson   ds - the ds corresponding to the first labeled point
4645f15299fSJeremy L Thompson 
4655f15299fSJeremy L Thompson   Level: developer
466817da375SSatish Balay @*/
467*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
468*d71ae5a4SJacob Faibussowitsch {
469a6e0b375SMatthew G. Knepley   DM              plex;
470a6e0b375SMatthew G. Knepley   DMEnclosureType enc;
471e9f0ba4eSJed Brown   PetscInt        ls = -1;
472e5e52638SMatthew G. Knepley 
473e5e52638SMatthew G. Knepley   PetscFunctionBegin;
4745f15299fSJeremy L Thompson   if (point) *point = -1;
475e5e52638SMatthew G. Knepley   if (!label) PetscFunctionReturn(0);
4769566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
4779566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
478e9f0ba4eSJed Brown   for (PetscInt i = 0; i < numIds; ++i) {
479e9f0ba4eSJed Brown     IS       labelIS;
480e9f0ba4eSJed Brown     PetscInt num_points, pStart, pEnd;
4819566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
482e9f0ba4eSJed Brown     if (!labelIS) continue; /* No points with that id on this process */
4839566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
4849566063dSJacob Faibussowitsch     PetscCall(ISGetSize(labelIS, &num_points));
485e9f0ba4eSJed Brown     if (num_points) {
486e5e52638SMatthew G. Knepley       const PetscInt *points;
4879566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(labelIS, &points));
488e9f0ba4eSJed Brown       for (PetscInt i = 0; i < num_points; i++) {
489e9f0ba4eSJed Brown         PetscInt point;
4909566063dSJacob Faibussowitsch         PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
491e9f0ba4eSJed Brown         if (pStart <= point && point < pEnd) {
492a6e0b375SMatthew G. Knepley           ls = point;
4939566063dSJacob Faibussowitsch           if (ds) PetscCall(DMGetCellDS(dm, ls, ds));
494e5e52638SMatthew G. Knepley         }
495e9f0ba4eSJed Brown       }
4969566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(labelIS, &points));
497e9f0ba4eSJed Brown     }
4989566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&labelIS));
499e5e52638SMatthew G. Knepley     if (ls >= 0) break;
500e5e52638SMatthew G. Knepley   }
5015f15299fSJeremy L Thompson   if (point) *point = ls;
5029566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
503e5e52638SMatthew G. Knepley   PetscFunctionReturn(0);
504e5e52638SMatthew G. Knepley }
505e5e52638SMatthew G. Knepley 
5060de51b56SMatthew G. Knepley /*
5070de51b56SMatthew G. Knepley   This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
5080de51b56SMatthew G. Knepley 
5090de51b56SMatthew G. Knepley   There are several different scenarios:
5100de51b56SMatthew G. Knepley 
5110de51b56SMatthew G. Knepley   1) Volumetric mesh with volumetric auxiliary data
5120de51b56SMatthew G. Knepley 
5130de51b56SMatthew G. Knepley      Here minHeight=0 since we loop over cells.
5140de51b56SMatthew G. Knepley 
5150de51b56SMatthew G. Knepley   2) Boundary mesh with boundary auxiliary data
5160de51b56SMatthew G. Knepley 
5170de51b56SMatthew 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.
5180de51b56SMatthew G. Knepley 
5190de51b56SMatthew G. Knepley   3) Volumetric mesh with boundary auxiliary data
5200de51b56SMatthew G. Knepley 
5210de51b56SMatthew 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.
5220de51b56SMatthew G. Knepley 
523636b9322SMatthew G. Knepley   4) Volumetric input mesh with boundary output mesh
524636b9322SMatthew G. Knepley 
525636b9322SMatthew G. Knepley      Here we must get a subspace for the input DS
526636b9322SMatthew G. Knepley 
5270de51b56SMatthew G. Knepley   The maxHeight is used to support enforcement of constraints in DMForest.
5280de51b56SMatthew G. Knepley 
5290de51b56SMatthew G. Knepley   If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
5300de51b56SMatthew G. Knepley 
5310de51b56SMatthew 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.
5320de51b56SMatthew 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
5335f790a90SMatthew 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
5340de51b56SMatthew G. Knepley       dual spaces for the boundary from our input spaces.
5350de51b56SMatthew G. Knepley     - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
5360de51b56SMatthew G. Knepley 
5370de51b56SMatthew G. Knepley   We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
5380de51b56SMatthew G. Knepley 
5390de51b56SMatthew 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.
5400de51b56SMatthew G. Knepley */
541*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, InsertMode mode, Vec localX)
542*d71ae5a4SJacob Faibussowitsch {
543a6e0b375SMatthew G. Knepley   DM               plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
544a6e0b375SMatthew G. Knepley   DMEnclosureType  encIn, encAux;
545a6e0b375SMatthew G. Knepley   PetscDS          ds = NULL, dsIn = NULL, dsAux = NULL;
546ca3d3a14SMatthew G. Knepley   Vec              localA = NULL, tv;
547aa0eca99SMatthew G. Knepley   IS               fieldIS;
54847923291SMatthew G. Knepley   PetscSection     section;
549636b9322SMatthew G. Knepley   PetscDualSpace  *sp, *cellsp, *spIn, *cellspIn;
550ef0bb6c7SMatthew G. Knepley   PetscTabulation *T = NULL, *TAux = NULL;
5518c6c5593SMatthew G. Knepley   PetscInt        *Nc;
552636b9322SMatthew G. Knepley   PetscInt         dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
5535fedec97SMatthew G. Knepley   PetscBool       *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, transform;
554c330f8ffSToby Isaac   DMField          coordField;
555c330f8ffSToby Isaac   DMLabel          depthLabel;
556c330f8ffSToby Isaac   PetscQuadrature  allPoints = NULL;
55747923291SMatthew G. Knepley 
55847923291SMatthew G. Knepley   PetscFunctionBegin;
5599566063dSJacob Faibussowitsch   if (localU) PetscCall(VecGetDM(localU, &dmIn));
560ad540459SPierre Jolivet   else dmIn = dm;
5619566063dSJacob Faibussowitsch   PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
5629371c9d4SSatish Balay   if (localA) PetscCall(VecGetDM(localA, &dmAux));
563ad540459SPierre Jolivet   else dmAux = NULL;
5649566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
5659566063dSJacob Faibussowitsch   PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
5669566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
5679566063dSJacob Faibussowitsch   PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
5689566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5699566063dSJacob Faibussowitsch   PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
5709566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
5719566063dSJacob Faibussowitsch   PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
5729566063dSJacob Faibussowitsch   PetscCall(DMHasBasisTransform(dm, &transform));
5730de51b56SMatthew G. Knepley   /* Auxiliary information can only be used with interpolation of field functions */
574a6e0b375SMatthew G. Knepley   if (dmAux) {
5759566063dSJacob Faibussowitsch     PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
576ad540459SPierre Jolivet     if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) PetscCheck(localA, PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector");
577636b9322SMatthew G. Knepley   }
578e04ae0b4SMatthew G. Knepley   if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL));
579e04ae0b4SMatthew G. Knepley   PetscCall(DMGetCoordinateField(dm, &coordField));
580e04ae0b4SMatthew G. Knepley   /**** No collective calls below this point ****/
581636b9322SMatthew G. Knepley   /* Determine height for iteration of all meshes */
582636b9322SMatthew G. Knepley   {
583636b9322SMatthew G. Knepley     DMPolytopeType ct, ctIn, ctAux;
58419552267SMatthew G. Knepley     PetscInt       minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux;
58519552267SMatthew G. Knepley     PetscInt       dim = -1, dimIn = -1, dimAux = -1;
58688aca1feSMatthew G. Knepley 
5879566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
588636b9322SMatthew G. Knepley     if (pEnd > pStart) {
5899566063dSJacob Faibussowitsch       PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
590636b9322SMatthew G. Knepley       p = lStart < 0 ? pStart : lStart;
5919566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plex, p, &ct));
592636b9322SMatthew G. Knepley       dim = DMPolytopeTypeGetDim(ct);
5939566063dSJacob Faibussowitsch       PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
5949566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
5959566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
596636b9322SMatthew G. Knepley       dimIn = DMPolytopeTypeGetDim(ctIn);
597636b9322SMatthew G. Knepley       if (dmAux) {
5989566063dSJacob Faibussowitsch         PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
59919552267SMatthew G. Knepley         PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
60019552267SMatthew G. Knepley         if (pStartAux < pEndAux) {
6019566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
602636b9322SMatthew G. Knepley           dimAux = DMPolytopeTypeGetDim(ctAux);
60319552267SMatthew G. Knepley         }
604636b9322SMatthew G. Knepley       } else dimAux = dim;
605e04ae0b4SMatthew G. Knepley     } else {
606e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plex));
607e04ae0b4SMatthew G. Knepley       PetscCall(DMDestroy(&plexIn));
608e04ae0b4SMatthew G. Knepley       if (dmAux) PetscCall(DMDestroy(&plexAux));
609e04ae0b4SMatthew G. Knepley       PetscFunctionReturn(0);
610e39fcb4eSStefano Zampini     }
611636b9322SMatthew G. Knepley     if (dim < 0) {
612636b9322SMatthew G. Knepley       DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
613636b9322SMatthew G. Knepley 
614636b9322SMatthew G. Knepley       /* Fall back to determination based on being a submesh */
6159566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plex, &spmap));
6169566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
6179566063dSJacob Faibussowitsch       if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
618636b9322SMatthew G. Knepley       dim    = spmap ? 1 : 0;
619636b9322SMatthew G. Knepley       dimIn  = spmapIn ? 1 : 0;
620636b9322SMatthew G. Knepley       dimAux = spmapAux ? 1 : 0;
62188aca1feSMatthew G. Knepley     }
622636b9322SMatthew G. Knepley     {
62319552267SMatthew G. Knepley       PetscInt dimProj   = PetscMin(PetscMin(dim, dimIn), (dimAux < 0 ? PETSC_MAX_INT : dimAux));
62419552267SMatthew G. Knepley       PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;
625636b9322SMatthew G. Knepley 
62619552267SMatthew G. Knepley       PetscCheck(PetscAbsInt(dimProj - dim) <= 1 && PetscAbsInt(dimProj - dimIn) <= 1 && PetscAbsInt(dimProj - dimAuxEff) <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not currently support differences of more than 1 in dimension");
627636b9322SMatthew G. Knepley       if (dimProj < dim) minHeight = 1;
628636b9322SMatthew G. Knepley       htInc    = dim - dimProj;
629636b9322SMatthew G. Knepley       htIncIn  = dimIn - dimProj;
63019552267SMatthew G. Knepley       htIncAux = dimAuxEff - dimProj;
6310de51b56SMatthew G. Knepley     }
632a6e0b375SMatthew G. Knepley   }
6339566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
6349566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
6359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
6362af58022SMatthew G. Knepley   maxHeight = PetscMax(maxHeight, minHeight);
63763a3b9bcSJacob Faibussowitsch   PetscCheck(maxHeight >= 0 && maxHeight <= dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", maxHeight, dim);
6389566063dSJacob Faibussowitsch   PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, 0, NULL, &ds));
6399566063dSJacob Faibussowitsch   if (!ds) PetscCall(DMGetDS(dm, &ds));
6409566063dSJacob Faibussowitsch   PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn));
6419566063dSJacob Faibussowitsch   if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
6429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6439566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
6449566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &NfTot));
6459566063dSJacob Faibussowitsch   PetscCall(DMFindRegionNum(dm, ds, &regionNum));
6469566063dSJacob Faibussowitsch   PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL));
6479566063dSJacob Faibussowitsch   PetscCall(PetscDSIsCohesive(ds, &isCohesive));
6489566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
6499566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
6500c898c18SMatthew G. Knepley   if (dmAux) {
6519566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dmAux, &dsAux));
6529566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6530c898c18SMatthew G. Knepley   }
6549566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponents(ds, &Nc));
6559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
6569566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
6579371c9d4SSatish Balay   else {
6589371c9d4SSatish Balay     cellsp   = sp;
6599371c9d4SSatish Balay     cellspIn = spIn;
6609371c9d4SSatish Balay   }
6618c6c5593SMatthew G. Knepley   /* Get cell dual spaces */
66247923291SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
663665f567fSMatthew G. Knepley     PetscDiscType disctype;
66447923291SMatthew G. Knepley 
6659566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
666665f567fSMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
667665f567fSMatthew G. Knepley       PetscFE fe;
66847923291SMatthew G. Knepley 
66947923291SMatthew G. Knepley       isFE[f] = PETSC_TRUE;
670665f567fSMatthew G. Knepley       hasFE   = PETSC_TRUE;
6719566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
6729566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
673665f567fSMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
674665f567fSMatthew G. Knepley       PetscFV fv;
67547923291SMatthew G. Knepley 
67647923291SMatthew G. Knepley       isFE[f] = PETSC_FALSE;
677665f567fSMatthew G. Knepley       hasFV   = PETSC_TRUE;
6789566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv));
6799566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
680665f567fSMatthew G. Knepley     } else {
681665f567fSMatthew G. Knepley       isFE[f]   = PETSC_FALSE;
682665f567fSMatthew G. Knepley       cellsp[f] = NULL;
683665f567fSMatthew G. Knepley     }
68447923291SMatthew G. Knepley   }
685636b9322SMatthew G. Knepley   for (f = 0; f < NfIn; ++f) {
686636b9322SMatthew G. Knepley     PetscDiscType disctype;
687636b9322SMatthew G. Knepley 
6889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
689636b9322SMatthew G. Knepley     if (disctype == PETSC_DISC_FE) {
690636b9322SMatthew G. Knepley       PetscFE fe;
691636b9322SMatthew G. Knepley 
6929566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
6939566063dSJacob Faibussowitsch       PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
694636b9322SMatthew G. Knepley     } else if (disctype == PETSC_DISC_FV) {
695636b9322SMatthew G. Knepley       PetscFV fv;
696636b9322SMatthew G. Knepley 
6979566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv));
6989566063dSJacob Faibussowitsch       PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
699636b9322SMatthew G. Knepley     } else {
700636b9322SMatthew G. Knepley       cellspIn[f] = NULL;
701636b9322SMatthew G. Knepley     }
702636b9322SMatthew G. Knepley   }
703636b9322SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
7049371c9d4SSatish Balay     if (!htInc) {
7059371c9d4SSatish Balay       sp[f] = cellsp[f];
7069371c9d4SSatish Balay     } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
707636b9322SMatthew G. Knepley   }
708ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
70974fc349aSMatthew G. Knepley     PetscFE          fem, subfem;
710665f567fSMatthew G. Knepley     PetscDiscType    disctype;
7114a3e9fdbSToby Isaac     const PetscReal *points;
7124a3e9fdbSToby Isaac     PetscInt         numPoints;
7130c898c18SMatthew G. Knepley 
71408401ef6SPierre Jolivet     PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
7159566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints));
7169566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
7179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
718a6e0b375SMatthew G. Knepley     for (f = 0; f < NfIn; ++f) {
7199371c9d4SSatish Balay       if (!htIncIn) {
7209371c9d4SSatish Balay         spIn[f] = cellspIn[f];
7219371c9d4SSatish Balay       } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
722636b9322SMatthew G. Knepley 
7239566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
724665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
7259566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem));
7269371c9d4SSatish Balay       if (!htIncIn) {
7279371c9d4SSatish Balay         subfem = fem;
7289371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
7299566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]));
7300c898c18SMatthew G. Knepley     }
7310c898c18SMatthew G. Knepley     for (f = 0; f < NfAux; ++f) {
7329566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
733665f567fSMatthew G. Knepley       if (disctype != PETSC_DISC_FE) continue;
7349566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem));
7359371c9d4SSatish Balay       if (!htIncAux) {
7369371c9d4SSatish Balay         subfem = fem;
7379371c9d4SSatish Balay       } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
7389566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]));
7390c898c18SMatthew G. Knepley     }
7400c898c18SMatthew G. Knepley   }
74147923291SMatthew G. Knepley   /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
7422af58022SMatthew G. Knepley   for (h = minHeight; h <= maxHeight; h++) {
743636b9322SMatthew G. Knepley     PetscInt     hEff     = h - minHeight + htInc;
744636b9322SMatthew G. Knepley     PetscInt     hEffIn   = h - minHeight + htIncIn;
745636b9322SMatthew G. Knepley     PetscInt     hEffAux  = h - minHeight + htIncAux;
746a6e0b375SMatthew G. Knepley     PetscDS      dsEff    = ds;
747636b9322SMatthew G. Knepley     PetscDS      dsEffIn  = dsIn;
748636b9322SMatthew G. Knepley     PetscDS      dsEffAux = dsAux;
7498c6c5593SMatthew G. Knepley     PetscScalar *values;
750b7260050SToby Isaac     PetscBool   *fieldActive;
751b7260050SToby Isaac     PetscInt     maxDegree;
752e5e52638SMatthew G. Knepley     PetscInt     pStart, pEnd, p, lStart, spDim, totDim, numValues;
753c330f8ffSToby Isaac     IS           heightIS;
7548c6c5593SMatthew G. Knepley 
755636b9322SMatthew G. Knepley     if (h > minHeight) {
7569566063dSJacob Faibussowitsch       for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
757636b9322SMatthew G. Knepley     }
7589566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
7599566063dSJacob Faibussowitsch     PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
7609566063dSJacob Faibussowitsch     PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
761c330f8ffSToby Isaac     if (pEnd <= pStart) {
7629566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
763c330f8ffSToby Isaac       continue;
764c330f8ffSToby Isaac     }
76547923291SMatthew G. Knepley     /* Compute totDim, the number of dofs in the closure of a point at this height */
76647923291SMatthew G. Knepley     totDim = 0;
76747923291SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
7685fedec97SMatthew G. Knepley       PetscBool cohesive;
7695fedec97SMatthew G. Knepley 
770665f567fSMatthew G. Knepley       if (!sp[f]) continue;
7719566063dSJacob Faibussowitsch       PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
7729566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
7739c3cf19fSMatthew G. Knepley       totDim += spDim;
7745fedec97SMatthew G. Knepley       if (isCohesive && !cohesive) totDim += spDim;
77547923291SMatthew G. Knepley     }
776636b9322SMatthew G. Knepley     p = lStart < 0 ? pStart : lStart;
7779566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
77863a3b9bcSJacob Faibussowitsch     PetscCheck(numValues == totDim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The output section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT " in [%" PetscInt_FMT ", %" PetscInt_FMT "]", p, numValues, totDim, h, minHeight, maxHeight);
779c330f8ffSToby Isaac     if (!totDim) {
7809566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&heightIS));
781c330f8ffSToby Isaac       continue;
782c330f8ffSToby Isaac     }
7839566063dSJacob Faibussowitsch     if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
784636b9322SMatthew G. Knepley     /* Compute totDimIn, the number of dofs in the closure of a point at this height */
785636b9322SMatthew G. Knepley     if (localU) {
786636b9322SMatthew G. Knepley       PetscInt totDimIn, pIn, numValuesIn;
787636b9322SMatthew G. Knepley 
788636b9322SMatthew G. Knepley       totDimIn = 0;
789636b9322SMatthew G. Knepley       for (f = 0; f < NfIn; ++f) {
7905fedec97SMatthew G. Knepley         PetscBool cohesive;
7915fedec97SMatthew G. Knepley 
792636b9322SMatthew G. Knepley         if (!spIn[f]) continue;
7939566063dSJacob Faibussowitsch         PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
7949566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
795636b9322SMatthew G. Knepley         totDimIn += spDim;
7965fedec97SMatthew G. Knepley         if (isCohesive && !cohesive) totDimIn += spDim;
797636b9322SMatthew G. Knepley       }
7989566063dSJacob Faibussowitsch       PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
7999566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
80063a3b9bcSJacob Faibussowitsch       PetscCheck(numValuesIn == totDimIn, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The input section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT, pIn, numValuesIn, totDimIn, htIncIn);
8019566063dSJacob Faibussowitsch       if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
802636b9322SMatthew G. Knepley     }
8039566063dSJacob Faibussowitsch     if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
80447923291SMatthew G. Knepley     /* Loop over points at this height */
8059566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
8069566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
807aa0eca99SMatthew G. Knepley     {
808aa0eca99SMatthew G. Knepley       const PetscInt *fields;
809aa0eca99SMatthew G. Knepley 
8109566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
811ad540459SPierre Jolivet       for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE;
812ad540459SPierre Jolivet       for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
8139566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
814aa0eca99SMatthew G. Knepley     }
8158c6c5593SMatthew G. Knepley     if (label) {
8168c6c5593SMatthew G. Knepley       PetscInt i;
81747923291SMatthew G. Knepley 
81847923291SMatthew G. Knepley       for (i = 0; i < numIds; ++i) {
819c330f8ffSToby Isaac         IS              pointIS, isectIS;
82047923291SMatthew G. Knepley         const PetscInt *points;
8218c6c5593SMatthew G. Knepley         PetscInt        n;
822c330f8ffSToby Isaac         PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
823c330f8ffSToby Isaac         PetscQuadrature quad = NULL;
82447923291SMatthew G. Knepley 
8259566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
82647923291SMatthew G. Knepley         if (!pointIS) continue; /* No points with that id on this process */
8279566063dSJacob Faibussowitsch         PetscCall(ISIntersect(pointIS, heightIS, &isectIS));
8289566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&pointIS));
829c330f8ffSToby Isaac         if (!isectIS) continue;
8309566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(isectIS, &n));
8319566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(isectIS, &points));
8329566063dSJacob Faibussowitsch         PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree));
83348a46eb9SPierre Jolivet         if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad));
834c330f8ffSToby Isaac         if (!quad) {
835c330f8ffSToby Isaac           if (!h && allPoints) {
836c330f8ffSToby Isaac             quad      = allPoints;
837c330f8ffSToby Isaac             allPoints = NULL;
838c330f8ffSToby Isaac           } else {
8399566063dSJacob Faibussowitsch             PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad));
840c330f8ffSToby Isaac           }
841c330f8ffSToby Isaac         }
8429566063dSJacob Faibussowitsch         PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
84347923291SMatthew G. Knepley         for (p = 0; p < n; ++p) {
84447923291SMatthew G. Knepley           const PetscInt point = points[p];
84547923291SMatthew G. Knepley 
8469566063dSJacob Faibussowitsch           PetscCall(PetscArrayzero(values, numValues));
8479566063dSJacob Faibussowitsch           PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom));
8489566063dSJacob Faibussowitsch           PetscCall(DMPlexSetActivePoint(dm, point));
849d0609cedSBarry Smith           PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values));
8509566063dSJacob Faibussowitsch           if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
8519566063dSJacob Faibussowitsch           PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
85247923291SMatthew G. Knepley         }
8539566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom));
8549566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomDestroy(&fegeom));
8559566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureDestroy(&quad));
8569566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(isectIS, &points));
8579566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isectIS));
85847923291SMatthew G. Knepley       }
8598c6c5593SMatthew G. Knepley     } else {
860c330f8ffSToby Isaac       PetscFEGeom    *fegeom = NULL, *chunkgeom = NULL;
861c330f8ffSToby Isaac       PetscQuadrature quad = NULL;
862c330f8ffSToby Isaac       IS              pointIS;
863c330f8ffSToby Isaac 
8649566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS));
8659566063dSJacob Faibussowitsch       PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
86648a46eb9SPierre Jolivet       if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
867c330f8ffSToby Isaac       if (!quad) {
868c330f8ffSToby Isaac         if (!h && allPoints) {
869c330f8ffSToby Isaac           quad      = allPoints;
870c330f8ffSToby Isaac           allPoints = NULL;
871c330f8ffSToby Isaac         } else {
8729566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad));
873c330f8ffSToby Isaac         }
874c330f8ffSToby Isaac       }
8759566063dSJacob Faibussowitsch       PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom));
8768c6c5593SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
8779566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(values, numValues));
8789566063dSJacob Faibussowitsch         PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom));
8799566063dSJacob Faibussowitsch         PetscCall(DMPlexSetActivePoint(dm, p));
880d0609cedSBarry Smith         PetscCall(DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsEffIn, plexAux, encAux, dsEffAux, chunkgeom, htInc, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values));
8819566063dSJacob Faibussowitsch         if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
8829566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
8838c6c5593SMatthew G. Knepley       }
8849566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom));
8859566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomDestroy(&fegeom));
8869566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureDestroy(&quad));
8879566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
8888c6c5593SMatthew G. Knepley     }
8899566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&heightIS));
8909566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
8919566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
89247923291SMatthew G. Knepley   }
8938c6c5593SMatthew G. Knepley   /* Cleanup */
894ece3a9fcSMatthew G. Knepley   if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
8959566063dSJacob Faibussowitsch     for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f]));
8969566063dSJacob Faibussowitsch     for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
8979566063dSJacob Faibussowitsch     PetscCall(PetscFree2(T, TAux));
8980c898c18SMatthew G. Knepley   }
8999566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&allPoints));
9009566063dSJacob Faibussowitsch   PetscCall(PetscFree3(isFE, sp, spIn));
9019566063dSJacob Faibussowitsch   if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
9029566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
9039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plexIn));
9049566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMDestroy(&plexAux));
9058c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
90647923291SMatthew G. Knepley }
9078c6c5593SMatthew G. Knepley 
908*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
909*d71ae5a4SJacob Faibussowitsch {
9108c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9119566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
9128c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
9138c6c5593SMatthew G. Knepley }
9148c6c5593SMatthew G. Knepley 
915*d71ae5a4SJacob Faibussowitsch 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)
916*d71ae5a4SJacob Faibussowitsch {
9178c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9189566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX));
91947923291SMatthew G. Knepley   PetscFunctionReturn(0);
92047923291SMatthew G. Knepley }
92147923291SMatthew G. Knepley 
922*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
923*d71ae5a4SJacob Faibussowitsch {
92447923291SMatthew G. Knepley   PetscFunctionBegin;
9259566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
9268c6c5593SMatthew G. Knepley   PetscFunctionReturn(0);
92747923291SMatthew G. Knepley }
92847923291SMatthew G. Knepley 
929*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
930*d71ae5a4SJacob Faibussowitsch {
9318c6c5593SMatthew G. Knepley   PetscFunctionBegin;
9329566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX));
93347923291SMatthew G. Knepley   PetscFunctionReturn(0);
93447923291SMatthew G. Knepley }
935ece3a9fcSMatthew G. Knepley 
936*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec localX)
937*d71ae5a4SJacob Faibussowitsch {
938ece3a9fcSMatthew G. Knepley   PetscFunctionBegin;
9399566063dSJacob Faibussowitsch   PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX));
940ece3a9fcSMatthew G. Knepley   PetscFunctionReturn(0);
941ece3a9fcSMatthew G. Knepley }
942