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 10*4165533cSJose E. Roman Input Parameter: 111b32699bSMatthew G. Knepley . dm - the DM 121b32699bSMatthew G. Knepley 13*4165533cSJose 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 181b32699bSMatthew G. Knepley .seealso: DMPlexSetActivePoint() 191b32699bSMatthew G. Knepley @*/ 20bdb10af2SPierre Jolivet PetscErrorCode DMPlexGetActivePoint(DM dm, PetscInt *point) 21bdb10af2SPierre Jolivet { 221b32699bSMatthew G. Knepley PetscFunctionBeginHot; 231b32699bSMatthew G. Knepley *point = ((DM_Plex *) dm->data)->activePoint; 241b32699bSMatthew G. Knepley PetscFunctionReturn(0); 251b32699bSMatthew G. Knepley } 261b32699bSMatthew G. Knepley 271b32699bSMatthew G. Knepley /*@ 281b32699bSMatthew G. Knepley DMPlexSetActivePoint - Set the point on which projection is currently working 291b32699bSMatthew G. Knepley 301b32699bSMatthew G. Knepley Not collective 311b32699bSMatthew G. Knepley 32*4165533cSJose E. Roman Input Parameters: 331b32699bSMatthew G. Knepley + dm - the DM 341b32699bSMatthew G. Knepley - point - The mesh point involved in the current projection 351b32699bSMatthew G. Knepley 361b32699bSMatthew G. Knepley Level: developer 371b32699bSMatthew G. Knepley 381b32699bSMatthew G. Knepley .seealso: DMPlexGetActivePoint() 391b32699bSMatthew G. Knepley @*/ 40bdb10af2SPierre Jolivet PetscErrorCode DMPlexSetActivePoint(DM dm, PetscInt point) 41bdb10af2SPierre Jolivet { 421b32699bSMatthew G. Knepley PetscFunctionBeginHot; 431b32699bSMatthew G. Knepley ((DM_Plex *) dm->data)->activePoint = point; 441b32699bSMatthew G. Knepley PetscFunctionReturn(0); 451b32699bSMatthew G. Knepley } 461b32699bSMatthew G. Knepley 47a6e0b375SMatthew G. Knepley /* 48a6e0b375SMatthew G. Knepley DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point 49a6e0b375SMatthew G. Knepley 50a6e0b375SMatthew G. Knepley Input Parameters: 51a6e0b375SMatthew G. Knepley + dm - The output DM 52a6e0b375SMatthew G. Knepley . ds - The output DS 53a6e0b375SMatthew G. Knepley . dmIn - The input DM 54a6e0b375SMatthew G. Knepley . dsIn - The input DS 55a6e0b375SMatthew G. Knepley . time - The time for this evaluation 56a6e0b375SMatthew G. Knepley . fegeom - The FE geometry for this point 57a6e0b375SMatthew G. Knepley . fvgeom - The FV geometry for this point 58a6e0b375SMatthew G. Knepley . isFE - Flag indicating whether each output field has an FE discretization 59a6e0b375SMatthew G. Knepley . sp - The output PetscDualSpace for each field 60a6e0b375SMatthew G. Knepley . funcs - The evaluation function for each field 61a6e0b375SMatthew G. Knepley - ctxs - The user context for each field 62a6e0b375SMatthew G. Knepley 63a6e0b375SMatthew G. Knepley Output Parameter: 64a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space 65a6e0b375SMatthew G. Knepley 66a6e0b375SMatthew G. Knepley Level: developer 67a6e0b375SMatthew G. Knepley 68a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private() 69a6e0b375SMatthew G. Knepley */ 70a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Func_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, PetscReal time, PetscFEGeom *fegeom, PetscFVCellGeom *fvgeom, PetscBool isFE[], PetscDualSpace sp[], 718c6c5593SMatthew G. Knepley PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, 728c6c5593SMatthew G. Knepley PetscScalar values[]) 7347923291SMatthew G. Knepley { 74a6e0b375SMatthew G. Knepley PetscInt coordDim, Nf, *Nc, f, spDim, d, v, tp; 7527f02ce8SMatthew G. Knepley PetscBool isAffine, isHybrid, transform; 768c6c5593SMatthew G. Knepley PetscErrorCode ierr; 778c6c5593SMatthew G. Knepley 788c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 79a6e0b375SMatthew G. Knepley ierr = DMGetCoordinateDim(dmIn, &coordDim);CHKERRQ(ierr); 80a6e0b375SMatthew G. Knepley ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr); 81a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 82a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 8327f02ce8SMatthew G. Knepley ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr); 848c6c5593SMatthew G. Knepley /* Get values for closure */ 85c330f8ffSToby Isaac isAffine = fegeom->isAffine; 86c330f8ffSToby Isaac for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 878c6c5593SMatthew G. Knepley void * const ctx = ctxs ? ctxs[f] : NULL; 888c6c5593SMatthew G. Knepley 898c6c5593SMatthew G. Knepley if (!sp[f]) continue; 908c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 918c6c5593SMatthew G. Knepley if (funcs[f]) { 92c330f8ffSToby Isaac if (isFE[f]) { 93c330f8ffSToby Isaac PetscQuadrature allPoints; 94c330f8ffSToby Isaac PetscInt q, dim, numPoints; 95c330f8ffSToby Isaac const PetscReal *points; 96c330f8ffSToby Isaac PetscScalar *pointEval; 97c330f8ffSToby Isaac PetscReal *x; 98ca3d3a14SMatthew G. Knepley DM rdm; 99c330f8ffSToby Isaac 100ca3d3a14SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr); 101b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr); 102c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 103ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 104ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 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) { 113665f567fSMatthew G. Knepley if (isHybrid) { 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; 117665f567fSMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %D != %D 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 } 124a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;} 125c330f8ffSToby Isaac ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr); 126c330f8ffSToby Isaac } 1274bee2e38SMatthew G. Knepley /* Transform point evaluations pointEval[q,c] */ 1282edcad52SToby Isaac ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr); 129c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 130ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 131ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 132c330f8ffSToby Isaac v += spDim; 13327f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 13427f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 13527f02ce8SMatthew G. Knepley } 136c330f8ffSToby Isaac } else { 137c330f8ffSToby Isaac for (d = 0; d < spDim; ++d, ++v) { 138c330f8ffSToby Isaac ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 139c330f8ffSToby Isaac } 140c330f8ffSToby Isaac } 141c330f8ffSToby Isaac } else { 14227f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 14327f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 14427f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 14527f02ce8SMatthew G. Knepley } 1468c6c5593SMatthew G. Knepley } 1479c3cf19fSMatthew G. Knepley } 1488c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 1498c6c5593SMatthew G. Knepley } 1508c6c5593SMatthew G. Knepley 151a6e0b375SMatthew G. Knepley /* 152a6e0b375SMatthew G. Knepley DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point 153a6e0b375SMatthew G. Knepley 154a6e0b375SMatthew G. Knepley Input Parameters: 155a6e0b375SMatthew G. Knepley + dm - The output DM 156a6e0b375SMatthew G. Knepley . ds - The output DS 157a6e0b375SMatthew G. Knepley . dmIn - The input DM 158a6e0b375SMatthew G. Knepley . dsIn - The input DS 159a6e0b375SMatthew G. Knepley . dmAux - The auxiliary DM, which is always for the input space 160a6e0b375SMatthew G. Knepley . dsAux - The auxiliary DS, which is always for the input space 161a6e0b375SMatthew G. Knepley . time - The time for this evaluation 162a6e0b375SMatthew G. Knepley . localU - The local solution 163a6e0b375SMatthew G. Knepley . localA - The local auziliary fields 164a6e0b375SMatthew G. Knepley . cgeom - The FE geometry for this point 165a6e0b375SMatthew G. Knepley . sp - The output PetscDualSpace for each field 166a6e0b375SMatthew G. Knepley . p - The point in the output DM 167a5b23f4aSJose E. Roman . T - Input basis and derivatives for each field tabulated on the quadrature points 168ef0bb6c7SMatthew G. Knepley . TAux - Auxiliary basis and derivatives for each aux field tabulated on the quadrature points 169a6e0b375SMatthew G. Knepley . funcs - The evaluation function for each field 170a6e0b375SMatthew G. Knepley - ctxs - The user context for each field 171a6e0b375SMatthew G. Knepley 172a6e0b375SMatthew G. Knepley Output Parameter: 173a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space 174a6e0b375SMatthew G. Knepley 175a6e0b375SMatthew G. Knepley Note: Not supported for FV 176a6e0b375SMatthew G. Knepley 177a6e0b375SMatthew G. Knepley Level: developer 178a6e0b375SMatthew G. Knepley 179a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private() 180a6e0b375SMatthew G. Knepley */ 181a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Field_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *cgeom, PetscDualSpace sp[], PetscInt p, 182ef0bb6c7SMatthew G. Knepley PetscTabulation *T, PetscTabulation *TAux, 1838c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 1848c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1858c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 186191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 1878c6c5593SMatthew G. Knepley PetscScalar values[]) 1888c6c5593SMatthew G. Knepley { 1898c6c5593SMatthew G. Knepley PetscSection section, sectionAux = NULL; 1904bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 1918c6c5593SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 1928c6c5593SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 193191494d9SMatthew G. Knepley const PetscScalar *constants; 1948c6c5593SMatthew G. Knepley PetscReal *x; 195ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc; 1964bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 1974bee2e38SMatthew G. Knepley const PetscInt dE = cgeom->dimEmbed; 198ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0; 19927f02ce8SMatthew G. Knepley PetscBool isAffine, isHybrid, transform; 2008c6c5593SMatthew G. Knepley PetscErrorCode ierr; 2018c6c5593SMatthew G. Knepley 2028c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 203a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 204a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 20527f02ce8SMatthew G. Knepley ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr); 206a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr); 207a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsIn, &uOff);CHKERRQ(ierr); 208a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x);CHKERRQ(ierr); 209a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 210a6e0b375SMatthew G. Knepley ierr = PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 211a6e0b375SMatthew G. Knepley ierr = PetscDSGetConstants(dsIn, &numConstants, &constants);CHKERRQ(ierr); 212a6e0b375SMatthew G. Knepley ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr); 213a6e0b375SMatthew G. Knepley ierr = DMGetLocalSection(dmIn, §ion);CHKERRQ(ierr); 214a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp);CHKERRQ(ierr); 215a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr); 2168c6c5593SMatthew G. Knepley if (dmAux) { 21744171101SMatthew G. Knepley PetscInt subp; 2181ba34526SMatthew G. Knepley 219a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr); 220a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 22192fd8e1eSJed Brown ierr = DMGetLocalSection(dmAux, §ionAux);CHKERRQ(ierr); 222a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 223a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 224a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 2251ba34526SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 2268c6c5593SMatthew G. Knepley } 2278c6c5593SMatthew G. Knepley /* Get values for closure */ 2284bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 22927f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 23027f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 2314bee2e38SMatthew G. Knepley if (isAffine) { 2324bee2e38SMatthew G. Knepley fegeom.v = x; 2334bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2344bee2e38SMatthew G. Knepley fegeom.J = cgeom->J; 2354bee2e38SMatthew G. Knepley fegeom.invJ = cgeom->invJ; 2364bee2e38SMatthew G. Knepley fegeom.detJ = cgeom->detJ; 2374bee2e38SMatthew G. Knepley } 2388c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 239c330f8ffSToby Isaac PetscQuadrature allPoints; 240c330f8ffSToby Isaac PetscInt q, dim, numPoints; 241c330f8ffSToby Isaac const PetscReal *points; 242c330f8ffSToby Isaac PetscScalar *pointEval; 243c330f8ffSToby Isaac DM dm; 244c330f8ffSToby Isaac 2458c6c5593SMatthew G. Knepley if (!sp[f]) continue; 2468c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 247c330f8ffSToby Isaac if (!funcs[f]) { 248be1504a2SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 24927f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 25027f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 25127f02ce8SMatthew G. Knepley } 252c330f8ffSToby Isaac continue; 253c330f8ffSToby Isaac } 254c330f8ffSToby Isaac ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 255b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr); 256c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 257c330f8ffSToby Isaac ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 2580c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 259c330f8ffSToby Isaac if (isAffine) { 260ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x); 2611c531cf8SMatthew G. Knepley } else { 2624bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[tp*dE]; 2634bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[tp*dE*dE]; 2644bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[tp*dE*dE]; 2654bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[tp]; 2668c6c5593SMatthew G. Knepley } 267ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 268ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 269a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);} 270a6e0b375SMatthew 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]); 2711c531cf8SMatthew G. Knepley } 272c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 273c330f8ffSToby Isaac ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 274c330f8ffSToby Isaac v += spDim; 27527f02ce8SMatthew G. Knepley /* TODO: For now, set both sides equal, but this should use info from other support cell */ 27627f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 27727f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 27827f02ce8SMatthew G. Knepley } 2798c6c5593SMatthew G. Knepley } 280a6e0b375SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr); 2818c6c5593SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 2828c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 2838c6c5593SMatthew G. Knepley } 2848c6c5593SMatthew G. Knepley 285a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_BdField_Private(DM dm, PetscDS ds, DM dmIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscReal time, Vec localU, Vec localA, PetscFEGeom *fgeom, PetscDualSpace sp[], PetscInt p, 286ef0bb6c7SMatthew G. Knepley PetscTabulation *T, PetscTabulation *TAux, 287b18799e0SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 288b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 289b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 290b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 291b18799e0SMatthew G. Knepley PetscScalar values[]) 292b18799e0SMatthew G. Knepley { 293b18799e0SMatthew G. Knepley PetscSection section, sectionAux = NULL; 2944bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 295b18799e0SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 296b18799e0SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 297b18799e0SMatthew G. Knepley const PetscScalar *constants; 298b18799e0SMatthew G. Knepley PetscReal *x; 299ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc; 3009f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 3014bee2e38SMatthew G. Knepley const PetscInt dE = fgeom->dimEmbed; 302ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 303b18799e0SMatthew G. Knepley PetscBool isAffine; 304b18799e0SMatthew G. Knepley PetscErrorCode ierr; 305b18799e0SMatthew G. Knepley 306b18799e0SMatthew G. Knepley PetscFunctionBeginHot; 307a6e0b375SMatthew G. Knepley if (dm != dmIn) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM"); 308a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 309a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 310a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 311a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 312a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 313a6e0b375SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 314a6e0b375SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 31592fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 316a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 317b18799e0SMatthew G. Knepley if (dmAux) { 318a6e0b375SMatthew G. Knepley PetscInt subp; 319b18799e0SMatthew G. Knepley 320a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr); 321a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 32292fd8e1eSJed Brown ierr = DMGetLocalSection(dmAux, §ionAux);CHKERRQ(ierr); 323a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 324a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 325a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 326b18799e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 327b18799e0SMatthew G. Knepley } 328b18799e0SMatthew G. Knepley /* Get values for closure */ 3294bee2e38SMatthew G. Knepley isAffine = fgeom->isAffine; 330ea78f98cSLisandro Dalcin fegeom.n = NULL; 331ea78f98cSLisandro Dalcin fegeom.J = NULL; 332ea78f98cSLisandro Dalcin fegeom.v = NULL; 333ea78f98cSLisandro Dalcin fegeom.xi = NULL; 334a6e0b375SMatthew G. Knepley cgeom.dim = fgeom->dim; 335a6e0b375SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3364bee2e38SMatthew G. Knepley if (isAffine) { 3374bee2e38SMatthew G. Knepley fegeom.v = x; 3384bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3394bee2e38SMatthew G. Knepley fegeom.J = fgeom->J; 3404bee2e38SMatthew G. Knepley fegeom.invJ = fgeom->invJ; 3414bee2e38SMatthew G. Knepley fegeom.detJ = fgeom->detJ; 3424bee2e38SMatthew G. Knepley fegeom.n = fgeom->n; 3439f209ee4SMatthew G. Knepley 3449f209ee4SMatthew G. Knepley cgeom.J = fgeom->suppJ[0]; 3459f209ee4SMatthew G. Knepley cgeom.invJ = fgeom->suppInvJ[0]; 3469f209ee4SMatthew G. Knepley cgeom.detJ = fgeom->suppDetJ[0]; 3474bee2e38SMatthew G. Knepley } 348b18799e0SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 349b18799e0SMatthew G. Knepley PetscQuadrature allPoints; 350b18799e0SMatthew G. Knepley PetscInt q, dim, numPoints; 351b18799e0SMatthew G. Knepley const PetscReal *points; 352b18799e0SMatthew G. Knepley PetscScalar *pointEval; 353b18799e0SMatthew G. Knepley DM dm; 354b18799e0SMatthew G. Knepley 355b18799e0SMatthew G. Knepley if (!sp[f]) continue; 356b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 357b18799e0SMatthew G. Knepley if (!funcs[f]) { 358b18799e0SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 359b18799e0SMatthew G. Knepley continue; 360b18799e0SMatthew G. Knepley } 361b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 362b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr); 363b18799e0SMatthew G. Knepley ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 364b18799e0SMatthew G. Knepley ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 365b18799e0SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 366b18799e0SMatthew G. Knepley if (isAffine) { 367ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x); 368b18799e0SMatthew G. Knepley } else { 3693fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[tp*dE]; 3709f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[tp*dE*dE]; 3719f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[tp*dE*dE]; 3724bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[tp]; 3734bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[tp*dE]; 3749f209ee4SMatthew G. Knepley 3759f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][tp*dE*dE]; 3769f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][tp*dE*dE]; 3779f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][tp]; 378b18799e0SMatthew G. Knepley } 379a6e0b375SMatthew 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 */ 380ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 381ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 3824bee2e38SMatthew 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]); 383b18799e0SMatthew G. Knepley } 384b18799e0SMatthew G. Knepley ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 385b18799e0SMatthew G. Knepley ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 386b18799e0SMatthew G. Knepley v += spDim; 387b18799e0SMatthew G. Knepley } 388a6e0b375SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 389b18799e0SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 390b18799e0SMatthew G. Knepley PetscFunctionReturn(0); 391b18799e0SMatthew G. Knepley } 392b18799e0SMatthew G. Knepley 393a6e0b375SMatthew G. Knepley static PetscErrorCode DMProjectPoint_Private(DM dm, PetscDS ds, DM dmIn, DMEnclosureType encIn, PetscDS dsIn, DM dmAux, DMEnclosureType encAux, PetscDS dsAux, PetscFEGeom *fegeom, PetscInt effectiveHeight, PetscReal time, Vec localU, Vec localA, PetscBool hasFE, PetscBool hasFV, PetscBool isFE[], 394ef0bb6c7SMatthew G. Knepley PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, 3958c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 3968c6c5593SMatthew G. Knepley { 3978c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom; 3988c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed; 3998c6c5593SMatthew G. Knepley PetscErrorCode ierr; 4008c6c5593SMatthew G. Knepley 4018c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 4028c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 4038c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 4048c6c5593SMatthew G. Knepley if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 4058c6c5593SMatthew G. Knepley switch (type) { 4068c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL: 4078c6c5593SMatthew G. Knepley case DM_BC_NATURAL: 408a6e0b375SMatthew G. Knepley ierr = DMProjectPoint_Func_Private(dm, ds, dmIn, dsIn, time, fegeom, &fvgeom, isFE, sp, (PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *)) funcs, ctxs, values);CHKERRQ(ierr);break; 4098c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 4108c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD: 411ef0bb6c7SMatthew G. Knepley ierr = DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, 4120c898c18SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 4130c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4140c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 415191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 416b18799e0SMatthew G. Knepley case DM_BC_ESSENTIAL_BD_FIELD: 417ef0bb6c7SMatthew G. Knepley ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, 418b18799e0SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 419b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 420b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 421b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 4228c6c5593SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 4238c6c5593SMatthew G. Knepley } 4248c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 4258c6c5593SMatthew G. Knepley } 4268c6c5593SMatthew G. Knepley 427df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 428df5c1128SToby Isaac { 429df5c1128SToby Isaac PetscReal *points; 430df5c1128SToby Isaac PetscInt f, numPoints; 431df5c1128SToby Isaac PetscErrorCode ierr; 432df5c1128SToby Isaac 433df5c1128SToby Isaac PetscFunctionBegin; 434df5c1128SToby Isaac numPoints = 0; 435df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 436df5c1128SToby Isaac if (funcs[f]) { 437df5c1128SToby Isaac PetscQuadrature fAllPoints; 438df5c1128SToby Isaac PetscInt fNumPoints; 439df5c1128SToby Isaac 440b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr); 441df5c1128SToby Isaac ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr); 442df5c1128SToby Isaac numPoints += fNumPoints; 443df5c1128SToby Isaac } 444df5c1128SToby Isaac } 445df5c1128SToby Isaac ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 446df5c1128SToby Isaac numPoints = 0; 447df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 448df5c1128SToby Isaac if (funcs[f]) { 449df5c1128SToby Isaac PetscQuadrature fAllPoints; 45054ee0cceSMatthew G. Knepley PetscInt qdim, fNumPoints, q; 451df5c1128SToby Isaac const PetscReal *fPoints; 452df5c1128SToby Isaac 453b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr); 45454ee0cceSMatthew G. Knepley ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr); 45554ee0cceSMatthew G. Knepley if (qdim != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Spatial dimension %D for dual basis does not match input dimension %D", qdim, dim); 456df5c1128SToby Isaac for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q]; 457df5c1128SToby Isaac numPoints += fNumPoints; 458df5c1128SToby Isaac } 459df5c1128SToby Isaac } 460df5c1128SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 461df5c1128SToby Isaac ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 462df5c1128SToby Isaac PetscFunctionReturn(0); 463df5c1128SToby Isaac } 464df5c1128SToby Isaac 465a6e0b375SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *ds) 466e5e52638SMatthew G. Knepley { 467a6e0b375SMatthew G. Knepley DM plex; 468a6e0b375SMatthew G. Knepley DMEnclosureType enc; 469e5e52638SMatthew G. Knepley DMLabel depthLabel; 470e5e52638SMatthew G. Knepley PetscInt dim, cdepth, ls = -1, i; 471e5e52638SMatthew G. Knepley PetscErrorCode ierr; 472e5e52638SMatthew G. Knepley 473e5e52638SMatthew G. Knepley PetscFunctionBegin; 474e5e52638SMatthew G. Knepley if (lStart) *lStart = -1; 475e5e52638SMatthew G. Knepley if (!label) PetscFunctionReturn(0); 476a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dm, odm, &enc);CHKERRQ(ierr); 477e5e52638SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 478a6e0b375SMatthew G. Knepley ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 479a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr); 480e5e52638SMatthew G. Knepley cdepth = dim - height; 481e5e52638SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 482e5e52638SMatthew G. Knepley IS pointIS; 483e5e52638SMatthew G. Knepley const PetscInt *points; 484a6e0b375SMatthew G. Knepley PetscInt pdepth, point; 485e5e52638SMatthew G. Knepley 486e5e52638SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 487e5e52638SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 488e5e52638SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 489a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dm, odm, enc, points[0], &point);CHKERRQ(ierr); 490a6e0b375SMatthew G. Knepley ierr = DMLabelGetValue(depthLabel, point, &pdepth);CHKERRQ(ierr); 491e5e52638SMatthew G. Knepley if (pdepth == cdepth) { 492a6e0b375SMatthew G. Knepley ls = point; 493a6e0b375SMatthew G. Knepley if (ds) {ierr = DMGetCellDS(dm, ls, ds);CHKERRQ(ierr);} 494e5e52638SMatthew G. Knepley } 495e5e52638SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 496e5e52638SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 497e5e52638SMatthew G. Knepley if (ls >= 0) break; 498e5e52638SMatthew G. Knepley } 499e5e52638SMatthew G. Knepley if (lStart) *lStart = ls; 500a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 501e5e52638SMatthew G. Knepley PetscFunctionReturn(0); 502e5e52638SMatthew G. Knepley } 503e5e52638SMatthew G. Knepley 5040de51b56SMatthew G. Knepley /* 5050de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 5060de51b56SMatthew G. Knepley 5070de51b56SMatthew G. Knepley There are several different scenarios: 5080de51b56SMatthew G. Knepley 5090de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data 5100de51b56SMatthew G. Knepley 5110de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells. 5120de51b56SMatthew G. Knepley 5130de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data 5140de51b56SMatthew G. Knepley 5150de51b56SMatthew 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. 5160de51b56SMatthew G. Knepley 5170de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data 5180de51b56SMatthew G. Knepley 5190de51b56SMatthew 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. 5200de51b56SMatthew G. Knepley 521636b9322SMatthew G. Knepley 4) Volumetric input mesh with boundary output mesh 522636b9322SMatthew G. Knepley 523636b9322SMatthew G. Knepley Here we must get a subspace for the input DS 524636b9322SMatthew G. Knepley 5250de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest. 5260de51b56SMatthew G. Knepley 5270de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 5280de51b56SMatthew G. Knepley 5290de51b56SMatthew 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. 5300de51b56SMatthew 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 5315f790a90SMatthew 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 5320de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces. 5330de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 5340de51b56SMatthew G. Knepley 5350de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 5360de51b56SMatthew G. Knepley 5370de51b56SMatthew 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. 5380de51b56SMatthew G. Knepley */ 53946fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 5401c531cf8SMatthew G. Knepley PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], 5418c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 5428c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 5438c6c5593SMatthew G. Knepley { 544a6e0b375SMatthew G. Knepley DM plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm; 545a6e0b375SMatthew G. Knepley DMEnclosureType encIn, encAux; 546a6e0b375SMatthew G. Knepley PetscDS ds = NULL, dsIn = NULL, dsAux = NULL; 547ca3d3a14SMatthew G. Knepley Vec localA = NULL, tv; 548aa0eca99SMatthew G. Knepley IS fieldIS; 54947923291SMatthew G. Knepley PetscSection section; 550636b9322SMatthew G. Knepley PetscDualSpace *sp, *cellsp, *spIn, *cellspIn; 551ef0bb6c7SMatthew G. Knepley PetscTabulation *T = NULL, *TAux = NULL; 5528c6c5593SMatthew G. Knepley PetscInt *Nc; 553636b9322SMatthew G. Knepley PetscInt dim, dimEmbed, depth, htInc = 0, htIncIn = 0, htIncAux = 0, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f; 554636b9322SMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, isHybrid = PETSC_FALSE, transform; 555c330f8ffSToby Isaac DMField coordField; 556c330f8ffSToby Isaac DMLabel depthLabel; 557c330f8ffSToby Isaac PetscQuadrature allPoints = NULL; 55847923291SMatthew G. Knepley PetscErrorCode ierr; 55947923291SMatthew G. Knepley 56047923291SMatthew G. Knepley PetscFunctionBegin; 561a6e0b375SMatthew G. Knepley if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);} 562a6e0b375SMatthew G. Knepley else {dmIn = dm;} 56304c51a94SMatthew G. Knepley ierr = DMGetAuxiliaryVec(dm, label, numIds ? ids[0] : 0, &localA);CHKERRQ(ierr); 5649a2a23afSMatthew G. Knepley if (localA) {ierr = VecGetDM(localA, &dmAux);CHKERRQ(ierr);} else {dmAux = NULL;} 565a6e0b375SMatthew G. Knepley ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 566a6e0b375SMatthew G. Knepley ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr); 567a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr); 568a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr); 5698c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 570a6e0b375SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr); 571ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr); 572ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr); 573ca3d3a14SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 5740de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */ 575a6e0b375SMatthew G. Knepley if (dmAux) { 576a6e0b375SMatthew G. Knepley ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr); 577a6e0b375SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 578a6e0b375SMatthew G. Knepley if (!localA) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Missing localA vector"); 579636b9322SMatthew G. Knepley } 580636b9322SMatthew G. Knepley } 581636b9322SMatthew G. Knepley /* Determine height for iteration of all meshes */ 582636b9322SMatthew G. Knepley { 583636b9322SMatthew G. Knepley DMPolytopeType ct, ctIn, ctAux; 584636b9322SMatthew G. Knepley PetscInt minHeightIn, minHeightAux, lStart, pStart, pEnd, p, pStartIn, pStartAux; 585636b9322SMatthew G. Knepley PetscInt dim = -1, dimIn, dimAux; 58688aca1feSMatthew G. Knepley 587636b9322SMatthew G. Knepley ierr = DMPlexGetSimplexOrBoxCells(plex, minHeight, &pStart, &pEnd);CHKERRQ(ierr); 588636b9322SMatthew G. Knepley if (pEnd > pStart) { 589636b9322SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, minHeight, &lStart, NULL);CHKERRQ(ierr); 590636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart; 591636b9322SMatthew G. Knepley ierr = DMPlexGetCellType(plex, p, &ct);CHKERRQ(ierr); 592636b9322SMatthew G. Knepley dim = DMPolytopeTypeGetDim(ct); 593636b9322SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(plexIn, &minHeightIn);CHKERRQ(ierr); 594636b9322SMatthew G. Knepley ierr = DMPlexGetSimplexOrBoxCells(plexIn, minHeightIn, &pStartIn, NULL);CHKERRQ(ierr); 595636b9322SMatthew G. Knepley ierr = DMPlexGetCellType(plexIn, pStartIn, &ctIn);CHKERRQ(ierr); 596636b9322SMatthew G. Knepley dimIn = DMPolytopeTypeGetDim(ctIn); 597636b9322SMatthew G. Knepley if (dmAux) { 598636b9322SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(plexAux, &minHeightAux);CHKERRQ(ierr); 599636b9322SMatthew G. Knepley ierr = DMPlexGetSimplexOrBoxCells(plexAux, minHeightAux, &pStartAux, NULL);CHKERRQ(ierr); 600636b9322SMatthew G. Knepley ierr = DMPlexGetCellType(plexAux, pStartAux, &ctAux);CHKERRQ(ierr); 601636b9322SMatthew G. Knepley dimAux = DMPolytopeTypeGetDim(ctAux); 602636b9322SMatthew G. Knepley } else dimAux = dim; 603e39fcb4eSStefano Zampini } 604636b9322SMatthew G. Knepley if (dim < 0) { 605636b9322SMatthew G. Knepley DMLabel spmap = NULL, spmapIn = NULL, spmapAux = NULL; 606636b9322SMatthew G. Knepley 607636b9322SMatthew G. Knepley /* Fall back to determination based on being a submesh */ 608636b9322SMatthew G. Knepley ierr = DMPlexGetSubpointMap(plex, &spmap);CHKERRQ(ierr); 609636b9322SMatthew G. Knepley ierr = DMPlexGetSubpointMap(plexIn, &spmapIn);CHKERRQ(ierr); 610636b9322SMatthew G. Knepley if (plexAux) {ierr = DMPlexGetSubpointMap(plexAux, &spmapAux);CHKERRQ(ierr);} 611636b9322SMatthew G. Knepley dim = spmap ? 1 : 0; 612636b9322SMatthew G. Knepley dimIn = spmapIn ? 1 : 0; 613636b9322SMatthew G. Knepley dimAux = spmapAux ? 1 : 0; 61488aca1feSMatthew G. Knepley } 615636b9322SMatthew G. Knepley { 616636b9322SMatthew G. Knepley PetscInt dimProj = PetscMin(PetscMin(dim, dimIn), dimAux); 617636b9322SMatthew G. Knepley 618636b9322SMatthew G. Knepley if (PetscAbsInt(dimProj - dim) > 1 || PetscAbsInt(dimProj - dimIn) > 1 || PetscAbsInt(dimProj - dimAux) > 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not currently support differences of more than 1 in dimension"); 619636b9322SMatthew G. Knepley if (dimProj < dim) minHeight = 1; 620636b9322SMatthew G. Knepley htInc = dim - dimProj; 621636b9322SMatthew G. Knepley htIncIn = dimIn - dimProj; 622636b9322SMatthew G. Knepley htIncAux = dimAux - dimProj; 6230de51b56SMatthew G. Knepley } 624a6e0b375SMatthew G. Knepley } 625a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 626a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr); 627a6e0b375SMatthew G. Knepley ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr); 6282af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight); 629e39fcb4eSStefano Zampini if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim); 630a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr); 631a6e0b375SMatthew G. Knepley if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);} 632a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr); 633a6e0b375SMatthew G. Knepley if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);} 634a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 635a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr); 636aa0eca99SMatthew G. Knepley ierr = DMGetNumFields(dm, &NfTot);CHKERRQ(ierr); 637aa0eca99SMatthew G. Knepley ierr = DMFindRegionNum(dm, ds, ®ionNum);CHKERRQ(ierr); 638aa0eca99SMatthew G. Knepley ierr = DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL);CHKERRQ(ierr); 63927f02ce8SMatthew G. Knepley ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr); 6408c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 64192fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 6420c898c18SMatthew G. Knepley if (dmAux) { 643a6e0b375SMatthew G. Knepley ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr); 644a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 6450c898c18SMatthew G. Knepley } 646a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 647636b9322SMatthew G. Knepley ierr = PetscMalloc3(Nf, &isFE, Nf, &sp, NfIn, &spIn);CHKERRQ(ierr); 648636b9322SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscMalloc2(Nf, &cellsp, NfIn, &cellspIn);CHKERRQ(ierr);} 649636b9322SMatthew G. Knepley else {cellsp = sp; cellspIn = spIn;} 650a6e0b375SMatthew G. Knepley if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 6518c6c5593SMatthew G. Knepley /* Get cell dual spaces */ 65247923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 653665f567fSMatthew G. Knepley PetscDiscType disctype; 65447923291SMatthew G. Knepley 655665f567fSMatthew G. Knepley ierr = PetscDSGetDiscType_Internal(ds, f, &disctype);CHKERRQ(ierr); 656665f567fSMatthew G. Knepley if (disctype == PETSC_DISC_FE) { 657665f567fSMatthew G. Knepley PetscFE fe; 65847923291SMatthew G. Knepley 65947923291SMatthew G. Knepley isFE[f] = PETSC_TRUE; 660665f567fSMatthew G. Knepley hasFE = PETSC_TRUE; 661665f567fSMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr); 66247923291SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 663665f567fSMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) { 664665f567fSMatthew G. Knepley PetscFV fv; 66547923291SMatthew G. Knepley 66647923291SMatthew G. Knepley isFE[f] = PETSC_FALSE; 667665f567fSMatthew G. Knepley hasFV = PETSC_TRUE; 668665f567fSMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fv);CHKERRQ(ierr); 66947923291SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 670665f567fSMatthew G. Knepley } else { 671665f567fSMatthew G. Knepley isFE[f] = PETSC_FALSE; 672665f567fSMatthew G. Knepley cellsp[f] = NULL; 673665f567fSMatthew G. Knepley } 67447923291SMatthew G. Knepley } 675636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 676636b9322SMatthew G. Knepley PetscDiscType disctype; 677636b9322SMatthew G. Knepley 678636b9322SMatthew G. Knepley ierr = PetscDSGetDiscType_Internal(dsIn, f, &disctype);CHKERRQ(ierr); 679636b9322SMatthew G. Knepley if (disctype == PETSC_DISC_FE) { 680636b9322SMatthew G. Knepley PetscFE fe; 681636b9322SMatthew G. Knepley 682636b9322SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fe);CHKERRQ(ierr); 683636b9322SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &cellspIn[f]);CHKERRQ(ierr); 684636b9322SMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) { 685636b9322SMatthew G. Knepley PetscFV fv; 686636b9322SMatthew G. Knepley 687636b9322SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fv);CHKERRQ(ierr); 688636b9322SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &cellspIn[f]);CHKERRQ(ierr); 689636b9322SMatthew G. Knepley } else { 690636b9322SMatthew G. Knepley cellspIn[f] = NULL; 691636b9322SMatthew G. Knepley } 692636b9322SMatthew G. Knepley } 693c330f8ffSToby Isaac ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr); 694636b9322SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 695636b9322SMatthew G. Knepley if (!htInc) {sp[f] = cellsp[f];} 696636b9322SMatthew G. Knepley else {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], htInc, &sp[f]);CHKERRQ(ierr);} 697636b9322SMatthew G. Knepley } 698ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 69974fc349aSMatthew G. Knepley PetscFE fem, subfem; 700665f567fSMatthew G. Knepley PetscDiscType disctype; 7014a3e9fdbSToby Isaac const PetscReal *points; 7024a3e9fdbSToby Isaac PetscInt numPoints; 7030c898c18SMatthew G. Knepley 7042af58022SMatthew G. Knepley if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 705636b9322SMatthew G. Knepley ierr = PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &allPoints);CHKERRQ(ierr); 706c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints, NULL, NULL, &numPoints, &points, NULL);CHKERRQ(ierr); 707ef0bb6c7SMatthew G. Knepley ierr = PetscMalloc2(NfIn, &T, NfAux, &TAux);CHKERRQ(ierr); 708a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 709636b9322SMatthew G. Knepley if (!htIncIn) {spIn[f] = cellspIn[f];} 710636b9322SMatthew G. Knepley else {ierr = PetscDualSpaceGetHeightSubspace(cellspIn[f], htIncIn, &spIn[f]);CHKERRQ(ierr);} 711636b9322SMatthew G. Knepley 712665f567fSMatthew G. Knepley ierr = PetscDSGetDiscType_Internal(dsIn, f, &disctype);CHKERRQ(ierr); 713665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 714a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr); 715636b9322SMatthew G. Knepley if (!htIncIn) {subfem = fem;} 716636b9322SMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, htIncIn, &subfem);CHKERRQ(ierr);} 717ef0bb6c7SMatthew G. Knepley ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]);CHKERRQ(ierr); 7180c898c18SMatthew G. Knepley } 7190c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 720665f567fSMatthew G. Knepley ierr = PetscDSGetDiscType_Internal(dsAux, f, &disctype);CHKERRQ(ierr); 721665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 722a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 723636b9322SMatthew G. Knepley if (!htIncAux) {subfem = fem;} 724636b9322SMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, htIncAux, &subfem);CHKERRQ(ierr);} 725ef0bb6c7SMatthew G. Knepley ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]);CHKERRQ(ierr); 7260c898c18SMatthew G. Knepley } 7270c898c18SMatthew G. Knepley } 72847923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 7292af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) { 730636b9322SMatthew G. Knepley PetscInt hEff = h - minHeight + htInc; 731636b9322SMatthew G. Knepley PetscInt hEffIn = h - minHeight + htIncIn; 732636b9322SMatthew G. Knepley PetscInt hEffAux = h - minHeight + htIncAux; 733a6e0b375SMatthew G. Knepley PetscDS dsEff = ds; 734636b9322SMatthew G. Knepley PetscDS dsEffIn = dsIn; 735636b9322SMatthew G. Knepley PetscDS dsEffAux = dsAux; 7368c6c5593SMatthew G. Knepley PetscScalar *values; 737b7260050SToby Isaac PetscBool *fieldActive; 738b7260050SToby Isaac PetscInt maxDegree; 739e5e52638SMatthew G. Knepley PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 740c330f8ffSToby Isaac IS heightIS; 7418c6c5593SMatthew G. Knepley 742636b9322SMatthew G. Knepley if (h > minHeight) { 743636b9322SMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], hEff, &sp[f]);CHKERRQ(ierr);} 744636b9322SMatthew G. Knepley } 745412e9a14SMatthew G. Knepley ierr = DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd);CHKERRQ(ierr); 746a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr); 747c330f8ffSToby Isaac ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr); 748c330f8ffSToby Isaac if (pEnd <= pStart) { 749c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 750c330f8ffSToby Isaac continue; 751c330f8ffSToby Isaac } 75247923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */ 75347923291SMatthew G. Knepley totDim = 0; 75447923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 755665f567fSMatthew G. Knepley if (!sp[f]) continue; 75647923291SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 7579c3cf19fSMatthew G. Knepley totDim += spDim; 75827f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) totDim += spDim; 75947923291SMatthew G. Knepley } 760636b9322SMatthew G. Knepley p = lStart < 0 ? pStart : lStart; 761636b9322SMatthew G. Knepley ierr = DMPlexVecGetClosure(plex, section, localX, p, &numValues, NULL);CHKERRQ(ierr); 762636b9322SMatthew G. Knepley if (numValues != totDim) SETERRQ6(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The output section point (%D) closure size %D != dual space dimension %D at height %D in [%D, %D]", p, numValues, totDim, h, minHeight, maxHeight); 763c330f8ffSToby Isaac if (!totDim) { 764c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 765c330f8ffSToby Isaac continue; 766c330f8ffSToby Isaac } 767636b9322SMatthew G. Knepley if (htInc) {ierr = PetscDSGetHeightSubspace(ds, hEff, &dsEff);CHKERRQ(ierr);} 768636b9322SMatthew G. Knepley /* Compute totDimIn, the number of dofs in the closure of a point at this height */ 769636b9322SMatthew G. Knepley if (localU) { 770636b9322SMatthew G. Knepley PetscInt totDimIn, pIn, numValuesIn; 771636b9322SMatthew G. Knepley 772636b9322SMatthew G. Knepley totDimIn = 0; 773636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 774636b9322SMatthew G. Knepley if (!spIn[f]) continue; 775636b9322SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(spIn[f], &spDim);CHKERRQ(ierr); 776636b9322SMatthew G. Knepley totDimIn += spDim; 777636b9322SMatthew G. Knepley if (isHybrid && (f < Nf-1)) totDimIn += spDim; 778636b9322SMatthew G. Knepley } 779636b9322SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmIn, dm, encIn, lStart < 0 ? pStart : lStart, &pIn);CHKERRQ(ierr); 780636b9322SMatthew G. Knepley ierr = DMPlexVecGetClosure(plexIn, NULL, localU, pIn, &numValuesIn, NULL);CHKERRQ(ierr); 781636b9322SMatthew G. Knepley if (numValuesIn != totDimIn) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The input section point (%D) closure size %D != dual space dimension %D at height %D", pIn, numValuesIn, totDimIn, htIncIn); 782636b9322SMatthew G. Knepley if (htIncIn) {ierr = PetscDSGetHeightSubspace(dsIn, hEffIn, &dsEffIn);CHKERRQ(ierr);} 783636b9322SMatthew G. Knepley } 784636b9322SMatthew G. Knepley if (htIncAux) {ierr = PetscDSGetHeightSubspace(dsAux, hEffAux, &dsEffAux);CHKERRQ(ierr);} 78547923291SMatthew G. Knepley /* Loop over points at this height */ 78669291d52SBarry Smith ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 787aa0eca99SMatthew G. Knepley ierr = DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive);CHKERRQ(ierr); 788aa0eca99SMatthew G. Knepley { 789aa0eca99SMatthew G. Knepley const PetscInt *fields; 790aa0eca99SMatthew G. Knepley 791aa0eca99SMatthew G. Knepley ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr); 792aa0eca99SMatthew G. Knepley for (f = 0; f < NfTot; ++f) {fieldActive[f] = PETSC_FALSE;} 793aa0eca99SMatthew G. Knepley for (f = 0; f < Nf; ++f) {fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;} 794aa0eca99SMatthew G. Knepley ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr); 795aa0eca99SMatthew G. Knepley } 7968c6c5593SMatthew G. Knepley if (label) { 7978c6c5593SMatthew G. Knepley PetscInt i; 79847923291SMatthew G. Knepley 79947923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 800c330f8ffSToby Isaac IS pointIS, isectIS; 80147923291SMatthew G. Knepley const PetscInt *points; 8028c6c5593SMatthew G. Knepley PetscInt n; 803c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 804c330f8ffSToby Isaac PetscQuadrature quad = NULL; 80547923291SMatthew G. Knepley 80647923291SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 80747923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 808c330f8ffSToby Isaac ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr); 809c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 810c330f8ffSToby Isaac if (!isectIS) continue; 811c330f8ffSToby Isaac ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr); 812c330f8ffSToby Isaac ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr); 813b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr); 814b7260050SToby Isaac if (maxDegree <= 1) { 815c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr); 816c330f8ffSToby Isaac } 817c330f8ffSToby Isaac if (!quad) { 818c330f8ffSToby Isaac if (!h && allPoints) { 819c330f8ffSToby Isaac quad = allPoints; 820c330f8ffSToby Isaac allPoints = NULL; 821c330f8ffSToby Isaac } else { 822636b9322SMatthew G. Knepley ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,isHybrid ? dim-htInc-1 : dim-htInc,funcs,&quad);CHKERRQ(ierr); 823c330f8ffSToby Isaac } 824c330f8ffSToby Isaac } 825636b9322SMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField, isectIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom);CHKERRQ(ierr); 82647923291SMatthew G. Knepley for (p = 0; p < n; ++p) { 82747923291SMatthew G. Knepley const PetscInt point = points[p]; 82847923291SMatthew G. Knepley 829580bdb30SBarry Smith ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr); 830c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 8311b32699bSMatthew G. Knepley ierr = DMPlexSetActivePoint(dm, point);CHKERRQ(ierr); 832636b9322SMatthew G. Knepley ierr = 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); 83347923291SMatthew G. Knepley if (ierr) { 83447923291SMatthew G. Knepley PetscErrorCode ierr2; 83569291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 83669291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 83747923291SMatthew G. Knepley CHKERRQ(ierr); 83847923291SMatthew G. Knepley } 839a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 8405f790a90SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode);CHKERRQ(ierr); 84147923291SMatthew G. Knepley } 842c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 843c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 844c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 845c330f8ffSToby Isaac ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr); 846c330f8ffSToby Isaac ierr = ISDestroy(&isectIS);CHKERRQ(ierr); 84747923291SMatthew G. Knepley } 8488c6c5593SMatthew G. Knepley } else { 849c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 850c330f8ffSToby Isaac PetscQuadrature quad = NULL; 851c330f8ffSToby Isaac IS pointIS; 852c330f8ffSToby Isaac 853c330f8ffSToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr); 854b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 855b7260050SToby Isaac if (maxDegree <= 1) { 856c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr); 857c330f8ffSToby Isaac } 858c330f8ffSToby Isaac if (!quad) { 859c330f8ffSToby Isaac if (!h && allPoints) { 860c330f8ffSToby Isaac quad = allPoints; 861c330f8ffSToby Isaac allPoints = NULL; 862c330f8ffSToby Isaac } else { 863636b9322SMatthew G. Knepley ierr = PetscDualSpaceGetAllPointsUnion(Nf, sp, dim-htInc, funcs, &quad);CHKERRQ(ierr); 864c330f8ffSToby Isaac } 865c330f8ffSToby Isaac } 866636b9322SMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField, pointIS, quad, (htInc && h == minHeight) ? PETSC_TRUE : PETSC_FALSE, &fegeom);CHKERRQ(ierr); 8678c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 868580bdb30SBarry Smith ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr); 869c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr); 8701b32699bSMatthew G. Knepley ierr = DMPlexSetActivePoint(dm, p);CHKERRQ(ierr); 871636b9322SMatthew G. Knepley ierr = 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); 8728c6c5593SMatthew G. Knepley if (ierr) { 8738c6c5593SMatthew G. Knepley PetscErrorCode ierr2; 87469291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 87569291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 8768c6c5593SMatthew G. Knepley CHKERRQ(ierr); 8778c6c5593SMatthew G. Knepley } 878a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 8795f790a90SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode);CHKERRQ(ierr); 8808c6c5593SMatthew G. Knepley } 881c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr); 882c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 883c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 884c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 8858c6c5593SMatthew G. Knepley } 886c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 88769291d52SBarry Smith ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 88869291d52SBarry Smith ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 88947923291SMatthew G. Knepley } 8908c6c5593SMatthew G. Knepley /* Cleanup */ 891ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 892636b9322SMatthew G. Knepley for (f = 0; f < NfIn; ++f) {ierr = PetscTabulationDestroy(&T[f]);CHKERRQ(ierr);} 893636b9322SMatthew G. Knepley for (f = 0; f < NfAux; ++f) {ierr = PetscTabulationDestroy(&TAux[f]);CHKERRQ(ierr);} 894ef0bb6c7SMatthew G. Knepley ierr = PetscFree2(T, TAux);CHKERRQ(ierr); 8950c898c18SMatthew G. Knepley } 896c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr); 897636b9322SMatthew G. Knepley ierr = PetscFree3(isFE, sp, spIn);CHKERRQ(ierr); 898636b9322SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscFree2(cellsp, cellspIn);CHKERRQ(ierr);} 899a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 900a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plexIn);CHKERRQ(ierr); 901a6e0b375SMatthew G. Knepley if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);} 9028c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 90347923291SMatthew G. Knepley } 9048c6c5593SMatthew G. Knepley 9058c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 9068c6c5593SMatthew G. Knepley { 9078c6c5593SMatthew G. Knepley PetscErrorCode ierr; 9088c6c5593SMatthew G. Knepley 9098c6c5593SMatthew G. Knepley PetscFunctionBegin; 910636b9322SMatthew G. Knepley ierr = DMProjectLocal_Generic_Plex(dm, time, NULL, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 9118c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 9128c6c5593SMatthew G. Knepley } 9138c6c5593SMatthew G. Knepley 9141c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 9158c6c5593SMatthew G. Knepley { 9168c6c5593SMatthew G. Knepley PetscErrorCode ierr; 9178c6c5593SMatthew G. Knepley 9188c6c5593SMatthew G. Knepley PetscFunctionBegin; 919636b9322SMatthew G. Knepley ierr = DMProjectLocal_Generic_Plex(dm, time, NULL, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 92047923291SMatthew G. Knepley PetscFunctionReturn(0); 92147923291SMatthew G. Knepley } 92247923291SMatthew G. Knepley 9238c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 92447923291SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 92547923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 92647923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 927191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 92847923291SMatthew G. Knepley InsertMode mode, Vec localX) 92947923291SMatthew G. Knepley { 93047923291SMatthew G. Knepley PetscErrorCode ierr; 93147923291SMatthew G. Knepley 93247923291SMatthew G. Knepley PetscFunctionBegin; 9331c531cf8SMatthew G. Knepley ierr = DMProjectLocal_Generic_Plex(dm, time, localU, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 9348c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 93547923291SMatthew G. Knepley } 93647923291SMatthew G. Knepley 9371c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 9388c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 9398c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 9408c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 941191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 9428c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 9438c6c5593SMatthew G. Knepley { 9448c6c5593SMatthew G. Knepley PetscErrorCode ierr; 94547923291SMatthew G. Knepley 9468c6c5593SMatthew G. Knepley PetscFunctionBegin; 9471c531cf8SMatthew G. Knepley ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 94847923291SMatthew G. Knepley PetscFunctionReturn(0); 94947923291SMatthew G. Knepley } 950ece3a9fcSMatthew G. Knepley 951ece3a9fcSMatthew G. Knepley PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 952ece3a9fcSMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 953ece3a9fcSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 954ece3a9fcSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 955ece3a9fcSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 956ece3a9fcSMatthew G. Knepley InsertMode mode, Vec localX) 957ece3a9fcSMatthew G. Knepley { 958ece3a9fcSMatthew G. Knepley PetscErrorCode ierr; 959ece3a9fcSMatthew G. Knepley 960ece3a9fcSMatthew G. Knepley PetscFunctionBegin; 961ece3a9fcSMatthew G. Knepley ierr = DMProjectLocal_Generic_Plex(dm, time, localU, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL_BD_FIELD, (void (**)(void)) funcs, NULL, mode, localX);CHKERRQ(ierr); 962ece3a9fcSMatthew G. Knepley PetscFunctionReturn(0); 963ece3a9fcSMatthew G. Knepley } 964