147923291SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 247923291SMatthew G. Knepley 38c6c5593SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> 48c6c5593SMatthew G. Knepley 51b32699bSMatthew G. Knepley /*@ 61b32699bSMatthew G. Knepley DMPlexGetActivePoint - Get the point on which projection is currently working 71b32699bSMatthew G. Knepley 81b32699bSMatthew G. Knepley Not collective 91b32699bSMatthew G. Knepley 104165533cSJose E. Roman Input Parameter: 111b32699bSMatthew G. Knepley . dm - the DM 121b32699bSMatthew G. Knepley 134165533cSJose E. Roman Output Parameter: 141b32699bSMatthew G. Knepley . point - The mesh point involved in the current projection 151b32699bSMatthew G. Knepley 161b32699bSMatthew G. Knepley Level: developer 171b32699bSMatthew G. Knepley 18db781477SPatrick Sanan .seealso: `DMPlexSetActivePoint()` 191b32699bSMatthew G. Knepley @*/ 209371c9d4SSatish Balay PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point) { 211b32699bSMatthew G. Knepley PetscFunctionBeginHot; 221b32699bSMatthew G. Knepley *point = ((DM_Plex *)dm->data)->activePoint; 231b32699bSMatthew G. Knepley PetscFunctionReturn(0); 241b32699bSMatthew G. Knepley } 251b32699bSMatthew G. Knepley 261b32699bSMatthew G. Knepley /*@ 271b32699bSMatthew G. Knepley DMPlexSetActivePoint - Set the point on which projection is currently working 281b32699bSMatthew G. Knepley 291b32699bSMatthew G. Knepley Not collective 301b32699bSMatthew G. Knepley 314165533cSJose E. Roman Input Parameters: 321b32699bSMatthew G. Knepley + dm - the DM 331b32699bSMatthew G. Knepley - point - The mesh point involved in the current projection 341b32699bSMatthew G. Knepley 351b32699bSMatthew G. Knepley Level: developer 361b32699bSMatthew G. Knepley 37db781477SPatrick Sanan .seealso: `DMPlexGetActivePoint()` 381b32699bSMatthew G. Knepley @*/ 399371c9d4SSatish Balay PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point) { 401b32699bSMatthew G. Knepley PetscFunctionBeginHot; 411b32699bSMatthew G. Knepley ((DM_Plex *)dm->data)->activePoint = point; 421b32699bSMatthew G. Knepley PetscFunctionReturn(0); 431b32699bSMatthew G. Knepley } 441b32699bSMatthew G. Knepley 45a6e0b375SMatthew G. Knepley /* 46a6e0b375SMatthew G. Knepley DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point 47a6e0b375SMatthew G. Knepley 48a6e0b375SMatthew G. Knepley Input Parameters: 49a6e0b375SMatthew G. Knepley + dm - The output DM 50a6e0b375SMatthew G. Knepley . ds - The output DS 51a6e0b375SMatthew G. Knepley . dmIn - The input DM 52a6e0b375SMatthew G. Knepley . dsIn - The input DS 53a6e0b375SMatthew G. Knepley . time - The time for this evaluation 54a6e0b375SMatthew G. Knepley . fegeom - The FE geometry for this point 55a6e0b375SMatthew G. Knepley . fvgeom - The FV geometry for this point 56a6e0b375SMatthew G. Knepley . isFE - Flag indicating whether each output field has an FE discretization 57a6e0b375SMatthew G. Knepley . sp - The output PetscDualSpace for each field 58a6e0b375SMatthew G. Knepley . funcs - The evaluation function for each field 59a6e0b375SMatthew G. Knepley - ctxs - The user context for each field 60a6e0b375SMatthew G. Knepley 61a6e0b375SMatthew G. Knepley Output Parameter: 62a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space 63a6e0b375SMatthew G. Knepley 64a6e0b375SMatthew G. Knepley Level: developer 65a6e0b375SMatthew G. Knepley 66db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()` 67a6e0b375SMatthew G. Knepley */ 689371c9d4SSatish Balay 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[]) { 69a6e0b375SMatthew G. Knepley PetscInt coordDim, Nf, *Nc, f, spDim, d, v, tp; 705fedec97SMatthew G. Knepley PetscBool isAffine, isCohesive, transform; 718c6c5593SMatthew G. Knepley 728c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 739566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dmIn, &coordDim)); 749566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dmIn, &transform)); 759566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 769566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc)); 779566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 788c6c5593SMatthew G. Knepley /* Get values for closure */ 79c330f8ffSToby Isaac isAffine = fegeom->isAffine; 80c330f8ffSToby Isaac for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 818c6c5593SMatthew G. Knepley void *const ctx = ctxs ? ctxs[f] : NULL; 825fedec97SMatthew G. Knepley PetscBool cohesive; 838c6c5593SMatthew G. Knepley 848c6c5593SMatthew G. Knepley if (!sp[f]) continue; 859566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 869566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 878c6c5593SMatthew G. Knepley if (funcs[f]) { 88c330f8ffSToby Isaac if (isFE[f]) { 89c330f8ffSToby Isaac PetscQuadrature allPoints; 90c330f8ffSToby Isaac PetscInt q, dim, numPoints; 91c330f8ffSToby Isaac const PetscReal *points; 92c330f8ffSToby Isaac PetscScalar *pointEval; 93c330f8ffSToby Isaac PetscReal *x; 94ca3d3a14SMatthew G. Knepley DM rdm; 95c330f8ffSToby Isaac 969566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &rdm)); 979566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 989566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 999566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 1009566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(rdm, coordDim, MPIU_REAL, &x)); 101737e23dcSMatthew G. Knepley PetscCall(PetscArrayzero(pointEval, numPoints * Nc[f])); 102c330f8ffSToby Isaac for (q = 0; q < numPoints; q++, tp++) { 103c330f8ffSToby Isaac const PetscReal *v0; 104c330f8ffSToby Isaac 105c330f8ffSToby Isaac if (isAffine) { 106665f567fSMatthew G. Knepley const PetscReal *refpoint = &points[q * dim]; 107665f567fSMatthew G. Knepley PetscReal injpoint[3] = {0., 0., 0.}; 108665f567fSMatthew G. Knepley 109665f567fSMatthew G. Knepley if (dim != fegeom->dim) { 1105fedec97SMatthew G. Knepley if (isCohesive) { 111665f567fSMatthew G. Knepley /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */ 112665f567fSMatthew G. Knepley for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d]; 113665f567fSMatthew G. Knepley refpoint = injpoint; 11463a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %" PetscInt_FMT " != %" PetscInt_FMT " dual basis spatial dimension", fegeom->dim, dim); 115665f567fSMatthew G. Knepley } 116665f567fSMatthew G. Knepley CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x); 117c330f8ffSToby Isaac v0 = x; 1188c6c5593SMatthew G. Knepley } else { 119c330f8ffSToby Isaac v0 = &fegeom->v[tp * coordDim]; 1208c6c5593SMatthew G. Knepley } 1219371c9d4SSatish Balay if (transform) { 1229371c9d4SSatish Balay PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx)); 1239371c9d4SSatish Balay v0 = x; 1249371c9d4SSatish Balay } 1259566063dSJacob Faibussowitsch PetscCall((*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f] * q], ctx)); 126c330f8ffSToby Isaac } 1274bee2e38SMatthew G. Knepley /* Transform point evaluations pointEval[q,c] */ 1289566063dSJacob Faibussowitsch PetscCall(PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval)); 1299566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 1309566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(rdm, coordDim, MPIU_REAL, &x)); 1319566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(rdm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 132c330f8ffSToby Isaac v += spDim; 1335fedec97SMatthew G. Knepley if (isCohesive && !cohesive) { 13427f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 13527f02ce8SMatthew G. Knepley } 136c330f8ffSToby Isaac } else { 137*48a46eb9SPierre Jolivet for (d = 0; d < spDim; ++d, ++v) PetscCall(PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v])); 138c330f8ffSToby Isaac } 139c330f8ffSToby Isaac } else { 14027f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 1415fedec97SMatthew G. Knepley if (isCohesive && !cohesive) { 14227f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 14327f02ce8SMatthew G. Knepley } 1448c6c5593SMatthew G. Knepley } 1459c3cf19fSMatthew G. Knepley } 1468c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 1478c6c5593SMatthew G. Knepley } 1488c6c5593SMatthew G. Knepley 149a6e0b375SMatthew G. Knepley /* 150a6e0b375SMatthew G. Knepley DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point 151a6e0b375SMatthew G. Knepley 152a6e0b375SMatthew G. Knepley Input Parameters: 153a6e0b375SMatthew G. Knepley + dm - The output DM 154a6e0b375SMatthew G. Knepley . ds - The output DS 155a6e0b375SMatthew G. Knepley . dmIn - The input DM 156a6e0b375SMatthew G. Knepley . dsIn - The input DS 157a6e0b375SMatthew G. Knepley . dmAux - The auxiliary DM, which is always for the input space 158a6e0b375SMatthew G. Knepley . dsAux - The auxiliary DS, which is always for the input space 159a6e0b375SMatthew G. Knepley . time - The time for this evaluation 160a6e0b375SMatthew G. Knepley . localU - The local solution 161a6e0b375SMatthew G. Knepley . localA - The local auziliary fields 162a6e0b375SMatthew G. Knepley . cgeom - The FE geometry for this point 163a6e0b375SMatthew G. Knepley . sp - The output PetscDualSpace for each field 164a6e0b375SMatthew G. Knepley . p - The point in the output DM 165a5b23f4aSJose E. Roman . T - Input basis and derivatives for each field tabulated on the quadrature points 166ef0bb6c7SMatthew G. Knepley . TAux - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points 167a6e0b375SMatthew G. Knepley . funcs - The evaluation function for each field 168a6e0b375SMatthew G. Knepley - ctxs - The user context for each field 169a6e0b375SMatthew G. Knepley 170a6e0b375SMatthew G. Knepley Output Parameter: 171a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space 172a6e0b375SMatthew G. Knepley 173a6e0b375SMatthew G. Knepley Note: Not supported for FV 174a6e0b375SMatthew G. Knepley 175a6e0b375SMatthew G. Knepley Level: developer 176a6e0b375SMatthew G. Knepley 177db781477SPatrick Sanan .seealso: `DMProjectPoint_Field_Private()` 178a6e0b375SMatthew G. Knepley */ 1799371c9d4SSatish Balay 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[]) { 1808c6c5593SMatthew G. Knepley PetscSection section, sectionAux = NULL; 1814bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 1828c6c5593SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 1838c6c5593SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 184191494d9SMatthew G. Knepley const PetscScalar *constants; 1858c6c5593SMatthew G. Knepley PetscReal *x; 186ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc; 1874bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 1884bee2e38SMatthew G. Knepley const PetscInt dE = cgeom->dimEmbed; 189ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0; 1905fedec97SMatthew G. Knepley PetscBool isAffine, isCohesive, transform; 1918c6c5593SMatthew G. Knepley 1928c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 1939566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1949566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc)); 1959566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 1969566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsIn, &NfIn)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsIn, &uOff)); 1989566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x)); 1999566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x)); 2009566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL)); 2019566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(dsIn, &numConstants, &constants)); 2029566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dmIn, &transform)); 2039566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmIn, §ion)); 2049566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp)); 2059566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients)); 2068c6c5593SMatthew G. Knepley if (dmAux) { 20744171101SMatthew G. Knepley PetscInt subp; 2081ba34526SMatthew G. Knepley 2099566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp)); 2109566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 2119566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmAux, §ionAux)); 2129566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 2139566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 2149566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x)); 2159566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux)); 2168c6c5593SMatthew G. Knepley } 2178c6c5593SMatthew G. Knepley /* Get values for closure */ 2184bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 21927f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 22027f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 2214bee2e38SMatthew G. Knepley if (isAffine) { 2224bee2e38SMatthew G. Knepley fegeom.v = x; 2234bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2244bee2e38SMatthew G. Knepley fegeom.J = cgeom->J; 2254bee2e38SMatthew G. Knepley fegeom.invJ = cgeom->invJ; 2264bee2e38SMatthew G. Knepley fegeom.detJ = cgeom->detJ; 2274bee2e38SMatthew G. Knepley } 2288c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 229c330f8ffSToby Isaac PetscQuadrature allPoints; 230c330f8ffSToby Isaac PetscInt q, dim, numPoints; 231c330f8ffSToby Isaac const PetscReal *points; 232c330f8ffSToby Isaac PetscScalar *pointEval; 2335fedec97SMatthew G. Knepley PetscBool cohesive; 234c330f8ffSToby Isaac DM dm; 235c330f8ffSToby Isaac 2368c6c5593SMatthew G. Knepley if (!sp[f]) continue; 2379566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 2389566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 239c330f8ffSToby Isaac if (!funcs[f]) { 240be1504a2SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 2415fedec97SMatthew G. Knepley if (isCohesive && !cohesive) { 24227f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 24327f02ce8SMatthew G. Knepley } 244c330f8ffSToby Isaac continue; 245c330f8ffSToby Isaac } 2469566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &dm)); 2479566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 2489566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 2499566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 2500c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 251c330f8ffSToby Isaac if (isAffine) { 252ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q * dim], x); 2531c531cf8SMatthew G. Knepley } else { 2544bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[tp * dE]; 2554bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[tp * dE * dE]; 2564bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[tp * dE * dE]; 2574bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[tp]; 2588c6c5593SMatthew G. Knepley } 2599566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t)); 2609566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t)); 2619566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx)); 262a6e0b375SMatthew 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]); 2631c531cf8SMatthew G. Knepley } 2649566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 2659566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 266c330f8ffSToby Isaac v += spDim; 26727f02ce8SMatthew G. Knepley /* TODO: For now, set both sides equal, but this should use info from other support cell */ 2685fedec97SMatthew G. Knepley if (isCohesive && !cohesive) { 26927f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 27027f02ce8SMatthew G. Knepley } 2718c6c5593SMatthew G. Knepley } 2729566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients)); 2739566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux)); 2748c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 2758c6c5593SMatthew G. Knepley } 2768c6c5593SMatthew G. Knepley 2779371c9d4SSatish Balay 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[]) { 278b18799e0SMatthew G. Knepley PetscSection section, sectionAux = NULL; 2794bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 280b18799e0SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 281b18799e0SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 282b18799e0SMatthew G. Knepley const PetscScalar *constants; 283b18799e0SMatthew G. Knepley PetscReal *x; 284ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc; 2859f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 2864bee2e38SMatthew G. Knepley const PetscInt dE = fgeom->dimEmbed; 287ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 288b18799e0SMatthew G. Knepley PetscBool isAffine; 289b18799e0SMatthew G. Knepley 290b18799e0SMatthew G. Knepley PetscFunctionBeginHot; 29108401ef6SPierre Jolivet PetscCheck(dm == dmIn, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM"); 2929566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2939566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc)); 2949566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 2959566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 2969566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x)); 2979566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL)); 2989566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 2999566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 3009566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients)); 301b18799e0SMatthew G. Knepley if (dmAux) { 302a6e0b375SMatthew G. Knepley PetscInt subp; 303b18799e0SMatthew G. Knepley 3049566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp)); 3059566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 3069566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmAux, §ionAux)); 3079566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 3089566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 3099566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x)); 3109566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux)); 311b18799e0SMatthew G. Knepley } 312b18799e0SMatthew G. Knepley /* Get values for closure */ 3134bee2e38SMatthew G. Knepley isAffine = fgeom->isAffine; 314ea78f98cSLisandro Dalcin fegeom.n = NULL; 315ea78f98cSLisandro Dalcin fegeom.J = NULL; 316ea78f98cSLisandro Dalcin fegeom.v = NULL; 317ea78f98cSLisandro Dalcin fegeom.xi = NULL; 318a6e0b375SMatthew G. Knepley cgeom.dim = fgeom->dim; 319a6e0b375SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3204bee2e38SMatthew G. Knepley if (isAffine) { 3214bee2e38SMatthew G. Knepley fegeom.v = x; 3224bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3234bee2e38SMatthew G. Knepley fegeom.J = fgeom->J; 3244bee2e38SMatthew G. Knepley fegeom.invJ = fgeom->invJ; 3254bee2e38SMatthew G. Knepley fegeom.detJ = fgeom->detJ; 3264bee2e38SMatthew G. Knepley fegeom.n = fgeom->n; 3279f209ee4SMatthew G. Knepley 3289f209ee4SMatthew G. Knepley cgeom.J = fgeom->suppJ[0]; 3299f209ee4SMatthew G. Knepley cgeom.invJ = fgeom->suppInvJ[0]; 3309f209ee4SMatthew G. Knepley cgeom.detJ = fgeom->suppDetJ[0]; 3314bee2e38SMatthew G. Knepley } 332b18799e0SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 333b18799e0SMatthew G. Knepley PetscQuadrature allPoints; 334b18799e0SMatthew G. Knepley PetscInt q, dim, numPoints; 335b18799e0SMatthew G. Knepley const PetscReal *points; 336b18799e0SMatthew G. Knepley PetscScalar *pointEval; 337b18799e0SMatthew G. Knepley DM dm; 338b18799e0SMatthew G. Knepley 339b18799e0SMatthew G. Knepley if (!sp[f]) continue; 3409566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 341b18799e0SMatthew G. Knepley if (!funcs[f]) { 342b18799e0SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 343b18799e0SMatthew G. Knepley continue; 344b18799e0SMatthew G. Knepley } 3459566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp[f], &dm)); 3469566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &allPoints, NULL)); 3479566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, &dim, NULL, &numPoints, &points, NULL)); 3489566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 349b18799e0SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 350b18799e0SMatthew G. Knepley if (isAffine) { 351ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q * dim], x); 352b18799e0SMatthew G. Knepley } else { 3533fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[tp * dE]; 3549f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[tp * dE * dE]; 3559f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[tp * dE * dE]; 3564bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[tp]; 3574bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[tp * dE]; 3589f209ee4SMatthew G. Knepley 3599f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][tp * dE * dE]; 3609f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][tp * dE * dE]; 3619f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][tp]; 362b18799e0SMatthew G. Knepley } 363a6e0b375SMatthew 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 */ 3649566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t)); 3659566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t)); 3664bee2e38SMatthew 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]); 367b18799e0SMatthew G. Knepley } 3689566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceApplyAll(sp[f], pointEval, &values[v])); 3699566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numPoints * Nc[f], MPIU_SCALAR, &pointEval)); 370b18799e0SMatthew G. Knepley v += spDim; 371b18799e0SMatthew G. Knepley } 3729566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients)); 3739566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux)); 374b18799e0SMatthew G. Knepley PetscFunctionReturn(0); 375b18799e0SMatthew G. Knepley } 376b18799e0SMatthew G. Knepley 3779371c9d4SSatish Balay 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[]) { 3788c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom; 3798c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed; 3808c6c5593SMatthew G. Knepley 3818c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 3829566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3839566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 3849566063dSJacob Faibussowitsch if (hasFV) PetscCall(DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL)); 3858c6c5593SMatthew G. Knepley switch (type) { 3868c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL: 3879371c9d4SSatish Balay case DM_BC_NATURAL: PetscCall(DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode(**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *))funcs, ctxs, values)); break; 3888c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 3898c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD: 3909371c9d4SSatish 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)); 391d0609cedSBarry Smith break; 392b18799e0SMatthew G. Knepley case DM_BC_ESSENTIAL_BD_FIELD: 3939371c9d4SSatish 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)); 394d0609cedSBarry Smith break; 39598921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int)type); 3968c6c5593SMatthew G. Knepley } 3978c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 3988c6c5593SMatthew G. Knepley } 3998c6c5593SMatthew G. Knepley 4009371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) { 401df5c1128SToby Isaac PetscReal *points; 402df5c1128SToby Isaac PetscInt f, numPoints; 403df5c1128SToby Isaac 404df5c1128SToby Isaac PetscFunctionBegin; 40519552267SMatthew G. Knepley if (!dim) { 40619552267SMatthew G. Knepley PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints)); 40719552267SMatthew G. Knepley PetscFunctionReturn(0); 40819552267SMatthew G. Knepley } 409df5c1128SToby Isaac numPoints = 0; 410df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 411df5c1128SToby Isaac if (funcs[f]) { 412df5c1128SToby Isaac PetscQuadrature fAllPoints; 413df5c1128SToby Isaac PetscInt fNumPoints; 414df5c1128SToby Isaac 4159566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL)); 4169566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL)); 417df5c1128SToby Isaac numPoints += fNumPoints; 418df5c1128SToby Isaac } 419df5c1128SToby Isaac } 4209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * numPoints, &points)); 421df5c1128SToby Isaac numPoints = 0; 422df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 423df5c1128SToby Isaac if (funcs[f]) { 424df5c1128SToby Isaac PetscQuadrature fAllPoints; 42554ee0cceSMatthew G. Knepley PetscInt qdim, fNumPoints, q; 426df5c1128SToby Isaac const PetscReal *fPoints; 427df5c1128SToby Isaac 4289566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp[f], &fAllPoints, NULL)); 4299566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL)); 43063a3b9bcSJacob 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); 431df5c1128SToby Isaac for (q = 0; q < fNumPoints * dim; ++q) points[numPoints * dim + q] = fPoints[q]; 432df5c1128SToby Isaac numPoints += fNumPoints; 433df5c1128SToby Isaac } 434df5c1128SToby Isaac } 4359566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, allPoints)); 4369566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(*allPoints, dim, 0, numPoints, points, NULL)); 437df5c1128SToby Isaac PetscFunctionReturn(0); 438df5c1128SToby Isaac } 439df5c1128SToby Isaac 4405f15299fSJeremy L Thompson /*@C 4415f15299fSJeremy L Thompson DMGetFirstLabeledPoint - Find first labeled point p_o in odm such that the corresponding point p in dm has the specified height. Return p and the corresponding ds. 4425f15299fSJeremy L Thompson 4435f15299fSJeremy L Thompson Input Parameters: 4445f15299fSJeremy L Thompson dm - the DM 4455f15299fSJeremy L Thompson odm - the enclosing DM 4465f15299fSJeremy L Thompson label - label for DM domain, or NULL for whole domain 4475f15299fSJeremy L Thompson numIds - the number of ids 4485f15299fSJeremy L Thompson ids - An array of the label ids in sequence for the domain 4495f15299fSJeremy L Thompson height - Height of target cells in DMPlex topology 4505f15299fSJeremy L Thompson 4515f15299fSJeremy L Thompson Output Parameters: 4525f15299fSJeremy L Thompson point - the first labeled point 4535f15299fSJeremy L Thompson ds - the ds corresponding to the first labeled point 4545f15299fSJeremy L Thompson 4555f15299fSJeremy L Thompson Level: developer 456817da375SSatish Balay @*/ 4579371c9d4SSatish Balay PetscErrorCode DMGetFirstLabeledPoint(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *point, PetscDS *ds) { 458a6e0b375SMatthew G. Knepley DM plex; 459a6e0b375SMatthew G. Knepley DMEnclosureType enc; 460e9f0ba4eSJed Brown PetscInt ls = -1; 461e5e52638SMatthew G. Knepley 462e5e52638SMatthew G. Knepley PetscFunctionBegin; 4635f15299fSJeremy L Thompson if (point) *point = -1; 464e5e52638SMatthew G. Knepley if (!label) PetscFunctionReturn(0); 4659566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dm, odm, &enc)); 4669566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex)); 467e9f0ba4eSJed Brown for (PetscInt i = 0; i < numIds; ++i) { 468e9f0ba4eSJed Brown IS labelIS; 469e9f0ba4eSJed Brown PetscInt num_points, pStart, pEnd; 4709566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, ids[i], &labelIS)); 471e9f0ba4eSJed Brown if (!labelIS) continue; /* No points with that id on this process */ 4729566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(plex, height, &pStart, &pEnd)); 4739566063dSJacob Faibussowitsch PetscCall(ISGetSize(labelIS, &num_points)); 474e9f0ba4eSJed Brown if (num_points) { 475e5e52638SMatthew G. Knepley const PetscInt *points; 4769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(labelIS, &points)); 477e9f0ba4eSJed Brown for (PetscInt i = 0; i < num_points; i++) { 478e9f0ba4eSJed Brown PetscInt point; 4799566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dm, odm, enc, points[i], &point)); 480e9f0ba4eSJed Brown if (pStart <= point && point < pEnd) { 481a6e0b375SMatthew G. Knepley ls = point; 4829566063dSJacob Faibussowitsch if (ds) PetscCall(DMGetCellDS(dm, ls, ds)); 483e5e52638SMatthew G. Knepley } 484e9f0ba4eSJed Brown } 4859566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(labelIS, &points)); 486e9f0ba4eSJed Brown } 4879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&labelIS)); 488e5e52638SMatthew G. Knepley if (ls >= 0) break; 489e5e52638SMatthew G. Knepley } 4905f15299fSJeremy L Thompson if (point) *point = ls; 4919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 492e5e52638SMatthew G. Knepley PetscFunctionReturn(0); 493e5e52638SMatthew G. Knepley } 494e5e52638SMatthew G. Knepley 4950de51b56SMatthew G. Knepley /* 4960de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 4970de51b56SMatthew G. Knepley 4980de51b56SMatthew G. Knepley There are several different scenarios: 4990de51b56SMatthew G. Knepley 5000de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data 5010de51b56SMatthew G. Knepley 5020de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells. 5030de51b56SMatthew G. Knepley 5040de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data 5050de51b56SMatthew G. Knepley 5060de51b56SMatthew 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. 5070de51b56SMatthew G. Knepley 5080de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data 5090de51b56SMatthew G. Knepley 5100de51b56SMatthew 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. 5110de51b56SMatthew G. Knepley 512636b9322SMatthew G. Knepley 4) Volumetric input mesh with boundary output mesh 513636b9322SMatthew G. Knepley 514636b9322SMatthew G. Knepley Here we must get a subspace for the input DS 515636b9322SMatthew G. Knepley 5160de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest. 5170de51b56SMatthew G. Knepley 5180de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 5190de51b56SMatthew G. Knepley 5200de51b56SMatthew 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. 5210de51b56SMatthew 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 5225f790a90SMatthew 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 5230de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces. 5240de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 5250de51b56SMatthew G. Knepley 5260de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 5270de51b56SMatthew G. Knepley 5280de51b56SMatthew 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. 5290de51b56SMatthew G. Knepley */ 5309371c9d4SSatish Balay 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) { 531a6e0b375SMatthew G. Knepley DM plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm; 532a6e0b375SMatthew G. Knepley DMEnclosureType encIn, encAux; 533a6e0b375SMatthew G. Knepley PetscDS ds = NULL, dsIn = NULL, dsAux = NULL; 534ca3d3a14SMatthew G. Knepley Vec localA = NULL, tv; 535aa0eca99SMatthew G. Knepley IS fieldIS; 53647923291SMatthew G. Knepley PetscSection section; 537636b9322SMatthew G. Knepley PetscDualSpace *sp, *cellsp, *spIn, *cellspIn; 538ef0bb6c7SMatthew G. Knepley PetscTabulation *T = NULL, *TAux = NULL; 5398c6c5593SMatthew G. Knepley PetscInt *Nc; 540636b9322SMatthew G. Knepley PetscInt dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f; 5415fedec97SMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isCohesive = PETSC_FALSE, transform; 542c330f8ffSToby Isaac DMField coordField; 543c330f8ffSToby Isaac DMLabel depthLabel; 544c330f8ffSToby Isaac PetscQuadrature allPoints = NULL; 54547923291SMatthew G. Knepley 54647923291SMatthew G. Knepley PetscFunctionBegin; 5479566063dSJacob Faibussowitsch if (localU) PetscCall(VecGetDM(localU, &dmIn)); 548a6e0b375SMatthew G. Knepley else { dmIn = dm; } 5499566063dSJacob Faibussowitsch PetscCall(DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, 0, &localA)); 5509371c9d4SSatish Balay if (localA) PetscCall(VecGetDM(localA, &dmAux)); 5519371c9d4SSatish Balay else { dmAux = NULL; } 5529566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex)); 5539566063dSJacob Faibussowitsch PetscCall(DMConvert(dmIn, DMPLEX, &plexIn)); 5549566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dmIn, dm, &encIn)); 5559566063dSJacob Faibussowitsch PetscCall(DMGetEnclosureRelation(dmAux, dm, &encAux)); 5569566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5579566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plex, &minHeight)); 5589566063dSJacob Faibussowitsch PetscCall(DMGetBasisTransformDM_Internal(dm, &tdm)); 5599566063dSJacob Faibussowitsch PetscCall(DMGetBasisTransformVec_Internal(dm, &tv)); 5609566063dSJacob Faibussowitsch PetscCall(DMHasBasisTransform(dm, &transform)); 5610de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */ 562a6e0b375SMatthew G. Knepley if (dmAux) { 5639566063dSJacob Faibussowitsch PetscCall(DMConvert(dmAux, DMPLEX, &plexAux)); 5649371c9d4SSatish Balay 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"); } 565636b9322SMatthew G. Knepley } 566e04ae0b4SMatthew G. Knepley if (localU && localU != localX) PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL)); 567e04ae0b4SMatthew G. Knepley PetscCall(DMGetCoordinateField(dm, &coordField)); 568e04ae0b4SMatthew G. Knepley /**** No collective calls below this point ****/ 569636b9322SMatthew G. Knepley /* Determine height for iteration of all meshes */ 570636b9322SMatthew G. Knepley { 571636b9322SMatthew G. Knepley DMPolytopeType ct, ctIn, ctAux; 57219552267SMatthew G. Knepley PetscInt minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux, pEndAux; 57319552267SMatthew G. Knepley PetscInt dim = -1, dimIn = -1, dimAux = -1; 57488aca1feSMatthew G. Knepley 5759566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd)); 576636b9322SMatthew G. Knepley if (pEnd > pStart) { 5779566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, minHeight, &lStart, NULL)); 578636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart; 5799566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plex, p, &ct)); 580636b9322SMatthew G. Knepley dim = DMPolytopeTypeGetDim(ct); 5819566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plexIn, &minHeightIn)); 5829566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL)); 5839566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plexIn, pStartIn, &ctIn)); 584636b9322SMatthew G. Knepley dimIn = DMPolytopeTypeGetDim(ctIn); 585636b9322SMatthew G. Knepley if (dmAux) { 5869566063dSJacob Faibussowitsch PetscCall(DMPlexGetVTKCellHeight(plexAux, &minHeightAux)); 58719552267SMatthew G. Knepley PetscCall(DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, &pEndAux)); 58819552267SMatthew G. Knepley if (pStartAux < pEndAux) { 5899566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(plexAux, pStartAux, &ctAux)); 590636b9322SMatthew G. Knepley dimAux = DMPolytopeTypeGetDim(ctAux); 59119552267SMatthew G. Knepley } 592636b9322SMatthew G. Knepley } else dimAux = dim; 593e04ae0b4SMatthew G. Knepley } else { 594e04ae0b4SMatthew G. Knepley PetscCall(DMDestroy(&plex)); 595e04ae0b4SMatthew G. Knepley PetscCall(DMDestroy(&plexIn)); 596e04ae0b4SMatthew G. Knepley if (dmAux) PetscCall(DMDestroy(&plexAux)); 597e04ae0b4SMatthew G. Knepley PetscFunctionReturn(0); 598e39fcb4eSStefano Zampini } 599636b9322SMatthew G. Knepley if (dim < 0) { 600636b9322SMatthew G. Knepley DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL; 601636b9322SMatthew G. Knepley 602636b9322SMatthew G. Knepley /* Fall back to determination based on being a submesh */ 6039566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(plex, &spmap)); 6049566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(plexIn, &spmapIn)); 6059566063dSJacob Faibussowitsch if (plexAux) PetscCall(DMPlexGetSubpointMap(plexAux, &spmapAux)); 606636b9322SMatthew G. Knepley dim = spmap ? 1 : 0; 607636b9322SMatthew G. Knepley dimIn = spmapIn ? 1 : 0; 608636b9322SMatthew G. Knepley dimAux = spmapAux ? 1 : 0; 60988aca1feSMatthew G. Knepley } 610636b9322SMatthew G. Knepley { 61119552267SMatthew G. Knepley PetscInt dimProj = PetscMin(PetscMin(dim, dimIn), (dimAux < 0 ? PETSC_MAX_INT : dimAux)); 61219552267SMatthew G. Knepley PetscInt dimAuxEff = dimAux < 0 ? dimProj : dimAux; 613636b9322SMatthew G. Knepley 61419552267SMatthew 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"); 615636b9322SMatthew G. Knepley if (dimProj < dim) minHeight = 1; 616636b9322SMatthew G. Knepley htInc = dim - dimProj; 617636b9322SMatthew G. Knepley htIncIn = dimIn - dimProj; 61819552267SMatthew G. Knepley htIncAux = dimAuxEff - dimProj; 6190de51b56SMatthew G. Knepley } 620a6e0b375SMatthew G. Knepley } 6219566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(plex, &depth)); 6229566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(plex, &depthLabel)); 6239566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxProjectionHeight(plex, &maxHeight)); 6242af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight); 62563a3b9bcSJacob 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); 6269566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, 0, NULL, &ds)); 6279566063dSJacob Faibussowitsch if (!ds) PetscCall(DMGetDS(dm, &ds)); 6289566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn)); 6299566063dSJacob Faibussowitsch if (!dsIn) PetscCall(DMGetDS(dmIn, &dsIn)); 6309566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 6319566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsIn, &NfIn)); 6329566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &NfTot)); 6339566063dSJacob Faibussowitsch PetscCall(DMFindRegionNum(dm, ds, ®ionNum)); 6349566063dSJacob Faibussowitsch PetscCall(DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL)); 6359566063dSJacob Faibussowitsch PetscCall(PetscDSIsCohesive(ds, &isCohesive)); 6369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 6379566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 6380c898c18SMatthew G. Knepley if (dmAux) { 6399566063dSJacob Faibussowitsch PetscCall(DMGetDS(dmAux, &dsAux)); 6409566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 6410c898c18SMatthew G. Knepley } 6429566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponents(ds, &Nc)); 6439566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn)); 6449566063dSJacob Faibussowitsch if (maxHeight > 0) PetscCall(PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn)); 6459371c9d4SSatish Balay else { 6469371c9d4SSatish Balay cellsp = sp; 6479371c9d4SSatish Balay cellspIn = spIn; 6489371c9d4SSatish Balay } 6498c6c5593SMatthew G. Knepley /* Get cell dual spaces */ 65047923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 651665f567fSMatthew G. Knepley PetscDiscType disctype; 65247923291SMatthew G. Knepley 6539566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(ds, f, &disctype)); 654665f567fSMatthew G. Knepley if (disctype == PETSC_DISC_FE) { 655665f567fSMatthew G. Knepley PetscFE fe; 65647923291SMatthew G. Knepley 65747923291SMatthew G. Knepley isFE[f] = PETSC_TRUE; 658665f567fSMatthew G. Knepley hasFE = PETSC_TRUE; 6599566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 6609566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &cellsp[f])); 661665f567fSMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) { 662665f567fSMatthew G. Knepley PetscFV fv; 66347923291SMatthew G. Knepley 66447923291SMatthew G. Knepley isFE[f] = PETSC_FALSE; 665665f567fSMatthew G. Knepley hasFV = PETSC_TRUE; 6669566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fv)); 6679566063dSJacob Faibussowitsch PetscCall(PetscFVGetDualSpace(fv, &cellsp[f])); 668665f567fSMatthew G. Knepley } else { 669665f567fSMatthew G. Knepley isFE[f] = PETSC_FALSE; 670665f567fSMatthew G. Knepley cellsp[f] = NULL; 671665f567fSMatthew G. Knepley } 67247923291SMatthew G. Knepley } 673636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 674636b9322SMatthew G. Knepley PetscDiscType disctype; 675636b9322SMatthew G. Knepley 6769566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype)); 677636b9322SMatthew G. Knepley if (disctype == PETSC_DISC_FE) { 678636b9322SMatthew G. Knepley PetscFE fe; 679636b9322SMatthew G. Knepley 6809566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fe)); 6819566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &cellspIn[f])); 682636b9322SMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) { 683636b9322SMatthew G. Knepley PetscFV fv; 684636b9322SMatthew G. Knepley 6859566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fv)); 6869566063dSJacob Faibussowitsch PetscCall(PetscFVGetDualSpace(fv, &cellspIn[f])); 687636b9322SMatthew G. Knepley } else { 688636b9322SMatthew G. Knepley cellspIn[f] = NULL; 689636b9322SMatthew G. Knepley } 690636b9322SMatthew G. Knepley } 691636b9322SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 6929371c9d4SSatish Balay if (!htInc) { 6939371c9d4SSatish Balay sp[f] = cellsp[f]; 6949371c9d4SSatish Balay } else PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f])); 695636b9322SMatthew G. Knepley } 696ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 69774fc349aSMatthew G. Knepley PetscFE fem, subfem; 698665f567fSMatthew G. Knepley PetscDiscType disctype; 6994a3e9fdbSToby Isaac const PetscReal *points; 7004a3e9fdbSToby Isaac PetscInt numPoints; 7010c898c18SMatthew G. Knepley 70208401ef6SPierre Jolivet PetscCheck(maxHeight <= minHeight, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 7039566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &allPoints)); 7049566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL)); 7059566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(NfIn, &T, NfAux, &TAux)); 706a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 7079371c9d4SSatish Balay if (!htIncIn) { 7089371c9d4SSatish Balay spIn[f] = cellspIn[f]; 7099371c9d4SSatish Balay } else PetscCall(PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f])); 710636b9322SMatthew G. Knepley 7119566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsIn, f, &disctype)); 712665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 7139566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsIn, f, (PetscObject *)&fem)); 7149371c9d4SSatish Balay if (!htIncIn) { 7159371c9d4SSatish Balay subfem = fem; 7169371c9d4SSatish Balay } else PetscCall(PetscFEGetHeightSubspace(fem, htIncIn, &subfem)); 7179566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f])); 7180c898c18SMatthew G. Knepley } 7190c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 7209566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscType_Internal(dsAux, f, &disctype)); 721665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 7229566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(dsAux, f, (PetscObject *)&fem)); 7239371c9d4SSatish Balay if (!htIncAux) { 7249371c9d4SSatish Balay subfem = fem; 7259371c9d4SSatish Balay } else PetscCall(PetscFEGetHeightSubspace(fem, htIncAux, &subfem)); 7269566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f])); 7270c898c18SMatthew G. Knepley } 7280c898c18SMatthew G. Knepley } 72947923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 7302af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) { 731636b9322SMatthew G. Knepley PetscInt hEff = h - minHeight + htInc; 732636b9322SMatthew G. Knepley PetscInt hEffIn = h - minHeight + htIncIn; 733636b9322SMatthew G. Knepley PetscInt hEffAux = h - minHeight + htIncAux; 734a6e0b375SMatthew G. Knepley PetscDS dsEff = ds; 735636b9322SMatthew G. Knepley PetscDS dsEffIn = dsIn; 736636b9322SMatthew G. Knepley PetscDS dsEffAux = dsAux; 7378c6c5593SMatthew G. Knepley PetscScalar *values; 738b7260050SToby Isaac PetscBool *fieldActive; 739b7260050SToby Isaac PetscInt maxDegree; 740e5e52638SMatthew G. Knepley PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 741c330f8ffSToby Isaac IS heightIS; 7428c6c5593SMatthew G. Knepley 743636b9322SMatthew G. Knepley if (h > minHeight) { 7449566063dSJacob Faibussowitsch for (f = 0; f < Nf; ++f) PetscCall(PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f])); 745636b9322SMatthew G. Knepley } 7469566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd)); 7479566063dSJacob Faibussowitsch PetscCall(DMGetFirstLabeledPoint(dm, dm, label, numIds, ids, h, &lStart, NULL)); 7489566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(depthLabel, depth - h, &heightIS)); 749c330f8ffSToby Isaac if (pEnd <= pStart) { 7509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS)); 751c330f8ffSToby Isaac continue; 752c330f8ffSToby Isaac } 75347923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */ 75447923291SMatthew G. Knepley totDim = 0; 75547923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 7565fedec97SMatthew G. Knepley PetscBool cohesive; 7575fedec97SMatthew G. Knepley 758665f567fSMatthew G. Knepley if (!sp[f]) continue; 7599566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, f, &cohesive)); 7609566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(sp[f], &spDim)); 7619c3cf19fSMatthew G. Knepley totDim += spDim; 7625fedec97SMatthew G. Knepley if (isCohesive && !cohesive) totDim += spDim; 76347923291SMatthew G. Knepley } 764636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart; 7659566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL)); 76663a3b9bcSJacob 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); 767c330f8ffSToby Isaac if (!totDim) { 7689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS)); 769c330f8ffSToby Isaac continue; 770c330f8ffSToby Isaac } 7719566063dSJacob Faibussowitsch if (htInc) PetscCall(PetscDSGetHeightSubspace(ds, hEff, &dsEff)); 772636b9322SMatthew G. Knepley /* Compute totDimIn, the number of dofs in the closure of a point at this height */ 773636b9322SMatthew G. Knepley if (localU) { 774636b9322SMatthew G. Knepley PetscInt totDimIn, pIn, numValuesIn; 775636b9322SMatthew G. Knepley 776636b9322SMatthew G. Knepley totDimIn = 0; 777636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 7785fedec97SMatthew G. Knepley PetscBool cohesive; 7795fedec97SMatthew G. Knepley 780636b9322SMatthew G. Knepley if (!spIn[f]) continue; 7819566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(dsIn, f, &cohesive)); 7829566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(spIn[f], &spDim)); 783636b9322SMatthew G. Knepley totDimIn += spDim; 7845fedec97SMatthew G. Knepley if (isCohesive && !cohesive) totDimIn += spDim; 785636b9322SMatthew G. Knepley } 7869566063dSJacob Faibussowitsch PetscCall(DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn)); 7879566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL)); 78863a3b9bcSJacob Faibussowitsch PetscCheck(numValuesIn == totDimIn, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The input section point (%" PetscInt_FMT ") closure size %" PetscInt_FMT " != dual space dimension %" PetscInt_FMT " at height %" PetscInt_FMT, pIn, numValuesIn, totDimIn, htIncIn); 7899566063dSJacob Faibussowitsch if (htIncIn) PetscCall(PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn)); 790636b9322SMatthew G. Knepley } 7919566063dSJacob Faibussowitsch if (htIncAux) PetscCall(PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux)); 79247923291SMatthew G. Knepley /* Loop over points at this height */ 7939566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values)); 7949566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive)); 795aa0eca99SMatthew G. Knepley { 796aa0eca99SMatthew G. Knepley const PetscInt *fields; 797aa0eca99SMatthew G. Knepley 7989566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fieldIS, &fields)); 799aa0eca99SMatthew G. Knepley for (f = 0; f < NfTot; ++f) { fieldActive[f] = PETSC_FALSE; } 800aa0eca99SMatthew G. Knepley for (f = 0; f < Nf; ++f) { fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; } 8019566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fieldIS, &fields)); 802aa0eca99SMatthew G. Knepley } 8038c6c5593SMatthew G. Knepley if (label) { 8048c6c5593SMatthew G. Knepley PetscInt i; 80547923291SMatthew G. Knepley 80647923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 807c330f8ffSToby Isaac IS pointIS, isectIS; 80847923291SMatthew G. Knepley const PetscInt *points; 8098c6c5593SMatthew G. Knepley PetscInt n; 810c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 811c330f8ffSToby Isaac PetscQuadrature quad = NULL; 81247923291SMatthew G. Knepley 8139566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, ids[i], &pointIS)); 81447923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 8159566063dSJacob Faibussowitsch PetscCall(ISIntersect(pointIS, heightIS, &isectIS)); 8169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 817c330f8ffSToby Isaac if (!isectIS) continue; 8189566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isectIS, &n)); 8199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isectIS, &points)); 8209566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, isectIS, NULL, &maxDegree)); 821*48a46eb9SPierre Jolivet if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, isectIS, &quad)); 822c330f8ffSToby Isaac if (!quad) { 823c330f8ffSToby Isaac if (!h && allPoints) { 824c330f8ffSToby Isaac quad = allPoints; 825c330f8ffSToby Isaac allPoints = NULL; 826c330f8ffSToby Isaac } else { 8279566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, isCohesive ? dim - htInc - 1 : dim - htInc, funcs, &quad)); 828c330f8ffSToby Isaac } 829c330f8ffSToby Isaac } 8309566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom)); 83147923291SMatthew G. Knepley for (p = 0; p < n; ++p) { 83247923291SMatthew G. Knepley const PetscInt point = points[p]; 83347923291SMatthew G. Knepley 8349566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, numValues)); 8359566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(fegeom, p, p + 1, &chunkgeom)); 8369566063dSJacob Faibussowitsch PetscCall(DMPlexSetActivePoint(dm, point)); 837d0609cedSBarry 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)); 8389566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values)); 8399566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode)); 84047923291SMatthew G. Knepley } 8419566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(fegeom, p, p + 1, &chunkgeom)); 8429566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&fegeom)); 8439566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 8449566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isectIS, &points)); 8459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isectIS)); 84647923291SMatthew G. Knepley } 8478c6c5593SMatthew G. Knepley } else { 848c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 849c330f8ffSToby Isaac PetscQuadrature quad = NULL; 850c330f8ffSToby Isaac IS pointIS; 851c330f8ffSToby Isaac 8529566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, pEnd - pStart, pStart, 1, &pointIS)); 8539566063dSJacob Faibussowitsch PetscCall(DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree)); 854*48a46eb9SPierre Jolivet if (maxDegree <= 1) PetscCall(DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad)); 855c330f8ffSToby Isaac if (!quad) { 856c330f8ffSToby Isaac if (!h && allPoints) { 857c330f8ffSToby Isaac quad = allPoints; 858c330f8ffSToby Isaac allPoints = NULL; 859c330f8ffSToby Isaac } else { 8609566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllPointsUnion(Nf, sp, dim - htInc, funcs, &quad)); 861c330f8ffSToby Isaac } 862c330f8ffSToby Isaac } 8639566063dSJacob Faibussowitsch PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom)); 8648c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 8659566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(values, numValues)); 8669566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetChunk(fegeom, p - pStart, p - pStart + 1, &chunkgeom)); 8679566063dSJacob Faibussowitsch PetscCall(DMPlexSetActivePoint(dm, p)); 868d0609cedSBarry 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)); 8699566063dSJacob Faibussowitsch if (transform) PetscCall(DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values)); 8709566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode)); 8718c6c5593SMatthew G. Knepley } 8729566063dSJacob Faibussowitsch PetscCall(PetscFEGeomRestoreChunk(fegeom, p - pStart, pStart - p + 1, &chunkgeom)); 8739566063dSJacob Faibussowitsch PetscCall(PetscFEGeomDestroy(&fegeom)); 8749566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 8759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 8768c6c5593SMatthew G. Knepley } 8779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&heightIS)); 8789566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values)); 8799566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive)); 88047923291SMatthew G. Knepley } 8818c6c5593SMatthew G. Knepley /* Cleanup */ 882ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 8839566063dSJacob Faibussowitsch for (f = 0; f < NfIn; ++f) PetscCall(PetscTabulationDestroy(&T[f])); 8849566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(PetscTabulationDestroy(&TAux[f])); 8859566063dSJacob Faibussowitsch PetscCall(PetscFree2(T, TAux)); 8860c898c18SMatthew G. Knepley } 8879566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&allPoints)); 8889566063dSJacob Faibussowitsch PetscCall(PetscFree3(isFE, sp, spIn)); 8899566063dSJacob Faibussowitsch if (maxHeight > 0) PetscCall(PetscFree2(cellsp, cellspIn)); 8909566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 8919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plexIn)); 8929566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMDestroy(&plexAux)); 8938c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 89447923291SMatthew G. Knepley } 8958c6c5593SMatthew G. Knepley 8969371c9d4SSatish Balay PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) { 8978c6c5593SMatthew G. Knepley PetscFunctionBegin; 8989566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX)); 8998c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 9008c6c5593SMatthew G. Knepley } 9018c6c5593SMatthew G. Knepley 9029371c9d4SSatish Balay 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) { 9038c6c5593SMatthew G. Knepley PetscFunctionBegin; 9049566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void))funcs, ctxs, mode, localX)); 90547923291SMatthew G. Knepley PetscFunctionReturn(0); 90647923291SMatthew G. Knepley } 90747923291SMatthew G. Knepley 9089371c9d4SSatish Balay 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) { 90947923291SMatthew G. Knepley PetscFunctionBegin; 9109566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 9118c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 91247923291SMatthew G. Knepley } 91347923291SMatthew G. Knepley 9149371c9d4SSatish Balay 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) { 9158c6c5593SMatthew G. Knepley PetscFunctionBegin; 9169566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 91747923291SMatthew G. Knepley PetscFunctionReturn(0); 91847923291SMatthew G. Knepley } 919ece3a9fcSMatthew G. Knepley 9209371c9d4SSatish Balay 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) { 921ece3a9fcSMatthew G. Knepley PetscFunctionBegin; 9229566063dSJacob Faibussowitsch PetscCall(DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void))funcs, NULL, mode, localX)); 923ece3a9fcSMatthew G. Knepley PetscFunctionReturn(0); 924ece3a9fcSMatthew G. Knepley } 925