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 @*/ 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 @*/ 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 */ 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 { 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 } 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 Level: developer 177a6e0b375SMatthew G. Knepley 178a1cb98faSBarry Smith Note: 179a1cb98faSBarry Smith Not supported for FV 180a1cb98faSBarry Smith 181db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()` 182a6e0b375SMatthew G. Knepley */ 183d71ae5a4SJacob 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[]) 184d71ae5a4SJacob Faibussowitsch { 1858c6c5593SMatthew G. Knepley PetscSection section, sectionAux = NULL; 1864bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 1878c6c5593SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 1888c6c5593SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 189191494d9SMatthew G. Knepley const PetscScalar *constants; 1908c6c5593SMatthew G. Knepley PetscReal *x; 1918e6d238bSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, face[2]; 1924bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 1938e6d238bSMatthew G. Knepley const PetscInt dE = cgeom->dimEmbed, *cone, *ornt; 194ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0; 1958e6d238bSMatthew G. Knepley PetscBool isAffine, isCohesive, isCohesiveIn, transform; 1968e6d238bSMatthew G. Knepley DMPolytopeType qct; 1978c6c5593SMatthew G. Knepley 1988c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 1999566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2009566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc)); 2019566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 2029566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsIn, &NfIn)); 2038e6d238bSMatthew G. Knepley PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn)); 2049566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff)); 2059566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x)); 2069566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x)); 2079566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL)); 2089566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants)); 2099566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dmIn, &transform)); 2109566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmIn, §ion)); 2119566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp)); 2128e6d238bSMatthew G. Knepley // Get cohesive cell hanging off face 2138e6d238bSMatthew G. Knepley if (isCohesiveIn) { 2148e6d238bSMatthew G. Knepley PetscCall(DMPlexGetCellType(dmIn, inp, &qct)); 2158e6d238bSMatthew 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)) { 2168e6d238bSMatthew G. Knepley DMPolytopeType ct; 2178e6d238bSMatthew G. Knepley const PetscInt *support; 2188e6d238bSMatthew G. Knepley PetscInt Ns, s; 2198e6d238bSMatthew G. Knepley 2208e6d238bSMatthew G. Knepley PetscCall(DMPlexGetSupport(dmIn, inp, &support)); 2218e6d238bSMatthew G. Knepley PetscCall(DMPlexGetSupportSize(dmIn, inp, &Ns)); 2228e6d238bSMatthew G. Knepley for (s = 0; s < Ns; ++s) { 2238e6d238bSMatthew G. Knepley PetscCall(DMPlexGetCellType(dmIn, support[s], &ct)); 2248e6d238bSMatthew 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)) { 2258e6d238bSMatthew G. Knepley inp = support[s]; 2268e6d238bSMatthew G. Knepley break; 2278e6d238bSMatthew G. Knepley } 2288e6d238bSMatthew G. Knepley } 2298e6d238bSMatthew G. Knepley PetscCheck(s < Ns, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cohesive cell not found from face %" PetscInt_FMT, inp); 2308e6d238bSMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 2, &uOff)); 2318e6d238bSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dmIn, inp, &cone, &ornt)); 2328e6d238bSMatthew G. Knepley face[0] = 0; 2338e6d238bSMatthew G. Knepley face[1] = 0; 2348e6d238bSMatthew G. Knepley } 2358e6d238bSMatthew G. Knepley } 236eb8f539aSJed Brown if (localU) PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients)); 2378c6c5593SMatthew G. Knepley if (dmAux) { 23844171101SMatthew G. Knepley PetscInt subp; 2391ba34526SMatthew G. Knepley 2409566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp)); 2419566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 2429566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmAux, §ionAux)); 2439566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 2449566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 2459566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x)); 2469566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux)); 2478c6c5593SMatthew G. Knepley } 2488c6c5593SMatthew G. Knepley /* Get values for closure */ 2494bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 25027f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 25127f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 2524bee2e38SMatthew G. Knepley if (isAffine) { 2534bee2e38SMatthew G. Knepley fegeom.v = x; 2544bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2554bee2e38SMatthew G. Knepley fegeom.J = cgeom->J; 2564bee2e38SMatthew G. Knepley fegeom.invJ = cgeom->invJ; 2574bee2e38SMatthew G. Knepley fegeom.detJ = cgeom->detJ; 2584bee2e38SMatthew G. Knepley } 2598c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 260c330f8ffSToby Isaac PetscQuadrature allPoints; 261c330f8ffSToby Isaac PetscInt q, dim, numPoints; 262c330f8ffSToby Isaac const PetscReal *points; 263c330f8ffSToby Isaac PetscScalar *pointEval; 2645fedec97SMatthew G. Knepley PetscBool cohesive; 265c330f8ffSToby Isaac DM dm; 266c330f8ffSToby Isaac 2678c6c5593SMatthew G. Knepley if (!sp[f]) continue; 2689566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 2699566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 270c330f8ffSToby Isaac if (!funcs[f]) { 271be1504a2SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 2725fedec97SMatthew G. Knepley if (isCohesive && !cohesive) { 27327f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 27427f02ce8SMatthew G. Knepley } 275c330f8ffSToby Isaac continue; 276c330f8ffSToby Isaac } 2779566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &dm)); 2789566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 2799566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 2809566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 2810c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 2828e6d238bSMatthew G. Knepley PetscInt qpt[2]; 2838e6d238bSMatthew G. Knepley 2848e6d238bSMatthew G. Knepley if (isCohesiveIn) { 2858e6d238bSMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(dsIn, ornt[0], f, q, &qpt[0])); 2868e6d238bSMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(dsIn, DMPolytopeTypeComposeOrientationInv(qct, ornt[1], 0), f, q, &qpt[1])); 2878e6d238bSMatthew G. Knepley } 288c330f8ffSToby Isaac if (isAffine) { 289ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x); 2901c531cf8SMatthew G. Knepley } else { 2914bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[tp * dE]; 2924bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[tp * dE * dE]; 2934bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[tp * dE * dE]; 2944bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[tp]; 2958c6c5593SMatthew G. Knepley } 2968e6d238bSMatthew G. Knepley if (coefficients) { 2978e6d238bSMatthew G. Knepley if (isCohesiveIn) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, NfIn, 0, tp, T, face, qpt, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t)); 2988e6d238bSMatthew G. Knepley else PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t)); 2998e6d238bSMatthew G. Knepley } 3009566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t)); 3019566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx)); 302a6e0b375SMatthew 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]); 3031c531cf8SMatthew G. Knepley } 3049566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 3059566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 306c330f8ffSToby Isaac v += spDim; 30727f02ce8SMatthew G. Knepley /* TODO: For now, set both sides equal, but this should use info from other support cell */ 3085fedec97SMatthew G. Knepley if (isCohesive && !cohesive) { 30927f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 31027f02ce8SMatthew G. Knepley } 3118c6c5593SMatthew G. Knepley } 312eb8f539aSJed Brown if (localU) PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients)); 3139566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux)); 3148e6d238bSMatthew G. Knepley if (isCohesiveIn) PetscCall(DMPlexRestoreOrientedCone(dmIn, inp, &cone, &ornt)); 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3168c6c5593SMatthew G. Knepley } 3178c6c5593SMatthew G. Knepley 318d71ae5a4SJacob 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[]) 319d71ae5a4SJacob Faibussowitsch { 320b18799e0SMatthew G. Knepley PetscSection section, sectionAux = NULL; 3214bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 322b18799e0SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 323b18799e0SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 324b18799e0SMatthew G. Knepley const PetscScalar *constants; 325b18799e0SMatthew G. Knepley PetscReal *x; 326ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc; 3279f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 3284bee2e38SMatthew G. Knepley const PetscInt dE = fgeom->dimEmbed; 329ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 330b18799e0SMatthew G. Knepley PetscBool isAffine; 331b18799e0SMatthew G. Knepley 332b18799e0SMatthew G. Knepley PetscFunctionBeginHot; 33308401ef6SPierre Jolivet PetscCheck(dm == dmIn, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM"); 3349566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 3359566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc)); 3369566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 3379566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 3389566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x)); 3399566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL)); 3409566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 3419566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 3429566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients)); 343b18799e0SMatthew G. Knepley if (dmAux) { 344a6e0b375SMatthew G. Knepley PetscInt subp; 345b18799e0SMatthew G. Knepley 3469566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp)); 3479566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 3489566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmAux, §ionAux)); 3499566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 3509566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 3519566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x)); 3529566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux)); 353b18799e0SMatthew G. Knepley } 354b18799e0SMatthew G. Knepley /* Get values for closure */ 3554bee2e38SMatthew G. Knepley isAffine = fgeom->isAffine; 356ea78f98cSLisandro Dalcin fegeom.n = NULL; 357ea78f98cSLisandro Dalcin fegeom.J = NULL; 358ea78f98cSLisandro Dalcin fegeom.v = NULL; 359ea78f98cSLisandro Dalcin fegeom.xi = NULL; 360a6e0b375SMatthew G. Knepley cgeom.dim = fgeom->dim; 361a6e0b375SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3624bee2e38SMatthew G. Knepley if (isAffine) { 3634bee2e38SMatthew G. Knepley fegeom.v = x; 3644bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3654bee2e38SMatthew G. Knepley fegeom.J = fgeom->J; 3664bee2e38SMatthew G. Knepley fegeom.invJ = fgeom->invJ; 3674bee2e38SMatthew G. Knepley fegeom.detJ = fgeom->detJ; 3684bee2e38SMatthew G. Knepley fegeom.n = fgeom->n; 3699f209ee4SMatthew G. Knepley 3709f209ee4SMatthew G. Knepley cgeom.J = fgeom->suppJ[0]; 3719f209ee4SMatthew G. Knepley cgeom.invJ = fgeom->suppInvJ[0]; 3729f209ee4SMatthew G. Knepley cgeom.detJ = fgeom->suppDetJ[0]; 3734bee2e38SMatthew G. Knepley } 374b18799e0SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 375b18799e0SMatthew G. Knepley PetscQuadrature allPoints; 376b18799e0SMatthew G. Knepley PetscInt q, dim, numPoints; 377b18799e0SMatthew G. Knepley const PetscReal *points; 378b18799e0SMatthew G. Knepley PetscScalar *pointEval; 379b18799e0SMatthew G. Knepley DM dm; 380b18799e0SMatthew G. Knepley 381b18799e0SMatthew G. Knepley if (!sp[f]) continue; 3829566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 383b18799e0SMatthew G. Knepley if (!funcs[f]) { 384b18799e0SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 385b18799e0SMatthew G. Knepley continue; 386b18799e0SMatthew G. Knepley } 3879566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &dm)); 3889566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 3899566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 3909566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 391b18799e0SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 392b18799e0SMatthew G. Knepley if (isAffine) { 393ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x); 394b18799e0SMatthew G. Knepley } else { 3953fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[tp * dE]; 3969f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[tp * dE * dE]; 3979f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[tp * dE * dE]; 3984bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[tp]; 3994bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[tp * dE]; 4009f209ee4SMatthew G. Knepley 4019f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][tp * dE * dE]; 4029f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE]; 4039f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][tp]; 404b18799e0SMatthew G. Knepley } 405a6e0b375SMatthew 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 */ 4069566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t)); 4079566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t)); 4084bee2e38SMatthew 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]); 409b18799e0SMatthew G. Knepley } 4109566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 4119566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 412b18799e0SMatthew G. Knepley v += spDim; 413b18799e0SMatthew G. Knepley } 4149566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients)); 4159566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux)); 4163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 417b18799e0SMatthew G. Knepley } 418b18799e0SMatthew G. Knepley 419d71ae5a4SJacob 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[]) 420d71ae5a4SJacob Faibussowitsch { 4218c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom; 4228c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed; 4238c6c5593SMatthew G. Knepley 4248c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 4259566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4269566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 4279566063dSJacob Faibussowitsch if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL)); 4288c6c5593SMatthew G. Knepley switch (type) { 4298c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL: 430d71ae5a4SJacob Faibussowitsch case DM_BC_NATURAL: 431d71ae5a4SJacob 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)); 432d71ae5a4SJacob Faibussowitsch break; 4338c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 4348c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD: 4359371c9d4SSatish 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)); 436d0609cedSBarry Smith break; 437b18799e0SMatthew G. Knepley case DM_BC_ESSENTIAL_BD_FIELD: 4389371c9d4SSatish 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)); 439d0609cedSBarry Smith break; 440d71ae5a4SJacob Faibussowitsch default: 441d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type); 4428c6c5593SMatthew G. Knepley } 4433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4448c6c5593SMatthew G. Knepley } 4458c6c5593SMatthew G. Knepley 446d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 447d71ae5a4SJacob Faibussowitsch { 448df5c1128SToby Isaac PetscReal *points; 449df5c1128SToby Isaac PetscInt f, numPoints; 450df5c1128SToby Isaac 451df5c1128SToby Isaac PetscFunctionBegin; 45219552267SMatthew G. Knepley if (!dim) { 45319552267SMatthew G. Knepley PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints)); 4543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45519552267SMatthew G. Knepley } 456df5c1128SToby Isaac numPoints = 0; 457df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 458df5c1128SToby Isaac if (funcs[f]) { 459df5c1128SToby Isaac PetscQuadrature fAllPoints; 460df5c1128SToby Isaac PetscInt fNumPoints; 461df5c1128SToby Isaac 4629566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL)); 4639566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL)); 464df5c1128SToby Isaac numPoints += fNumPoints; 465df5c1128SToby Isaac } 466df5c1128SToby Isaac } 4679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * numPoints, &points)); 468df5c1128SToby Isaac numPoints = 0; 469df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 470df5c1128SToby Isaac if (funcs[f]) { 471df5c1128SToby Isaac PetscQuadrature fAllPoints; 47254ee0cceSMatthew G. Knepley PetscInt qdim, fNumPoints, q; 473df5c1128SToby Isaac const PetscReal *fPoints; 474df5c1128SToby Isaac 4759566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL)); 4769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL)); 47763a3b9bcSJacob 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); 478df5c1128SToby Isaac for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q]; 479df5c1128SToby Isaac numPoints += fNumPoints; 480df5c1128SToby Isaac } 481df5c1128SToby Isaac } 4829566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints)); 4839566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL)); 4843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 485df5c1128SToby Isaac } 486df5c1128SToby Isaac 4875f15299fSJeremy L Thompson /*@C 48820f4b53cSBarry 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`. 4895f15299fSJeremy L Thompson 4905f15299fSJeremy L Thompson Input Parameters: 4912fe279fdSBarry Smith + dm - the `DM` 4922fe279fdSBarry Smith . odm - the enclosing `DM` 4932fe279fdSBarry Smith . label - label for `DM` domain, or `NULL` for whole domain 4942fe279fdSBarry Smith . numIds - the number of `ids` 4952fe279fdSBarry Smith . ids - An array of the label ids in sequence for the domain 4962fe279fdSBarry Smith - height - Height of target cells in `DMPLEX` topology 4975f15299fSJeremy L Thompson 4985f15299fSJeremy L Thompson Output Parameters: 4992fe279fdSBarry Smith + point - the first labeled point 5002fe279fdSBarry Smith - ds - the ds corresponding to the first labeled point 5015f15299fSJeremy L Thompson 5025f15299fSJeremy L Thompson Level: developer 503a1cb98faSBarry Smith 5041cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexSetActivePoint()`, `DMLabel`, `PetscDS` 505817da375SSatish Balay @*/ 506d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds) 507d71ae5a4SJacob Faibussowitsch { 508a6e0b375SMatthew G. Knepley DM plex; 509a6e0b375SMatthew G. Knepley DMEnclosureType enc; 510e9f0ba4eSJed Brown PetscInt ls = -1; 511e5e52638SMatthew G. Knepley 512e5e52638SMatthew G. Knepley PetscFunctionBegin; 5135f15299fSJeremy L Thompson if (point) *point = -1; 5143ba16761SJacob Faibussowitsch if (!label) PetscFunctionReturn(PETSC_SUCCESS); 5159566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dm, odm, &enc)); 5169566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex)); 517e9f0ba4eSJed Brown for (PetscInt i = 0; i < numIds; ++i) { 518e9f0ba4eSJed Brown IS labelIS; 519e9f0ba4eSJed Brown PetscInt num_points, pStart, pEnd; 5209566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS)); 521e9f0ba4eSJed Brown if (!labelIS) continue; /* No points with that id on this process */ 5229566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd)); 5239566063dSJacob Faibussowitsch PetscCall(ISGetSize(labelIS, &num_points)); 524e9f0ba4eSJed Brown if (num_points) { 525e5e52638SMatthew G. Knepley const PetscInt *points; 5269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(labelIS, &points)); 527e9f0ba4eSJed Brown for (PetscInt i = 0; i < num_points; i++) { 528e9f0ba4eSJed Brown PetscInt point; 5299566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point)); 530e9f0ba4eSJed Brown if (pStart <= point && point < pEnd) { 531a6e0b375SMatthew G. Knepley ls = point; 53207218a29SMatthew G. Knepley if (ds) PetscCall(DMGetCellDS(dm, ls, ds, NULL)); 5338e6d238bSMatthew G. Knepley if (ls >= 0) break; 534e5e52638SMatthew G. Knepley } 535e9f0ba4eSJed Brown } 5369566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(labelIS, &points)); 537e9f0ba4eSJed Brown } 5389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&labelIS)); 539e5e52638SMatthew G. Knepley if (ls >= 0) break; 540e5e52638SMatthew G. Knepley } 5415f15299fSJeremy L Thompson if (point) *point = ls; 5429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 544e5e52638SMatthew G. Knepley } 545e5e52638SMatthew G. Knepley 5460de51b56SMatthew G. Knepley /* 5470de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 5480de51b56SMatthew G. Knepley 5490de51b56SMatthew G. Knepley There are several different scenarios: 5500de51b56SMatthew G. Knepley 5510de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data 5520de51b56SMatthew G. Knepley 5530de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells. 5540de51b56SMatthew G. Knepley 5550de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data 5560de51b56SMatthew G. Knepley 5570de51b56SMatthew 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. 5580de51b56SMatthew G. Knepley 5590de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data 5600de51b56SMatthew G. Knepley 5610de51b56SMatthew 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. 5620de51b56SMatthew G. Knepley 563636b9322SMatthew G. Knepley 4) Volumetric input mesh with boundary output mesh 564636b9322SMatthew G. Knepley 565636b9322SMatthew G. Knepley Here we must get a subspace for the input DS 566636b9322SMatthew G. Knepley 5670de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest. 5680de51b56SMatthew G. Knepley 5690de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 5700de51b56SMatthew G. Knepley 5710de51b56SMatthew 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. 5720de51b56SMatthew 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 5735f790a90SMatthew 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 5740de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces. 5750de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 5760de51b56SMatthew G. Knepley 5770de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 5780de51b56SMatthew G. Knepley 5790de51b56SMatthew 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. 5800de51b56SMatthew G. Knepley */ 581d71ae5a4SJacob 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) 582d71ae5a4SJacob Faibussowitsch { 583a6e0b375SMatthew G. Knepley DM plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm; 584a6e0b375SMatthew G. Knepley DMEnclosureType encIn, encAux; 585a6e0b375SMatthew G. Knepley PetscDS ds = NULL, dsIn = NULL, dsAux = NULL; 586ca3d3a14SMatthew G. Knepley Vec localA = NULL, tv; 587aa0eca99SMatthew G. Knepley IS fieldIS; 58847923291SMatthew G. Knepley PetscSection section; 589636b9322SMatthew G. Knepley PetscDualSpace *sp, *cellsp, *spIn, *cellspIn; 590ef0bb6c7SMatthew G. Knepley PetscTabulation *T = NULL, *TAux = NULL; 5918c6c5593SMatthew G. Knepley PetscInt *Nc; 592636b9322SMatthew G. Knepley PetscInt dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f; 5938e6d238bSMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, isCohesiveIn = PETSC_FALSE, transform; 594c330f8ffSToby Isaac DMField coordField; 595c330f8ffSToby Isaac DMLabel depthLabel; 596c330f8ffSToby Isaac PetscQuadrature allPoints = NULL; 59747923291SMatthew G. Knepley 59847923291SMatthew G. Knepley PetscFunctionBegin; 5999566063dSJacob Faibussowitsch if (localU) PetscCall(VecGetDM(localU, &dmIn)); 600ad540459SPierre Jolivet else dmIn = dm; 6019566063dSJacob Faibussowitsch PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA)); 6029371c9d4SSatish Balay if (localA) PetscCall(VecGetDM(localA, &dmAux)); 603ad540459SPierre Jolivet else dmAux = NULL; 6049566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex)); 6059566063dSJacob Faibussowitsch PetscCall(DMConvert(dmIn, DMPLEX, &plexIn)); 6069566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn)); 6079566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux)); 6089566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 6099566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight)); 6109566063dSJacob Faibussowitsch PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm)); 6119566063dSJacob Faibussowitsch PetscCall(DMGetBasisTransformVec_Internal(dm, &tv)); 6129566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dm, &transform)); 6130de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */ 614a6e0b375SMatthew G. Knepley if (dmAux) { 6159566063dSJacob Faibussowitsch PetscCall(DMConvert(dmAux, DMPLEX, &plexAux)); 616ad540459SPierre 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"); 617636b9322SMatthew G. Knepley } 618e04ae0b4SMatthew G. Knepley if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL)); 619e04ae0b4SMatthew G. Knepley PetscCall(DMGetCoordinateField(dm, &coordField)); 620*e44f6aebSMatthew G. Knepley PetscCheck(coordField, PETSC_COMM_SELF, PETSC_ERR_USER, "DM must have a coordinate field"); 621e04ae0b4SMatthew G. Knepley /**** No collective calls below this point ****/ 622636b9322SMatthew G. Knepley /* Determine height for iteration of all meshes */ 623636b9322SMatthew G. Knepley { 624636b9322SMatthew G. Knepley DMPolytopeType ct, ctIn, ctAux; 62519552267SMatthew G. Knepley PetscInt minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux; 62619552267SMatthew G. Knepley PetscInt dim = -1, dimIn = -1, dimAux = -1; 62788aca1feSMatthew G. Knepley 6289566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd)); 629636b9322SMatthew G. Knepley if (pEnd > pStart) { 6309566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL)); 631636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart; 6329566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plex, p, &ct)); 633636b9322SMatthew G. Knepley dim = DMPolytopeTypeGetDim(ct); 6349566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn)); 6359566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL)); 6369566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn)); 637636b9322SMatthew G. Knepley dimIn = DMPolytopeTypeGetDim(ctIn); 638636b9322SMatthew G. Knepley if (dmAux) { 6399566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux)); 64019552267SMatthew G. Knepley PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux)); 64119552267SMatthew G. Knepley if (pStartAux < pEndAux) { 6429566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux)); 643636b9322SMatthew G. Knepley dimAux = DMPolytopeTypeGetDim(ctAux); 64419552267SMatthew G. Knepley } 645636b9322SMatthew G. Knepley } else dimAux = dim; 646e04ae0b4SMatthew G. Knepley } else { 647e04ae0b4SMatthew G. Knepley PetscCall(DMDestroy(&plex)); 648e04ae0b4SMatthew G. Knepley PetscCall(DMDestroy(&plexIn)); 649e04ae0b4SMatthew G. Knepley if (dmAux) PetscCall(DMDestroy(&plexAux)); 6503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 651e39fcb4eSStefano Zampini } 652636b9322SMatthew G. Knepley if (dim < 0) { 653636b9322SMatthew G. Knepley DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL; 654636b9322SMatthew G. Knepley 655636b9322SMatthew G. Knepley /* Fall back to determination based on being a submesh */ 6569566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(plex, &spmap)); 6579566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn)); 6589566063dSJacob Faibussowitsch if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux)); 659636b9322SMatthew G. Knepley dim = spmap ? 1 : 0; 660636b9322SMatthew G. Knepley dimIn = spmapIn ? 1 : 0; 661636b9322SMatthew G. Knepley dimAux = spmapAux ? 1 : 0; 66288aca1feSMatthew G. Knepley } 663636b9322SMatthew G. Knepley { 66419552267SMatthew G. Knepley PetscInt dimProj = PetscMin(PetscMin(dim, dimIn), (dimAux < 0 ? PETSC_MAX_INT : dimAux)); 66519552267SMatthew G. Knepley PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux; 666636b9322SMatthew G. Knepley 66719552267SMatthew 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"); 668636b9322SMatthew G. Knepley if (dimProj < dim) minHeight = 1; 669636b9322SMatthew G. Knepley htInc = dim - dimProj; 670636b9322SMatthew G. Knepley htIncIn = dimIn - dimProj; 67119552267SMatthew G. Knepley htIncAux = dimAuxEff - dimProj; 6720de51b56SMatthew G. Knepley } 673a6e0b375SMatthew G. Knepley } 6749566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(plex, &depth)); 6759566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(plex, &depthLabel)); 6769566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight)); 6772af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight); 67863a3b9bcSJacob 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); 6799566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, 0, NULL, &ds)); 6809566063dSJacob Faibussowitsch if (!ds) PetscCall(DMGetDS(dm, &ds)); 6819566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn)); 6829566063dSJacob Faibussowitsch if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn)); 6839566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 6849566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsIn, &NfIn)); 6858e6d238bSMatthew G. Knepley PetscCall(PetscDSIsCohesive(dsIn, &isCohesiveIn)); 6868e6d238bSMatthew G. Knepley if (isCohesiveIn) --htIncIn; // Should be rearranged 6879566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &NfTot)); 6889566063dSJacob Faibussowitsch PetscCall(DMFindRegionNum(dm, ds, ®ionNum)); 68907218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL, NULL)); 6909566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 6919566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 6929566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 6930c898c18SMatthew G. Knepley if (dmAux) { 6949566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmAux, &dsAux)); 6959566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 6960c898c18SMatthew G. Knepley } 6979566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc)); 6989566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn)); 6999566063dSJacob Faibussowitsch if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn)); 7009371c9d4SSatish Balay else { 7019371c9d4SSatish Balay cellsp = sp; 7029371c9d4SSatish Balay cellspIn = spIn; 7039371c9d4SSatish Balay } 7048c6c5593SMatthew G. Knepley /* Get cell dual spaces */ 70547923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 706665f567fSMatthew G. Knepley PetscDiscType disctype; 70747923291SMatthew G. Knepley 7089566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype)); 709665f567fSMatthew G. Knepley if (disctype == PETSC_DISC_FE) { 710665f567fSMatthew G. Knepley PetscFE fe; 71147923291SMatthew G. Knepley 71247923291SMatthew G. Knepley isFE[f] = PETSC_TRUE; 713665f567fSMatthew G. Knepley hasFE = PETSC_TRUE; 7149566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 7159566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &cellsp[f])); 716665f567fSMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) { 717665f567fSMatthew G. Knepley PetscFV fv; 71847923291SMatthew G. Knepley 71947923291SMatthew G. Knepley isFE[f] = PETSC_FALSE; 720665f567fSMatthew G. Knepley hasFV = PETSC_TRUE; 7219566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv)); 7229566063dSJacob Faibussowitsch PetscCall(PetscFVGetDualSpace(fv, &cellsp[f])); 723665f567fSMatthew G. Knepley } else { 724665f567fSMatthew G. Knepley isFE[f] = PETSC_FALSE; 725665f567fSMatthew G. Knepley cellsp[f] = NULL; 726665f567fSMatthew G. Knepley } 72747923291SMatthew G. Knepley } 728636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 729636b9322SMatthew G. Knepley PetscDiscType disctype; 730636b9322SMatthew G. Knepley 7319566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype)); 732636b9322SMatthew G. Knepley if (disctype == PETSC_DISC_FE) { 733636b9322SMatthew G. Knepley PetscFE fe; 734636b9322SMatthew G. Knepley 7359566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe)); 7369566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f])); 737636b9322SMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) { 738636b9322SMatthew G. Knepley PetscFV fv; 739636b9322SMatthew G. Knepley 7409566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv)); 7419566063dSJacob Faibussowitsch PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f])); 742636b9322SMatthew G. Knepley } else { 743636b9322SMatthew G. Knepley cellspIn[f] = NULL; 744636b9322SMatthew G. Knepley } 745636b9322SMatthew G. Knepley } 746636b9322SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 7479371c9d4SSatish Balay if (!htInc) { 7489371c9d4SSatish Balay sp[f] = cellsp[f]; 7499371c9d4SSatish Balay } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f])); 750636b9322SMatthew G. Knepley } 751ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 75274fc349aSMatthew G. Knepley PetscFE fem, subfem; 753665f567fSMatthew G. Knepley PetscDiscType disctype; 7544a3e9fdbSToby Isaac const PetscReal *points; 7554a3e9fdbSToby Isaac PetscInt numPoints; 7560c898c18SMatthew G. Knepley 75708401ef6SPierre Jolivet PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 7589566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints)); 7599566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL)); 7609566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux)); 761a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 7629371c9d4SSatish Balay if (!htIncIn) { 7639371c9d4SSatish Balay spIn[f] = cellspIn[f]; 7649371c9d4SSatish Balay } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f])); 765636b9322SMatthew G. Knepley 7669566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype)); 767665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 7689566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem)); 7699371c9d4SSatish Balay if (!htIncIn) { 7709371c9d4SSatish Balay subfem = fem; 7719371c9d4SSatish Balay } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem)); 7729566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f])); 7730c898c18SMatthew G. Knepley } 7740c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 7759566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype)); 776665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 7779566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem)); 7789371c9d4SSatish Balay if (!htIncAux) { 7799371c9d4SSatish Balay subfem = fem; 7809371c9d4SSatish Balay } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem)); 7819566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f])); 7820c898c18SMatthew G. Knepley } 7830c898c18SMatthew G. Knepley } 78447923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 7852af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) { 786636b9322SMatthew G. Knepley PetscInt hEff = h - minHeight + htInc; 787636b9322SMatthew G. Knepley PetscInt hEffIn = h - minHeight + htIncIn; 788636b9322SMatthew G. Knepley PetscInt hEffAux = h - minHeight + htIncAux; 789a6e0b375SMatthew G. Knepley PetscDS dsEff = ds; 790636b9322SMatthew G. Knepley PetscDS dsEffIn = dsIn; 791636b9322SMatthew G. Knepley PetscDS dsEffAux = dsAux; 7928c6c5593SMatthew G. Knepley PetscScalar *values; 793b7260050SToby Isaac PetscBool *fieldActive; 794b7260050SToby Isaac PetscInt maxDegree; 795e5e52638SMatthew G. Knepley PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 796c330f8ffSToby Isaac IS heightIS; 7978c6c5593SMatthew G. Knepley 798636b9322SMatthew G. Knepley if (h > minHeight) { 7999566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f])); 800636b9322SMatthew G. Knepley } 8019566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd)); 8029566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL)); 8039566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS)); 804c330f8ffSToby Isaac if (pEnd <= pStart) { 8059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS)); 806c330f8ffSToby Isaac continue; 807c330f8ffSToby Isaac } 80847923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */ 80947923291SMatthew G. Knepley totDim = 0; 81047923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 8115fedec97SMatthew G. Knepley PetscBool cohesive; 8125fedec97SMatthew G. Knepley 813665f567fSMatthew G. Knepley if (!sp[f]) continue; 8149566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 8159566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 8169c3cf19fSMatthew G. Knepley totDim += spDim; 8175fedec97SMatthew G. Knepley if (isCohesive && !cohesive) totDim += spDim; 81847923291SMatthew G. Knepley } 819636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart; 8209566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL)); 82163a3b9bcSJacob 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); 822c330f8ffSToby Isaac if (!totDim) { 8239566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS)); 824c330f8ffSToby Isaac continue; 825c330f8ffSToby Isaac } 8269566063dSJacob Faibussowitsch if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff)); 827636b9322SMatthew G. Knepley /* Compute totDimIn, the number of dofs in the closure of a point at this height */ 828636b9322SMatthew G. Knepley if (localU) { 829636b9322SMatthew G. Knepley PetscInt totDimIn, pIn, numValuesIn; 830636b9322SMatthew G. Knepley 831636b9322SMatthew G. Knepley totDimIn = 0; 832636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 8335fedec97SMatthew G. Knepley PetscBool cohesive; 8345fedec97SMatthew G. Knepley 835636b9322SMatthew G. Knepley if (!spIn[f]) continue; 8369566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive)); 8379566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim)); 838636b9322SMatthew G. Knepley totDimIn += spDim; 8398e6d238bSMatthew G. Knepley if (isCohesiveIn && !cohesive) totDimIn += spDim; 840636b9322SMatthew G. Knepley } 8419566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn)); 8429566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL)); 8438e6d238bSMatthew G. Knepley // TODO We could check that pIn is a cohesive cell for this check 8448e6d238bSMatthew 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); 8459566063dSJacob Faibussowitsch if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn)); 846636b9322SMatthew G. Knepley } 8479566063dSJacob Faibussowitsch if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux)); 84847923291SMatthew G. Knepley /* Loop over points at this height */ 8499566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values)); 8509566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive)); 851aa0eca99SMatthew G. Knepley { 852aa0eca99SMatthew G. Knepley const PetscInt *fields; 853aa0eca99SMatthew G. Knepley 8549566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fieldIS, &fields)); 855ad540459SPierre Jolivet for (f = 0; f < NfTot; ++f) fieldActive[f] = PETSC_FALSE; 856ad540459SPierre Jolivet for (f = 0; f < Nf; ++f) fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 8579566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fieldIS, &fields)); 858aa0eca99SMatthew G. Knepley } 8598c6c5593SMatthew G. Knepley if (label) { 8608c6c5593SMatthew G. Knepley PetscInt i; 86147923291SMatthew G. Knepley 86247923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 863c330f8ffSToby Isaac IS pointIS, isectIS; 86447923291SMatthew G. Knepley const PetscInt *points; 8658c6c5593SMatthew G. Knepley PetscInt n; 866c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 867c330f8ffSToby Isaac PetscQuadrature quad = NULL; 86847923291SMatthew G. Knepley 8699566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS)); 87047923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 8719566063dSJacob Faibussowitsch PetscCall(ISIntersect(pointIS, heightIS, &isectIS)); 8729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 873c330f8ffSToby Isaac if (!isectIS) continue; 8749566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isectIS, &n)); 8759566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isectIS, &points)); 8769566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree)); 87748a46eb9SPierre Jolivet if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad)); 878c330f8ffSToby Isaac if (!quad) { 879c330f8ffSToby Isaac if (!h && allPoints) { 880c330f8ffSToby Isaac quad = allPoints; 881c330f8ffSToby Isaac allPoints = NULL; 882c330f8ffSToby Isaac } else { 8839566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad)); 884c330f8ffSToby Isaac } 885c330f8ffSToby Isaac } 8869566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom)); 88747923291SMatthew G. Knepley for (p = 0; p < n; ++p) { 88847923291SMatthew G. Knepley const PetscInt point = points[p]; 88947923291SMatthew G. Knepley 8909566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, numValues)); 8919566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom)); 8929566063dSJacob Faibussowitsch PetscCall(DMPlexSetActivePoint(dm, point)); 893d0609cedSBarry 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)); 8949566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values)); 8959566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode)); 89647923291SMatthew G. Knepley } 8979566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom)); 8989566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&fegeom)); 8999566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 9009566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isectIS, &points)); 9019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isectIS)); 90247923291SMatthew G. Knepley } 9038c6c5593SMatthew G. Knepley } else { 904c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 905c330f8ffSToby Isaac PetscQuadrature quad = NULL; 906c330f8ffSToby Isaac IS pointIS; 907c330f8ffSToby Isaac 9089566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS)); 9099566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree)); 91048a46eb9SPierre Jolivet if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad)); 911c330f8ffSToby Isaac if (!quad) { 912c330f8ffSToby Isaac if (!h && allPoints) { 913c330f8ffSToby Isaac quad = allPoints; 914c330f8ffSToby Isaac allPoints = NULL; 915c330f8ffSToby Isaac } else { 9169566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad)); 917c330f8ffSToby Isaac } 918c330f8ffSToby Isaac } 9199566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom)); 9208c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 9219566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, numValues)); 9229566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom)); 9239566063dSJacob Faibussowitsch PetscCall(DMPlexSetActivePoint(dm, p)); 924d0609cedSBarry 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)); 9259566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values)); 9269566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode)); 9278c6c5593SMatthew G. Knepley } 9289566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom)); 9299566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&fegeom)); 9309566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 9319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 9328c6c5593SMatthew G. Knepley } 9339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS)); 9349566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values)); 9359566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive)); 93647923291SMatthew G. Knepley } 9378c6c5593SMatthew G. Knepley /* Cleanup */ 938ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 9399566063dSJacob Faibussowitsch for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f])); 9409566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f])); 9419566063dSJacob Faibussowitsch PetscCall(PetscFree2(T, TAux)); 9420c898c18SMatthew G. Knepley } 9439566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&allPoints)); 9449566063dSJacob Faibussowitsch PetscCall(PetscFree3(isFE, sp, spIn)); 9459566063dSJacob Faibussowitsch if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn)); 9469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 9479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plexIn)); 9489566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMDestroy(&plexAux)); 9493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95047923291SMatthew G. Knepley } 9518c6c5593SMatthew G. Knepley 952d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 953d71ae5a4SJacob Faibussowitsch { 9548c6c5593SMatthew G. Knepley PetscFunctionBegin; 9559566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX)); 9563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9578c6c5593SMatthew G. Knepley } 9588c6c5593SMatthew G. Knepley 959d71ae5a4SJacob 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) 960d71ae5a4SJacob Faibussowitsch { 9618c6c5593SMatthew G. Knepley PetscFunctionBegin; 9629566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX)); 9633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96447923291SMatthew G. Knepley } 96547923291SMatthew G. Knepley 966d71ae5a4SJacob 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) 967d71ae5a4SJacob Faibussowitsch { 96847923291SMatthew G. Knepley PetscFunctionBegin; 9699566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 9703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97147923291SMatthew G. Knepley } 97247923291SMatthew G. Knepley 973d71ae5a4SJacob 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) 974d71ae5a4SJacob Faibussowitsch { 9758c6c5593SMatthew G. Knepley PetscFunctionBegin; 9769566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 9773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97847923291SMatthew G. Knepley } 979ece3a9fcSMatthew G. Knepley 980d71ae5a4SJacob 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) 981d71ae5a4SJacob Faibussowitsch { 982ece3a9fcSMatthew G. Knepley PetscFunctionBegin; 9839566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 9843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 985ece3a9fcSMatthew G. Knepley } 986