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
820f4b53cSBarry Smith Not Collective
91b32699bSMatthew G. Knepley
104165533cSJose E. Roman Input Parameter:
11a1cb98faSBarry Smith . 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
181cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`
191b32699bSMatthew G. Knepley @*/
DMPlexGetActivePoint(DM dm,PetscInt * point)20d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point)
21d71ae5a4SJacob Faibussowitsch {
221b32699bSMatthew G. Knepley PetscFunctionBeginHot;
231b32699bSMatthew G. Knepley *point = ((DM_Plex *)dm->data)->activePoint;
243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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
3020f4b53cSBarry Smith Not Collective
311b32699bSMatthew G. Knepley
324165533cSJose E. Roman Input Parameters:
33a1cb98faSBarry Smith + 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
381cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexGetActivePoint()`
391b32699bSMatthew G. Knepley @*/
DMPlexSetActivePoint(DM dm,PetscInt point)40d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point)
41d71ae5a4SJacob Faibussowitsch {
421b32699bSMatthew G. Knepley PetscFunctionBeginHot;
431b32699bSMatthew G. Knepley ((DM_Plex *)dm->data)->activePoint = point;
443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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:
51a1cb98faSBarry Smith + dm - The output `DM`
52a1cb98faSBarry Smith . ds - The output `DS`
53a1cb98faSBarry Smith . dmIn - The input `DM`
54a1cb98faSBarry Smith . 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
59a1cb98faSBarry Smith . 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
681cc06b55SBarry Smith .seealso:[](ch_unstructured), `DM`, `DMPLEX`, `PetscDS`, `PetscFEGeom`, `PetscFVCellGeom`, `PetscDualSpace`
69a6e0b375SMatthew G. Knepley */
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[])70d71ae5a4SJacob 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[])
71d71ae5a4SJacob Faibussowitsch {
72a77a5016SMatthew G. Knepley PetscInt debug = ((DM_Plex *)dm->data)->printProject;
73a6e0b375SMatthew G. Knepley PetscInt coordDim, Nf, *Nc, f, spDim, d, v, tp;
745fedec97SMatthew G. Knepley PetscBool isAffine, isCohesive, transform;
758c6c5593SMatthew G. Knepley
768c6c5593SMatthew G. Knepley PetscFunctionBeginHot;
779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dmIn, &coordDim));
789566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dmIn, &transform));
799566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
809566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc));
819566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive));
828c6c5593SMatthew G. Knepley /* Get values for closure */
83c330f8ffSToby Isaac isAffine = fegeom->isAffine;
84c330f8ffSToby Isaac for (f = 0, v = 0, tp = 0; f < Nf; ++f) {
858c6c5593SMatthew G. Knepley void *const ctx = ctxs ? ctxs[f] : NULL;
865fedec97SMatthew G. Knepley PetscBool cohesive;
878c6c5593SMatthew G. Knepley
888c6c5593SMatthew G. Knepley if (!sp[f]) continue;
899566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
909566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
918c6c5593SMatthew G. Knepley if (funcs[f]) {
92c330f8ffSToby Isaac if (isFE[f]) {
93c330f8ffSToby Isaac PetscQuadrature allPoints;
94c330f8ffSToby Isaac PetscInt q, dim, numPoints;
95c330f8ffSToby Isaac const PetscReal *points;
96c330f8ffSToby Isaac PetscScalar *pointEval;
97c330f8ffSToby Isaac PetscReal *x;
98ca3d3a14SMatthew G. Knepley DM rdm;
99c330f8ffSToby Isaac
1009566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &rdm));
1019566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
1029566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
1039566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
1049566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x));
105737e23dcSMatthew G. Knepley PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f]));
106c330f8ffSToby Isaac for (q = 0; q < numPoints; q++, tp++) {
107c330f8ffSToby Isaac const PetscReal *v0;
108c330f8ffSToby Isaac
109c330f8ffSToby Isaac if (isAffine) {
110665f567fSMatthew G. Knepley const PetscReal *refpoint = &points[q * dim];
111665f567fSMatthew G. Knepley PetscReal injpoint[3] = {0., 0., 0.};
112665f567fSMatthew G. Knepley
113665f567fSMatthew G. Knepley if (dim != fegeom->dim) {
114966bd95aSPierre Jolivet PetscCheck(isCohesive, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim);
115665f567fSMatthew G. Knepley /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */
116665f567fSMatthew G. Knepley for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d];
117665f567fSMatthew G. Knepley refpoint = injpoint;
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 }
128a77a5016SMatthew G. Knepley if (debug > 3) {
129a77a5016SMatthew G. Knepley PetscInt ap;
130a77a5016SMatthew G. Knepley PetscCall(DMPlexGetActivePoint(dm, &ap));
131a77a5016SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "Project point %" PetscInt_FMT ", analytic: ref (", ap));
132a77a5016SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) {
133a77a5016SMatthew G. Knepley if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
134a77a5016SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)points[q * dim + d]));
135a77a5016SMatthew G. Knepley }
136a77a5016SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, ") real ("));
137a77a5016SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) {
138a77a5016SMatthew G. Knepley if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
139a77a5016SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)v0[d]));
140a77a5016SMatthew G. Knepley }
141a77a5016SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
142a77a5016SMatthew G. Knepley }
1439566063dSJacob Faibussowitsch PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx));
144c330f8ffSToby Isaac }
1454bee2e38SMatthew G. Knepley /* Transform point evaluations pointEval[q,c] */
1469566063dSJacob Faibussowitsch PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval));
1479566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
1489566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x));
1499566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
150c330f8ffSToby Isaac v += spDim;
1515fedec97SMatthew G. Knepley if (isCohesive && !cohesive) {
15227f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
15327f02ce8SMatthew G. Knepley }
154c330f8ffSToby Isaac } else {
15548a46eb9SPierre Jolivet for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]));
156c330f8ffSToby Isaac }
157c330f8ffSToby Isaac } else {
15827f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.;
1595fedec97SMatthew G. Knepley if (isCohesive && !cohesive) {
16027f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.;
16127f02ce8SMatthew G. Knepley }
1628c6c5593SMatthew G. Knepley }
1639c3cf19fSMatthew G. Knepley }
1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1658c6c5593SMatthew G. Knepley }
1668c6c5593SMatthew G. Knepley
167a6e0b375SMatthew G. Knepley /*
168a6e0b375SMatthew G. Knepley DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point
169a6e0b375SMatthew G. Knepley
170a6e0b375SMatthew G. Knepley Input Parameters:
171a6e0b375SMatthew G. Knepley + dm - The output DM
172a6e0b375SMatthew G. Knepley . ds - The output DS
173a6e0b375SMatthew G. Knepley . dmIn - The input DM
174a6e0b375SMatthew G. Knepley . dsIn - The input DS
175a6e0b375SMatthew G. Knepley . dmAux - The auxiliary DM, which is always for the input space
176a6e0b375SMatthew G. Knepley . dsAux - The auxiliary DS, which is always for the input space
177a6e0b375SMatthew G. Knepley . time - The time for this evaluation
178a6e0b375SMatthew G. Knepley . localU - The local solution
179a6e0b375SMatthew G. Knepley . localA - The local auziliary fields
180a6e0b375SMatthew G. Knepley . cgeom - The FE geometry for this point
181a6e0b375SMatthew G. Knepley . sp - The output PetscDualSpace for each field
182a6e0b375SMatthew G. Knepley . p - The point in the output DM
183a5b23f4aSJose E. Roman . T - Input basis and derivatives for each field tabulated on the quadrature points
184ef0bb6c7SMatthew G. Knepley . TAux - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points
185a6e0b375SMatthew G. Knepley . funcs - The evaluation function for each field
186a6e0b375SMatthew G. Knepley - ctxs - The user context for each field
187a6e0b375SMatthew G. Knepley
188a6e0b375SMatthew G. Knepley Output Parameter:
189a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space
190a6e0b375SMatthew G. Knepley
191a6e0b375SMatthew G. Knepley Level: developer
192a6e0b375SMatthew G. Knepley
193a1cb98faSBarry Smith Note:
194a1cb98faSBarry Smith Not supported for FV
195a1cb98faSBarry Smith
196db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()`
197a6e0b375SMatthew G. Knepley */
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[])198d71ae5a4SJacob 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[])
199d71ae5a4SJacob Faibussowitsch {
2008c6c5593SMatthew G. Knepley PetscSection section, sectionAux = NULL;
2014bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
2028c6c5593SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
2038c6c5593SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL;
204191494d9SMatthew G. Knepley const PetscScalar *constants;
2058c6c5593SMatthew G. Knepley PetscReal *x;
2068e6d238bSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
207989fa639SBrad Aagaard PetscFEGeom fegeom, fgeomN[2];
2088e6d238bSMatthew G. Knepley const PetscInt dE = cgeom->dimEmbed, *cone, *ornt;
209ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
2108e6d238bSMatthew G. Knepley PetscBool isAffine, isCohesive, isCohesiveIn, transform;
2118e6d238bSMatthew G. Knepley DMPolytopeType qct;
2128c6c5593SMatthew G. Knepley
2138c6c5593SMatthew G. Knepley PetscFunctionBeginHot;
2149566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
2159566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc));
2169566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive));
2179566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
2188e6d238bSMatthew G. Knepley PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
2199566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
2209566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
2219566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
2229566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
2239566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
2249566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dmIn, &transform));
2259566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmIn, §ion));
2269566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
2278e6d238bSMatthew G. Knepley // Get cohesive cell hanging off face
2288e6d238bSMatthew G. Knepley if (isCohesiveIn) {
2298e6d238bSMatthew G. Knepley PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
2308e6d238bSMatthew G. Knepley if ((qct != DM_POLYTOPE_POINT_PRISM_TENSOR) && (qct != DM_POLYTOPE_SEG_PRISM_TENSOR) && (qct != DM_POLYTOPE_TRI_PRISM_TENSOR) && (qct != DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
2318e6d238bSMatthew G. Knepley DMPolytopeType ct;
2328e6d238bSMatthew G. Knepley const PetscInt *support;
2338e6d238bSMatthew G. Knepley PetscInt Ns, s;
2348e6d238bSMatthew G. Knepley
2358e6d238bSMatthew G. Knepley PetscCall(DMPlexGetSupport(dmIn, inp, &support));
2368e6d238bSMatthew G. Knepley PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
2378e6d238bSMatthew G. Knepley for (s = 0; s < Ns; ++s) {
2388e6d238bSMatthew G. Knepley PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
2398e6d238bSMatthew G. Knepley if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
2408e6d238bSMatthew G. Knepley inp = support[s];
2418e6d238bSMatthew G. Knepley break;
2428e6d238bSMatthew G. Knepley }
2438e6d238bSMatthew G. Knepley }
2448e6d238bSMatthew G. Knepley PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
2458e6d238bSMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
2468e6d238bSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
2478e6d238bSMatthew G. Knepley face[0] = 0;
2488e6d238bSMatthew G. Knepley face[1] = 0;
2498e6d238bSMatthew G. Knepley }
2508e6d238bSMatthew G. Knepley }
251eb8f539aSJed Brown if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
2528c6c5593SMatthew G. Knepley if (dmAux) {
25344171101SMatthew G. Knepley PetscInt subp;
2541ba34526SMatthew G. Knepley
2559566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
2569566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2579566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmAux, §ionAux));
2589566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2599566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2609566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
2619566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
2628c6c5593SMatthew G. Knepley }
2638c6c5593SMatthew G. Knepley /* Get values for closure */
2644bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine;
26527f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim;
26627f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed;
267989fa639SBrad Aagaard if (isCohesiveIn) {
268989fa639SBrad Aagaard fgeomN[0].dim = cgeom->dim;
269989fa639SBrad Aagaard fgeomN[0].dimEmbed = cgeom->dimEmbed;
270989fa639SBrad Aagaard fgeomN[1].dim = cgeom->dim;
271989fa639SBrad Aagaard fgeomN[1].dimEmbed = cgeom->dimEmbed;
272989fa639SBrad Aagaard }
2734bee2e38SMatthew G. Knepley if (isAffine) {
2744bee2e38SMatthew G. Knepley fegeom.v = x;
2754bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi;
2764bee2e38SMatthew G. Knepley fegeom.J = cgeom->J;
2774bee2e38SMatthew G. Knepley fegeom.invJ = cgeom->invJ;
2784bee2e38SMatthew G. Knepley fegeom.detJ = cgeom->detJ;
279989fa639SBrad Aagaard if (isCohesiveIn) {
280989fa639SBrad Aagaard fgeomN[0].J = cgeom->suppJ[0];
281989fa639SBrad Aagaard fgeomN[0].invJ = cgeom->suppInvJ[0];
282989fa639SBrad Aagaard fgeomN[0].detJ = cgeom->suppDetJ[0];
283989fa639SBrad Aagaard fgeomN[1].J = cgeom->suppJ[1];
284989fa639SBrad Aagaard fgeomN[1].invJ = cgeom->suppInvJ[1];
285989fa639SBrad Aagaard fgeomN[1].detJ = cgeom->suppDetJ[1];
286989fa639SBrad Aagaard }
2874bee2e38SMatthew G. Knepley }
2888c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) {
289c330f8ffSToby Isaac PetscQuadrature allPoints;
290c330f8ffSToby Isaac PetscInt q, dim, numPoints;
291c330f8ffSToby Isaac const PetscReal *points;
292c330f8ffSToby Isaac PetscScalar *pointEval;
2935fedec97SMatthew G. Knepley PetscBool cohesive;
294c330f8ffSToby Isaac DM dm;
295c330f8ffSToby Isaac
2968c6c5593SMatthew G. Knepley if (!sp[f]) continue;
2979566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
2989566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
299c330f8ffSToby Isaac if (!funcs[f]) {
300be1504a2SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.;
3015fedec97SMatthew G. Knepley if (isCohesive && !cohesive) {
30227f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.;
30327f02ce8SMatthew G. Knepley }
304c330f8ffSToby Isaac continue;
305c330f8ffSToby Isaac }
3066f0e67eaSMatthew G. Knepley const PetscInt ***perms;
3079566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
3086f0e67eaSMatthew G. Knepley PetscCall(PetscDualSpaceGetSymmetries(sp[f], &perms, NULL));
3099566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
3109566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
3119566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
3120c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) {
3138e6d238bSMatthew G. Knepley PetscInt qpt[2];
3148e6d238bSMatthew G. Knepley
3158e6d238bSMatthew G. Knepley if (isCohesiveIn) {
3166f0e67eaSMatthew G. Knepley qpt[0] = perms ? perms[0][ornt[0]][q] : q;
3176f0e67eaSMatthew G. Knepley qpt[1] = perms ? perms[0][DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0)][q] : q;
3188e6d238bSMatthew G. Knepley }
319c330f8ffSToby Isaac if (isAffine) {
320ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x);
3211c531cf8SMatthew G. Knepley } else {
3224bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[tp * dE];
3234bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[tp * dE * dE];
3244bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[tp * dE * dE];
3254bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[tp];
326989fa639SBrad Aagaard if (isCohesiveIn) {
327989fa639SBrad Aagaard fgeomN[0].J = &cgeom->suppJ[0][tp * dE * dE];
328989fa639SBrad Aagaard fgeomN[0].invJ = &cgeom->suppInvJ[0][tp * dE * dE];
329989fa639SBrad Aagaard fgeomN[0].detJ = &cgeom->suppDetJ[0][tp];
330989fa639SBrad Aagaard fgeomN[1].J = &cgeom->suppJ[1][tp * dE * dE];
331989fa639SBrad Aagaard fgeomN[1].invJ = &cgeom->suppInvJ[1][tp * dE * dE];
332989fa639SBrad Aagaard fgeomN[1].detJ = &cgeom->suppDetJ[1][tp];
333989fa639SBrad Aagaard }
3348c6c5593SMatthew G. Knepley }
3358e6d238bSMatthew G. Knepley if (coefficients) {
336989fa639SBrad Aagaard if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, fgeomN, coefficients, coefficients_t, u, u_x, u_t));
3378e6d238bSMatthew G. Knepley else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t));
3388e6d238bSMatthew G. Knepley }
3399566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
3409566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
341a6e0b375SMatthew 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]);
3421c531cf8SMatthew G. Knepley }
3439566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
3449566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
345c330f8ffSToby Isaac v += spDim;
34627f02ce8SMatthew G. Knepley /* TODO: For now, set both sides equal, but this should use info from other support cell */
3475fedec97SMatthew G. Knepley if (isCohesive && !cohesive) {
34827f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
34927f02ce8SMatthew G. Knepley }
3508c6c5593SMatthew G. Knepley }
351eb8f539aSJed Brown if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
3529566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
3538e6d238bSMatthew G. Knepley if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3558c6c5593SMatthew G. Knepley }
3568c6c5593SMatthew G. Knepley
DMProjectPoint_BdField_Private(DM dm,PetscDS ds,DM dmIn,DMEnclosureType encIn,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[])35779f2dbaeSMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, 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[])
358d71ae5a4SJacob Faibussowitsch {
359b18799e0SMatthew G. Knepley PetscSection section, sectionAux = NULL;
3604bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc;
361b18799e0SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
362b18799e0SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL;
363b18799e0SMatthew G. Knepley const PetscScalar *constants;
364b18799e0SMatthew G. Knepley PetscReal *x;
36579f2dbaeSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2];
366989fa639SBrad Aagaard PetscFEGeom fegeom, cgeom, fgeomN[2];
36779f2dbaeSMatthew G. Knepley const PetscInt dE = fgeom->dimEmbed, *cone, *ornt;
36879f2dbaeSMatthew G. Knepley PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0;
36979f2dbaeSMatthew G. Knepley PetscBool isAffine, isCohesive, isCohesiveIn, transform;
37079f2dbaeSMatthew G. Knepley DMPolytopeType qct;
371b18799e0SMatthew G. Knepley
372b18799e0SMatthew G. Knepley PetscFunctionBeginHot;
3739566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
3749566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc));
37579f2dbaeSMatthew G. Knepley PetscCall(PetscDSIsCohesive(ds, &isCohesive));
37679f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
37779f2dbaeSMatthew G. Knepley PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
37879f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff));
37979f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x));
38079f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x));
38179f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL));
38279f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants));
38379f2dbaeSMatthew G. Knepley PetscCall(DMHasBasisTransform(dmIn, &transform));
38479f2dbaeSMatthew G. Knepley PetscCall(DMGetLocalSection(dmIn, §ion));
38579f2dbaeSMatthew G. Knepley PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp));
38679f2dbaeSMatthew G. Knepley // Get cohesive cell hanging off face
38779f2dbaeSMatthew G. Knepley if (isCohesiveIn) {
38879f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetCellType(dmIn, inp, &qct));
38979f2dbaeSMatthew G. Knepley if ((qct != DM_POLYTOPE_POINT_PRISM_TENSOR) && (qct != DM_POLYTOPE_SEG_PRISM_TENSOR) && (qct != DM_POLYTOPE_TRI_PRISM_TENSOR) && (qct != DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
39079f2dbaeSMatthew G. Knepley DMPolytopeType ct;
39179f2dbaeSMatthew G. Knepley const PetscInt *support;
39279f2dbaeSMatthew G. Knepley PetscInt Ns, s;
39379f2dbaeSMatthew G. Knepley
39479f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetSupport(dmIn, inp, &support));
39579f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns));
39679f2dbaeSMatthew G. Knepley for (s = 0; s < Ns; ++s) {
39779f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetCellType(dmIn, support[s], &ct));
39879f2dbaeSMatthew G. Knepley if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
39979f2dbaeSMatthew G. Knepley inp = support[s];
40079f2dbaeSMatthew G. Knepley break;
40179f2dbaeSMatthew G. Knepley }
40279f2dbaeSMatthew G. Knepley }
40379f2dbaeSMatthew G. Knepley PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp);
40479f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff));
40579f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt));
40679f2dbaeSMatthew G. Knepley face[0] = 0;
40779f2dbaeSMatthew G. Knepley face[1] = 0;
40879f2dbaeSMatthew G. Knepley }
40979f2dbaeSMatthew G. Knepley }
41079f2dbaeSMatthew G. Knepley if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients));
411b18799e0SMatthew G. Knepley if (dmAux) {
412a6e0b375SMatthew G. Knepley PetscInt subp;
413b18799e0SMatthew G. Knepley
4149566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp));
4159566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4169566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmAux, §ionAux));
4179566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4189566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4199566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x));
4209566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux));
421b18799e0SMatthew G. Knepley }
422b18799e0SMatthew G. Knepley /* Get values for closure */
4234bee2e38SMatthew G. Knepley isAffine = fgeom->isAffine;
424989fa639SBrad Aagaard fegeom.dim = fgeom->dim;
425989fa639SBrad Aagaard fegeom.dimEmbed = fgeom->dimEmbed;
426ea78f98cSLisandro Dalcin fegeom.n = NULL;
427ea78f98cSLisandro Dalcin fegeom.J = NULL;
428ea78f98cSLisandro Dalcin fegeom.v = NULL;
429ea78f98cSLisandro Dalcin fegeom.xi = NULL;
430a6e0b375SMatthew G. Knepley cgeom.dim = fgeom->dim;
431a6e0b375SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed;
432989fa639SBrad Aagaard if (isCohesiveIn) {
433989fa639SBrad Aagaard fgeomN[0].dim = fgeom->dim;
434989fa639SBrad Aagaard fgeomN[0].dimEmbed = fgeom->dimEmbed;
435989fa639SBrad Aagaard fgeomN[1].dim = fgeom->dim;
436989fa639SBrad Aagaard fgeomN[1].dimEmbed = fgeom->dimEmbed;
437989fa639SBrad Aagaard }
4384bee2e38SMatthew G. Knepley if (isAffine) {
4394bee2e38SMatthew G. Knepley fegeom.v = x;
4404bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi;
4414bee2e38SMatthew G. Knepley fegeom.J = fgeom->J;
4424bee2e38SMatthew G. Knepley fegeom.invJ = fgeom->invJ;
4434bee2e38SMatthew G. Knepley fegeom.detJ = fgeom->detJ;
4444bee2e38SMatthew G. Knepley fegeom.n = fgeom->n;
4459f209ee4SMatthew G. Knepley
4469f209ee4SMatthew G. Knepley cgeom.J = fgeom->suppJ[0];
4479f209ee4SMatthew G. Knepley cgeom.invJ = fgeom->suppInvJ[0];
4489f209ee4SMatthew G. Knepley cgeom.detJ = fgeom->suppDetJ[0];
449989fa639SBrad Aagaard
450989fa639SBrad Aagaard if (isCohesiveIn) {
451989fa639SBrad Aagaard fgeomN[0].J = fgeom->suppJ[0];
452989fa639SBrad Aagaard fgeomN[0].invJ = fgeom->suppInvJ[0];
453989fa639SBrad Aagaard fgeomN[0].detJ = fgeom->suppDetJ[0];
454989fa639SBrad Aagaard fgeomN[1].J = fgeom->suppJ[1];
455989fa639SBrad Aagaard fgeomN[1].invJ = fgeom->suppInvJ[1];
456989fa639SBrad Aagaard fgeomN[1].detJ = fgeom->suppDetJ[1];
457989fa639SBrad Aagaard }
4584bee2e38SMatthew G. Knepley }
459b18799e0SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) {
460b18799e0SMatthew G. Knepley PetscQuadrature allPoints;
461b18799e0SMatthew G. Knepley PetscInt q, dim, numPoints;
462b18799e0SMatthew G. Knepley const PetscReal *points;
463b18799e0SMatthew G. Knepley PetscScalar *pointEval;
46479f2dbaeSMatthew G. Knepley PetscBool cohesive;
465b18799e0SMatthew G. Knepley DM dm;
466b18799e0SMatthew G. Knepley
467b18799e0SMatthew G. Knepley if (!sp[f]) continue;
46879f2dbaeSMatthew G. Knepley PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
4699566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
470b18799e0SMatthew G. Knepley if (!funcs[f]) {
471b18799e0SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.;
47279f2dbaeSMatthew G. Knepley if (isCohesive && !cohesive) {
47379f2dbaeSMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.;
47479f2dbaeSMatthew G. Knepley }
475b18799e0SMatthew G. Knepley continue;
476b18799e0SMatthew G. Knepley }
4779566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &dm));
4789566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL));
4799566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL));
4809566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
481b18799e0SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) {
48279f2dbaeSMatthew G. Knepley PetscInt qpt[2];
48379f2dbaeSMatthew G. Knepley
48479f2dbaeSMatthew G. Knepley if (isCohesiveIn) {
48579f2dbaeSMatthew G. Knepley // These points are not integration quadratures, but dual space quadratures
486d8b4a066SPierre Jolivet // If they had multiple points we should match them from both sides, similar to hybrid residual eval
48779f2dbaeSMatthew G. Knepley qpt[0] = qpt[1] = q;
48879f2dbaeSMatthew G. Knepley }
489b18799e0SMatthew G. Knepley if (isAffine) {
490ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x);
491b18799e0SMatthew G. Knepley } else {
4923fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[tp * dE];
4939f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[tp * dE * dE];
4949f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[tp * dE * dE];
4954bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[tp];
4964bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[tp * dE];
4979f209ee4SMatthew G. Knepley
4989f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][tp * dE * dE];
4999f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE];
5009f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][tp];
501989fa639SBrad Aagaard if (isCohesiveIn) {
502989fa639SBrad Aagaard fgeomN[0].J = &fgeom->suppJ[0][tp * dE * dE];
503989fa639SBrad Aagaard fgeomN[0].invJ = &fgeom->suppInvJ[0][tp * dE * dE];
504989fa639SBrad Aagaard fgeomN[0].detJ = &fgeom->suppDetJ[0][tp];
505989fa639SBrad Aagaard fgeomN[1].J = &fgeom->suppJ[1][tp * dE * dE];
506989fa639SBrad Aagaard fgeomN[1].invJ = &fgeom->suppInvJ[1][tp * dE * dE];
507989fa639SBrad Aagaard fgeomN[1].detJ = &fgeom->suppDetJ[1][tp];
508989fa639SBrad Aagaard }
509b18799e0SMatthew G. Knepley }
510a6e0b375SMatthew 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 */
51179f2dbaeSMatthew G. Knepley if (coefficients) {
512989fa639SBrad Aagaard if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, fgeomN, coefficients, coefficients_t, u, u_x, u_t));
51379f2dbaeSMatthew G. Knepley else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t));
51479f2dbaeSMatthew G. Knepley }
5159566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t));
51679f2dbaeSMatthew G. Knepley if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx));
51779f2dbaeSMatthew 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, fegeom.n, numConstants, constants, &pointEval[Nc[f] * q]);
518b18799e0SMatthew G. Knepley }
5199566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]));
5209566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval));
521b18799e0SMatthew G. Knepley v += spDim;
52279f2dbaeSMatthew G. Knepley /* TODO: For now, set both sides equal, but this should use info from other support cell */
52379f2dbaeSMatthew G. Knepley if (isCohesive && !cohesive) {
52479f2dbaeSMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim];
525b18799e0SMatthew G. Knepley }
52679f2dbaeSMatthew G. Knepley }
52779f2dbaeSMatthew G. Knepley if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients));
5289566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux));
52979f2dbaeSMatthew G. Knepley if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt));
5303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
531b18799e0SMatthew G. Knepley }
532b18799e0SMatthew G. Knepley
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,PetscVoidFn ** funcs,void ** ctxs,PetscBool fieldActive[],PetscScalar values[])533*57d50842SBarry Smith 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, PetscVoidFn **funcs, void **ctxs, PetscBool fieldActive[], PetscScalar values[])
534d71ae5a4SJacob Faibussowitsch {
5358c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom;
5368c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed;
5378c6c5593SMatthew G. Knepley
5388c6c5593SMatthew G. Knepley PetscFunctionBeginHot;
5399566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
5409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
5419566063dSJacob Faibussowitsch if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL));
5428c6c5593SMatthew G. Knepley switch (type) {
5438c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL:
544d71ae5a4SJacob Faibussowitsch case DM_BC_NATURAL:
545d71ae5a4SJacob 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));
546d71ae5a4SJacob Faibussowitsch break;
5478c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD:
5488c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD:
5499371c9d4SSatish 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));
550d0609cedSBarry Smith break;
551b18799e0SMatthew G. Knepley case DM_BC_ESSENTIAL_BD_FIELD:
55279f2dbaeSMatthew G. Knepley PetscCall(DMProjectPoint_BdField_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[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))funcs, ctxs, values));
553d0609cedSBarry Smith break;
554d71ae5a4SJacob Faibussowitsch default:
555d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type);
5568c6c5593SMatthew G. Knepley }
5573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5588c6c5593SMatthew G. Knepley }
5598c6c5593SMatthew G. Knepley
PetscDualSpaceGetAllPointsUnion(PetscInt Nf,PetscDualSpace * sp,PetscInt dim,PetscVoidFn ** funcs,PetscQuadrature * allPoints)560*57d50842SBarry Smith static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, PetscVoidFn **funcs, PetscQuadrature *allPoints)
561d71ae5a4SJacob Faibussowitsch {
562df5c1128SToby Isaac PetscReal *points;
563df5c1128SToby Isaac PetscInt f, numPoints;
564df5c1128SToby Isaac
565df5c1128SToby Isaac PetscFunctionBegin;
56619552267SMatthew G. Knepley if (!dim) {
56719552267SMatthew G. Knepley PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
56919552267SMatthew G. Knepley }
570df5c1128SToby Isaac numPoints = 0;
571df5c1128SToby Isaac for (f = 0; f < Nf; ++f) {
572df5c1128SToby Isaac if (funcs[f]) {
573df5c1128SToby Isaac PetscQuadrature fAllPoints;
574df5c1128SToby Isaac PetscInt fNumPoints;
575df5c1128SToby Isaac
5769566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5779566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL));
578df5c1128SToby Isaac numPoints += fNumPoints;
579df5c1128SToby Isaac }
580df5c1128SToby Isaac }
5819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * numPoints, &points));
582df5c1128SToby Isaac numPoints = 0;
583df5c1128SToby Isaac for (f = 0; f < Nf; ++f) {
584df5c1128SToby Isaac if (funcs[f]) {
585df5c1128SToby Isaac PetscQuadrature fAllPoints;
58654ee0cceSMatthew G. Knepley PetscInt qdim, fNumPoints, q;
587df5c1128SToby Isaac const PetscReal *fPoints;
588df5c1128SToby Isaac
5899566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL));
5909566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL));
59163a3b9bcSJacob 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);
592df5c1128SToby Isaac for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q];
593df5c1128SToby Isaac numPoints += fNumPoints;
594df5c1128SToby Isaac }
595df5c1128SToby Isaac }
5969566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints));
5979566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL));
5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
599df5c1128SToby Isaac }
600df5c1128SToby Isaac
6015f15299fSJeremy L Thompson /*@C
60220f4b53cSBarry Smith DMGetFirstLabeledPoint - Find first labeled `point` in `odm` such that the corresponding point in `dm` has the specified `height`. Return `point` and the corresponding `ds`.
6035f15299fSJeremy L Thompson
6045f15299fSJeremy L Thompson Input Parameters:
6052fe279fdSBarry Smith + dm - the `DM`
6062fe279fdSBarry Smith . odm - the enclosing `DM`
6072fe279fdSBarry Smith . label - label for `DM` domain, or `NULL` for whole domain
6082fe279fdSBarry Smith . numIds - the number of `ids`
6092fe279fdSBarry Smith . ids - An array of the label ids in sequence for the domain
6102fe279fdSBarry Smith - height - Height of target cells in `DMPLEX` topology
6115f15299fSJeremy L Thompson
6125f15299fSJeremy L Thompson Output Parameters:
6132fe279fdSBarry Smith + point - the first labeled point
614a3b724e8SBarry Smith - ds - the `PetscDS` corresponding to the first labeled point
6155f15299fSJeremy L Thompson
6165f15299fSJeremy L Thompson Level: developer
617a1cb98faSBarry Smith
6181cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS`
619817da375SSatish Balay @*/
DMGetFirstLabeledPoint(DM dm,DM odm,DMLabel label,PetscInt numIds,const PetscInt ids[],PetscInt height,PetscInt * point,PetscDS * ds)620d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds)
621d71ae5a4SJacob Faibussowitsch {
622a6e0b375SMatthew G. Knepley DM plex;
623a6e0b375SMatthew G. Knepley DMEnclosureType enc;
624e9f0ba4eSJed Brown PetscInt ls = -1;
625e5e52638SMatthew G. Knepley
626e5e52638SMatthew G. Knepley PetscFunctionBegin;
6275f15299fSJeremy L Thompson if (point) *point = -1;
6283ba16761SJacob Faibussowitsch if (!label) PetscFunctionReturn(PETSC_SUCCESS);
6299566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dm, odm, &enc));
6309566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex));
631e9f0ba4eSJed Brown for (PetscInt i = 0; i < numIds; ++i) {
632e9f0ba4eSJed Brown IS labelIS;
633e9f0ba4eSJed Brown PetscInt num_points, pStart, pEnd;
6349566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS));
635e9f0ba4eSJed Brown if (!labelIS) continue; /* No points with that id on this process */
6369566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd));
6379566063dSJacob Faibussowitsch PetscCall(ISGetSize(labelIS, &num_points));
638e9f0ba4eSJed Brown if (num_points) {
639e5e52638SMatthew G. Knepley const PetscInt *points;
6409566063dSJacob Faibussowitsch PetscCall(ISGetIndices(labelIS, &points));
641e9f0ba4eSJed Brown for (PetscInt i = 0; i < num_points; i++) {
642e9f0ba4eSJed Brown PetscInt point;
6439566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point));
644e9f0ba4eSJed Brown if (pStart <= point && point < pEnd) {
645a6e0b375SMatthew G. Knepley ls = point;
64679f2dbaeSMatthew G. Knepley if (ds) {
64779f2dbaeSMatthew G. Knepley // If this is a face of a cohesive cell, then prefer that DS
64879f2dbaeSMatthew G. Knepley if (height == 1) {
64979f2dbaeSMatthew G. Knepley const PetscInt *supp;
65079f2dbaeSMatthew G. Knepley PetscInt suppSize;
65179f2dbaeSMatthew G. Knepley DMPolytopeType ct;
65279f2dbaeSMatthew G. Knepley
65379f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetSupport(dm, ls, &supp));
65479f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetSupportSize(dm, ls, &suppSize));
65579f2dbaeSMatthew G. Knepley for (PetscInt s = 0; s < suppSize; ++s) {
65679f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, supp[s], &ct));
65779f2dbaeSMatthew G. Knepley if ((ct == DM_POLYTOPE_POINT_PRISM_TENSOR) || (ct == DM_POLYTOPE_SEG_PRISM_TENSOR) || (ct == DM_POLYTOPE_TRI_PRISM_TENSOR) || (ct == DM_POLYTOPE_QUAD_PRISM_TENSOR)) {
65879f2dbaeSMatthew G. Knepley ls = supp[s];
65979f2dbaeSMatthew G. Knepley break;
66079f2dbaeSMatthew G. Knepley }
66179f2dbaeSMatthew G. Knepley }
66279f2dbaeSMatthew G. Knepley }
66379f2dbaeSMatthew G. Knepley PetscCall(DMGetCellDS(dm, ls, ds, NULL));
66479f2dbaeSMatthew G. Knepley }
6658e6d238bSMatthew G. Knepley if (ls >= 0) break;
666e5e52638SMatthew G. Knepley }
667e9f0ba4eSJed Brown }
6689566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(labelIS, &points));
669e9f0ba4eSJed Brown }
6709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&labelIS));
671e5e52638SMatthew G. Knepley if (ls >= 0) break;
672e5e52638SMatthew G. Knepley }
6735f15299fSJeremy L Thompson if (point) *point = ls;
6749566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex));
6753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
676e5e52638SMatthew G. Knepley }
677e5e52638SMatthew G. Knepley
6780de51b56SMatthew G. Knepley /*
6790de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM
6800de51b56SMatthew G. Knepley
6810de51b56SMatthew G. Knepley There are several different scenarios:
6820de51b56SMatthew G. Knepley
6830de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data
6840de51b56SMatthew G. Knepley
6850de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells.
6860de51b56SMatthew G. Knepley
6870de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data
6880de51b56SMatthew G. Knepley
6890de51b56SMatthew 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.
6900de51b56SMatthew G. Knepley
6910de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data
6920de51b56SMatthew G. Knepley
6930de51b56SMatthew 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.
6940de51b56SMatthew G. Knepley
695636b9322SMatthew G. Knepley 4) Volumetric input mesh with boundary output mesh
696636b9322SMatthew G. Knepley
697636b9322SMatthew G. Knepley Here we must get a subspace for the input DS
698636b9322SMatthew G. Knepley
6990de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest.
7000de51b56SMatthew G. Knepley
7010de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it.
7020de51b56SMatthew G. Knepley
7030de51b56SMatthew 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.
7040de51b56SMatthew 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
7055f790a90SMatthew 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
7060de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces.
7070de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them.
7080de51b56SMatthew G. Knepley
7090de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration
7100de51b56SMatthew G. Knepley
7110de51b56SMatthew 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.
7120de51b56SMatthew G. Knepley */
DMProjectLocal_Generic_Plex(DM dm,PetscReal time,Vec localU,PetscInt Ncc,const PetscInt comps[],DMLabel label,PetscInt numIds,const PetscInt ids[],DMBoundaryConditionType type,PetscVoidFn ** funcs,void ** ctxs,InsertMode mode,Vec localX)713*57d50842SBarry Smith 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, PetscVoidFn **funcs, void **ctxs, InsertMode mode, Vec localX)
714d71ae5a4SJacob Faibussowitsch {
715a6e0b375SMatthew G. Knepley DM plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm;
716a6e0b375SMatthew G. Knepley DMEnclosureType encIn, encAux;
717a6e0b375SMatthew G. Knepley PetscDS ds = NULL, dsIn = NULL, dsAux = NULL;
718ca3d3a14SMatthew G. Knepley Vec localA = NULL, tv;
719aa0eca99SMatthew G. Knepley IS fieldIS;
72047923291SMatthew G. Knepley PetscSection section;
721636b9322SMatthew G. Knepley PetscDualSpace *sp, *cellsp, *spIn, *cellspIn;
722ef0bb6c7SMatthew G. Knepley PetscTabulation *T = NULL, *TAux = NULL;
7238c6c5593SMatthew G. Knepley PetscInt *Nc;
724989fa639SBrad Aagaard PetscInt dim, dimEmbed, depth, pStart, pEnd, lStart = PETSC_DETERMINE, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, minHeightIn, minHeightAux, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f;
7258e6d238bSMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform;
726c330f8ffSToby Isaac DMField coordField;
727c330f8ffSToby Isaac DMLabel depthLabel;
728c330f8ffSToby Isaac PetscQuadrature allPoints = NULL;
72947923291SMatthew G. Knepley
73047923291SMatthew G. Knepley PetscFunctionBegin;
7319566063dSJacob Faibussowitsch if (localU) PetscCall(VecGetDM(localU, &dmIn));
732ad540459SPierre Jolivet else dmIn = dm;
7339566063dSJacob Faibussowitsch PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA));
7349371c9d4SSatish Balay if (localA) PetscCall(VecGetDM(localA, &dmAux));
735ad540459SPierre Jolivet else dmAux = NULL;
7369566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex));
7379566063dSJacob Faibussowitsch PetscCall(DMConvert(dmIn, DMPLEX, &plexIn));
7389566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn));
7399566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux));
7409566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
7419566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight));
7429566063dSJacob Faibussowitsch PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm));
7439566063dSJacob Faibussowitsch PetscCall(DMGetBasisTransformVec_Internal(dm, &tv));
7449566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dm, &transform));
7450de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */
746a6e0b375SMatthew G. Knepley if (dmAux) {
7479566063dSJacob Faibussowitsch PetscCall(DMConvert(dmAux, DMPLEX, &plexAux));
748ad540459SPierre 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");
749636b9322SMatthew G. Knepley }
7505ee220baSJed Brown if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plexIn, PETSC_TRUE, localU, time, NULL, NULL, NULL));
751e04ae0b4SMatthew G. Knepley PetscCall(DMGetCoordinateField(dm, &coordField));
752e44f6aebSMatthew G. Knepley PetscCheck(coordField, PETSC_COMM_SELF, PETSC_ERR_USER, "DM must have a coordinate field");
753e04ae0b4SMatthew G. Knepley /**** No collective calls below this point ****/
754636b9322SMatthew G. Knepley /* Determine height for iteration of all meshes */
755636b9322SMatthew G. Knepley {
756636b9322SMatthew G. Knepley DMPolytopeType ct, ctIn, ctAux;
757989fa639SBrad Aagaard PetscInt p, pStartIn, pStartAux, pEndAux;
75819552267SMatthew G. Knepley PetscInt dim = -1, dimIn = -1, dimAux = -1;
75988aca1feSMatthew G. Knepley
7609566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd));
761636b9322SMatthew G. Knepley if (pEnd > pStart) {
7629566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL));
763636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart;
7649566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plex, p, &ct));
765636b9322SMatthew G. Knepley dim = DMPolytopeTypeGetDim(ct);
7669566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn));
7679566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL));
7689566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn));
769636b9322SMatthew G. Knepley dimIn = DMPolytopeTypeGetDim(ctIn);
770636b9322SMatthew G. Knepley if (dmAux) {
7719566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux));
77219552267SMatthew G. Knepley PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux));
77319552267SMatthew G. Knepley if (pStartAux < pEndAux) {
7749566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux));
775636b9322SMatthew G. Knepley dimAux = DMPolytopeTypeGetDim(ctAux);
77619552267SMatthew G. Knepley }
777636b9322SMatthew G. Knepley } else dimAux = dim;
778e04ae0b4SMatthew G. Knepley } else {
779e04ae0b4SMatthew G. Knepley PetscCall(DMDestroy(&plex));
780e04ae0b4SMatthew G. Knepley PetscCall(DMDestroy(&plexIn));
781e04ae0b4SMatthew G. Knepley if (dmAux) PetscCall(DMDestroy(&plexAux));
7823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
783e39fcb4eSStefano Zampini }
784636b9322SMatthew G. Knepley if (dim < 0) {
785636b9322SMatthew G. Knepley DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL;
786636b9322SMatthew G. Knepley
787636b9322SMatthew G. Knepley /* Fall back to determination based on being a submesh */
7889566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(plex, &spmap));
7899566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn));
7909566063dSJacob Faibussowitsch if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux));
791636b9322SMatthew G. Knepley dim = spmap ? 1 : 0;
792636b9322SMatthew G. Knepley dimIn = spmapIn ? 1 : 0;
793636b9322SMatthew G. Knepley dimAux = spmapAux ? 1 : 0;
79488aca1feSMatthew G. Knepley }
795636b9322SMatthew G. Knepley {
79657508eceSPierre Jolivet PetscInt dimProj = PetscMin(PetscMin(dim, dimIn), dimAux < 0 ? PETSC_INT_MAX : dimAux);
79719552267SMatthew G. Knepley PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux;
798636b9322SMatthew G. Knepley
79919552267SMatthew 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");
800636b9322SMatthew G. Knepley if (dimProj < dim) minHeight = 1;
801636b9322SMatthew G. Knepley htInc = dim - dimProj;
802636b9322SMatthew G. Knepley htIncIn = dimIn - dimProj;
80319552267SMatthew G. Knepley htIncAux = dimAuxEff - dimProj;
8040de51b56SMatthew G. Knepley }
805a6e0b375SMatthew G. Knepley }
8069566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(plex, &depth));
8079566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(plex, &depthLabel));
8089566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight));
8092af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight);
81063a3b9bcSJacob 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);
81179f2dbaeSMatthew G. Knepley PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, NULL, &ds));
8129566063dSJacob Faibussowitsch if (!ds) PetscCall(DMGetDS(dm, &ds));
81379f2dbaeSMatthew G. Knepley PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, minHeight, NULL, &dsIn));
8149566063dSJacob Faibussowitsch if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn));
815989fa639SBrad Aagaard if (!dsIn) {
816989fa639SBrad Aagaard if (encIn == DM_ENC_SUPERMESH) {
817989fa639SBrad Aagaard PetscInt p = pStart, pIn;
818989fa639SBrad Aagaard
819989fa639SBrad Aagaard PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? p : lStart, &pIn));
820989fa639SBrad Aagaard // If the input mesh is higher dimensional than the output mesh, get a cell from the output mesh
821989fa639SBrad Aagaard if (htIncIn) PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &p, NULL));
822989fa639SBrad Aagaard PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? p : lStart, &pIn));
823989fa639SBrad Aagaard PetscCall(DMGetCellDS(dmIn, pIn, &dsIn, NULL));
824989fa639SBrad Aagaard } else PetscCall(DMGetDS(dmIn, &dsIn));
825989fa639SBrad Aagaard }
8269566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf));
8279566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsIn, &NfIn));
8288e6d238bSMatthew G. Knepley PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn));
8298e6d238bSMatthew G. Knepley if (isCohesiveIn) --htIncIn; // Should be rearranged
8309566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &NfTot));
8319566063dSJacob Faibussowitsch PetscCall(DMFindRegionNum(dm, ds, ®ionNum));
83207218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL));
8339566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive));
8349566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
8359566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion));
8360c898c18SMatthew G. Knepley if (dmAux) {
837989fa639SBrad Aagaard if (encAux == DM_ENC_SUPERMESH) {
838989fa639SBrad Aagaard PetscInt p = pStart, pAux;
839989fa639SBrad Aagaard
840989fa639SBrad Aagaard // If the auxiliary mesh is higher dimensional than the output mesh, get a cell from the output mesh
841989fa639SBrad Aagaard if (htIncAux) PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &p, NULL));
842989fa639SBrad Aagaard PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, lStart < 0 ? p : lStart, &pAux));
843989fa639SBrad Aagaard PetscCall(DMGetCellDS(dmAux, pAux, &dsAux, NULL));
844989fa639SBrad Aagaard } else PetscCall(DMGetDS(dmAux, &dsAux));
8459566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
8460c898c18SMatthew G. Knepley }
8479566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc));
8489566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn));
8499566063dSJacob Faibussowitsch if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn));
8509371c9d4SSatish Balay else {
8519371c9d4SSatish Balay cellsp = sp;
8529371c9d4SSatish Balay cellspIn = spIn;
8539371c9d4SSatish Balay }
8548c6c5593SMatthew G. Knepley /* Get cell dual spaces */
85547923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
856665f567fSMatthew G. Knepley PetscDiscType disctype;
85747923291SMatthew G. Knepley
8589566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype));
859665f567fSMatthew G. Knepley if (disctype == PETSC_DISC_FE) {
860665f567fSMatthew G. Knepley PetscFE fe;
86147923291SMatthew G. Knepley
86247923291SMatthew G. Knepley isFE[f] = PETSC_TRUE;
863665f567fSMatthew G. Knepley hasFE = PETSC_TRUE;
8649566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe));
8659566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &cellsp[f]));
866665f567fSMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) {
867665f567fSMatthew G. Knepley PetscFV fv;
86847923291SMatthew G. Knepley
86947923291SMatthew G. Knepley isFE[f] = PETSC_FALSE;
870665f567fSMatthew G. Knepley hasFV = PETSC_TRUE;
8719566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv));
8729566063dSJacob Faibussowitsch PetscCall(PetscFVGetDualSpace(fv, &cellsp[f]));
873665f567fSMatthew G. Knepley } else {
874665f567fSMatthew G. Knepley isFE[f] = PETSC_FALSE;
875665f567fSMatthew G. Knepley cellsp[f] = NULL;
876665f567fSMatthew G. Knepley }
87747923291SMatthew G. Knepley }
878636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) {
879636b9322SMatthew G. Knepley PetscDiscType disctype;
880636b9322SMatthew G. Knepley
8819566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
882636b9322SMatthew G. Knepley if (disctype == PETSC_DISC_FE) {
883636b9322SMatthew G. Knepley PetscFE fe;
884636b9322SMatthew G. Knepley
8859566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe));
8869566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f]));
887636b9322SMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) {
888636b9322SMatthew G. Knepley PetscFV fv;
889636b9322SMatthew G. Knepley
8909566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv));
8919566063dSJacob Faibussowitsch PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f]));
892636b9322SMatthew G. Knepley } else {
893636b9322SMatthew G. Knepley cellspIn[f] = NULL;
894636b9322SMatthew G. Knepley }
895636b9322SMatthew G. Knepley }
896636b9322SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
8979371c9d4SSatish Balay if (!htInc) {
8989371c9d4SSatish Balay sp[f] = cellsp[f];
8999371c9d4SSatish Balay } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]));
900636b9322SMatthew G. Knepley }
901ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
90274fc349aSMatthew G. Knepley PetscFE fem, subfem;
903665f567fSMatthew G. Knepley PetscDiscType disctype;
9044a3e9fdbSToby Isaac const PetscReal *points;
905409bb69aSDarsh Nathawani PetscInt numPoints, k;
9060c898c18SMatthew G. Knepley
90708401ef6SPierre Jolivet PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation");
9089566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints));
9099566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL));
9109566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux));
911a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) {
9129371c9d4SSatish Balay if (!htIncIn) {
9139371c9d4SSatish Balay spIn[f] = cellspIn[f];
9149371c9d4SSatish Balay } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]));
915636b9322SMatthew G. Knepley
9169566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype));
917665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue;
9189566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem));
9199371c9d4SSatish Balay if (!htIncIn) {
9209371c9d4SSatish Balay subfem = fem;
9219371c9d4SSatish Balay } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem));
922409bb69aSDarsh Nathawani PetscCall(PetscDSGetJetDegree(dsIn, f, &k));
923409bb69aSDarsh Nathawani PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &T[f]));
9240c898c18SMatthew G. Knepley }
9250c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) {
9269566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype));
927665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue;
9289566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem));
9299371c9d4SSatish Balay if (!htIncAux) {
9309371c9d4SSatish Balay subfem = fem;
9319371c9d4SSatish Balay } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem));
932409bb69aSDarsh Nathawani PetscCall(PetscDSGetJetDegree(dsAux, f, &k));
933409bb69aSDarsh Nathawani PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, k, &TAux[f]));
9340c898c18SMatthew G. Knepley }
9350c898c18SMatthew G. Knepley }
93647923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */
9372af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) {
938636b9322SMatthew G. Knepley PetscInt hEff = h - minHeight + htInc;
939636b9322SMatthew G. Knepley PetscInt hEffIn = h - minHeight + htIncIn;
940636b9322SMatthew G. Knepley PetscInt hEffAux = h - minHeight + htIncAux;
941a6e0b375SMatthew G. Knepley PetscDS dsEff = ds;
942636b9322SMatthew G. Knepley PetscDS dsEffIn = dsIn;
943636b9322SMatthew G. Knepley PetscDS dsEffAux = dsAux;
9448c6c5593SMatthew G. Knepley PetscScalar *values;
945b7260050SToby Isaac PetscBool *fieldActive;
946b7260050SToby Isaac PetscInt maxDegree;
947989fa639SBrad Aagaard PetscInt p, spDim, totDim, numValues;
948c330f8ffSToby Isaac IS heightIS;
9498c6c5593SMatthew G. Knepley
950636b9322SMatthew G. Knepley if (h > minHeight) {
9519566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]));
952636b9322SMatthew G. Knepley }
9539566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd));
9549566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL));
9559566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS));
956c330f8ffSToby Isaac if (pEnd <= pStart) {
9579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS));
958c330f8ffSToby Isaac continue;
959c330f8ffSToby Isaac }
96047923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */
96147923291SMatthew G. Knepley totDim = 0;
96247923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) {
9635fedec97SMatthew G. Knepley PetscBool cohesive;
9645fedec97SMatthew G. Knepley
965665f567fSMatthew G. Knepley if (!sp[f]) continue;
9669566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive));
9679566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim));
9689c3cf19fSMatthew G. Knepley totDim += spDim;
9695fedec97SMatthew G. Knepley if (isCohesive && !cohesive) totDim += spDim;
97047923291SMatthew G. Knepley }
971636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart;
9729566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL));
97363a3b9bcSJacob 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);
974c330f8ffSToby Isaac if (!totDim) {
9759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS));
976c330f8ffSToby Isaac continue;
977c330f8ffSToby Isaac }
9789566063dSJacob Faibussowitsch if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff));
979636b9322SMatthew G. Knepley /* Compute totDimIn, the number of dofs in the closure of a point at this height */
980636b9322SMatthew G. Knepley if (localU) {
981636b9322SMatthew G. Knepley PetscInt totDimIn, pIn, numValuesIn;
982636b9322SMatthew G. Knepley
983636b9322SMatthew G. Knepley totDimIn = 0;
984636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) {
9855fedec97SMatthew G. Knepley PetscBool cohesive;
9865fedec97SMatthew G. Knepley
987636b9322SMatthew G. Knepley if (!spIn[f]) continue;
9889566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive));
9899566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim));
990636b9322SMatthew G. Knepley totDimIn += spDim;
9918e6d238bSMatthew G. Knepley if (isCohesiveIn && !cohesive) totDimIn += spDim;
992636b9322SMatthew G. Knepley }
9939566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn));
9949566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL));
9958e6d238bSMatthew G. Knepley // TODO We could check that pIn is a cohesive cell for this check
9968e6d238bSMatthew G. Knepley PetscCheck(isCohesiveIn || (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);
9979566063dSJacob Faibussowitsch if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn));
998636b9322SMatthew G. Knepley }
9999566063dSJacob Faibussowitsch if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux));
100047923291SMatthew G. Knepley /* Loop over points at this height */
10019566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values));
10029566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive));
1003aa0eca99SMatthew G. Knepley {
1004aa0eca99SMatthew G. Knepley const PetscInt *fields;
1005aa0eca99SMatthew G. Knepley
10069566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fieldIS, &fields));
1007ad540459SPierre Jolivet for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE;
1008ad540459SPierre Jolivet for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;
10099566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fieldIS, &fields));
1010aa0eca99SMatthew G. Knepley }
10118c6c5593SMatthew G. Knepley if (label) {
10128c6c5593SMatthew G. Knepley PetscInt i;
101347923291SMatthew G. Knepley
101447923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) {
1015c330f8ffSToby Isaac IS pointIS, isectIS;
101647923291SMatthew G. Knepley const PetscInt *points;
10178c6c5593SMatthew G. Knepley PetscInt n;
1018c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL;
1019c330f8ffSToby Isaac PetscQuadrature quad = NULL;
102047923291SMatthew G. Knepley
10219566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS));
102247923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */
10239566063dSJacob Faibussowitsch PetscCall(ISIntersect(pointIS, heightIS, &isectIS));
10249566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS));
1025c330f8ffSToby Isaac if (!isectIS) continue;
10269566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isectIS, &n));
10279566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isectIS, &points));
10289566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree));
102948a46eb9SPierre Jolivet if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad));
1030c330f8ffSToby Isaac if (!quad) {
1031c330f8ffSToby Isaac if (!h && allPoints) {
1032c330f8ffSToby Isaac quad = allPoints;
1033c330f8ffSToby Isaac allPoints = NULL;
1034c330f8ffSToby Isaac } else {
10359566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad));
1036c330f8ffSToby Isaac }
1037c330f8ffSToby Isaac }
1038ac9d17c7SMatthew G. Knepley PetscFEGeomMode geommode = htInc && h == minHeight ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC;
103979f2dbaeSMatthew G. Knepley
104079f2dbaeSMatthew G. Knepley if (n) {
104179f2dbaeSMatthew G. Knepley PetscInt depth, dep;
104279f2dbaeSMatthew G. Knepley
104379f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetDepth(dm, &depth));
104479f2dbaeSMatthew G. Knepley PetscCall(DMPlexGetPointDepth(dm, points[0], &dep));
1045ac9d17c7SMatthew G. Knepley if (dep < depth && h == minHeight) geommode = PETSC_FEGEOM_BOUNDARY;
104679f2dbaeSMatthew G. Knepley }
1047ac9d17c7SMatthew G. Knepley PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, geommode, &fegeom));
104847923291SMatthew G. Knepley for (p = 0; p < n; ++p) {
104947923291SMatthew G. Knepley const PetscInt point = points[p];
105047923291SMatthew G. Knepley
10519566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, numValues));
10529566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom));
10539566063dSJacob Faibussowitsch PetscCall(DMPlexSetActivePoint(dm, point));
1054d0609cedSBarry 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));
10559566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values));
10569566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode));
105747923291SMatthew G. Knepley }
10589566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom));
10599566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&fegeom));
10609566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad));
10619566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isectIS, &points));
10629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isectIS));
106347923291SMatthew G. Knepley }
10648c6c5593SMatthew G. Knepley } else {
1065c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL;
1066c330f8ffSToby Isaac PetscQuadrature quad = NULL;
1067c330f8ffSToby Isaac IS pointIS;
1068c330f8ffSToby Isaac
10699566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS));
10709566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree));
107148a46eb9SPierre Jolivet if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad));
1072c330f8ffSToby Isaac if (!quad) {
1073c330f8ffSToby Isaac if (!h && allPoints) {
1074c330f8ffSToby Isaac quad = allPoints;
1075c330f8ffSToby Isaac allPoints = NULL;
1076c330f8ffSToby Isaac } else {
10779566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad));
1078c330f8ffSToby Isaac }
1079c330f8ffSToby Isaac }
1080ac9d17c7SMatthew G. Knepley PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_FEGEOM_BOUNDARY : PETSC_FEGEOM_BASIC, &fegeom));
10818c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {
10829566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, numValues));
10839566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom));
10849566063dSJacob Faibussowitsch PetscCall(DMPlexSetActivePoint(dm, p));
1085d0609cedSBarry 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));
10869566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values));
10879566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode));
10888c6c5593SMatthew G. Knepley }
10899566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom));
10909566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&fegeom));
10919566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad));
10929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS));
10938c6c5593SMatthew G. Knepley }
10949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS));
10959566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values));
10969566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive));
109747923291SMatthew G. Knepley }
10988c6c5593SMatthew G. Knepley /* Cleanup */
1099ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) {
11009566063dSJacob Faibussowitsch for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f]));
11019566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f]));
11029566063dSJacob Faibussowitsch PetscCall(PetscFree2(T, TAux));
11030c898c18SMatthew G. Knepley }
11049566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&allPoints));
11059566063dSJacob Faibussowitsch PetscCall(PetscFree3(isFE, sp, spIn));
11069566063dSJacob Faibussowitsch if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn));
11079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex));
11089566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plexIn));
11099566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMDestroy(&plexAux));
11103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
111147923291SMatthew G. Knepley }
11128c6c5593SMatthew G. Knepley
DMProjectFunctionLocal_Plex(DM dm,PetscReal time,PetscErrorCode (** funcs)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void ** ctxs,InsertMode mode,Vec localX)1113d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
1114d71ae5a4SJacob Faibussowitsch {
11158c6c5593SMatthew G. Knepley PetscFunctionBegin;
1116*57d50842SBarry Smith PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (PetscVoidFn **)funcs, ctxs, mode, localX));
11173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11188c6c5593SMatthew G. Knepley }
11198c6c5593SMatthew G. Knepley
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)1120d71ae5a4SJacob 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)
1121d71ae5a4SJacob Faibussowitsch {
11228c6c5593SMatthew G. Knepley PetscFunctionBegin;
1123*57d50842SBarry Smith PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (PetscVoidFn **)funcs, ctxs, mode, localX));
11243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
112547923291SMatthew G. Knepley }
112647923291SMatthew G. Knepley
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)1127d71ae5a4SJacob 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)
1128d71ae5a4SJacob Faibussowitsch {
112947923291SMatthew G. Knepley PetscFunctionBegin;
1130*57d50842SBarry Smith PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (PetscVoidFn **)funcs, NULL, mode, localX));
11313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
113247923291SMatthew G. Knepley }
113347923291SMatthew G. Knepley
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)1134d71ae5a4SJacob 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)
1135d71ae5a4SJacob Faibussowitsch {
11368c6c5593SMatthew G. Knepley PetscFunctionBegin;
1137*57d50842SBarry Smith PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (PetscVoidFn **)funcs, NULL, mode, localX));
11383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
113947923291SMatthew G. Knepley }
1140ece3a9fcSMatthew G. Knepley
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)1141d71ae5a4SJacob 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)
1142d71ae5a4SJacob Faibussowitsch {
1143ece3a9fcSMatthew G. Knepley PetscFunctionBegin;
1144*57d50842SBarry Smith PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (PetscVoidFn **)funcs, NULL, mode, localX));
11453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1146ece3a9fcSMatthew G. Knepley }
1147