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 101b32699bSMatthew G. Knepley Input Argument: 111b32699bSMatthew G. Knepley . dm - the DM 121b32699bSMatthew G. Knepley 131b32699bSMatthew G. Knepley Output Argument: 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 @*/ 201b32699bSMatthew G. Knepley 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 311b32699bSMatthew G. Knepley Input Arguments: 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 371b32699bSMatthew G. Knepley .seealso: DMPlexGetActivePoint() 381b32699bSMatthew G. Knepley @*/ 391b32699bSMatthew G. Knepley 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 66a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private() 67a6e0b375SMatthew G. Knepley */ 68a6e0b375SMatthew 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[], 698c6c5593SMatthew G. Knepley PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, 708c6c5593SMatthew G. Knepley PetscScalar values[]) 7147923291SMatthew G. Knepley { 72a6e0b375SMatthew G. Knepley PetscInt coordDim, Nf, *Nc, f, spDim, d, v, tp; 7327f02ce8SMatthew G. Knepley PetscBool isAffine, isHybrid, transform; 748c6c5593SMatthew G. Knepley PetscErrorCode ierr; 758c6c5593SMatthew G. Knepley 768c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 77a6e0b375SMatthew G. Knepley ierr = DMGetCoordinateDim(dmIn, &coordDim);CHKERRQ(ierr); 78a6e0b375SMatthew G. Knepley ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr); 79a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 80a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 8127f02ce8SMatthew G. Knepley ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr); 828c6c5593SMatthew G. Knepley /* Get values for closure */ 83c330f8ffSToby Isaac isAffine = fegeom->isAffine; 84c330f8ffSToby Isaac for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 858c6c5593SMatthew G. Knepley void * const ctx = ctxs ? ctxs[f] : NULL; 868c6c5593SMatthew G. Knepley 878c6c5593SMatthew G. Knepley if (!sp[f]) continue; 888c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 898c6c5593SMatthew G. Knepley if (funcs[f]) { 90c330f8ffSToby Isaac if (isFE[f]) { 91c330f8ffSToby Isaac PetscQuadrature allPoints; 92c330f8ffSToby Isaac PetscInt q, dim, numPoints; 93c330f8ffSToby Isaac const PetscReal *points; 94c330f8ffSToby Isaac PetscScalar *pointEval; 95c330f8ffSToby Isaac PetscReal *x; 96ca3d3a14SMatthew G. Knepley DM rdm; 97c330f8ffSToby Isaac 98ca3d3a14SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr); 99b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr); 100c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 101ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 102ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 103c330f8ffSToby Isaac for (q = 0; q < numPoints; q++, tp++) { 104c330f8ffSToby Isaac const PetscReal *v0; 105c330f8ffSToby Isaac 106c330f8ffSToby Isaac if (isAffine) { 107665f567fSMatthew G. Knepley const PetscReal *refpoint = &points[q*dim]; 108665f567fSMatthew G. Knepley PetscReal injpoint[3] = {0., 0., 0.}; 109665f567fSMatthew G. Knepley 110665f567fSMatthew G. Knepley if (dim != fegeom->dim) { 111665f567fSMatthew G. Knepley if (isHybrid) { 112665f567fSMatthew G. Knepley /* We just need to inject into the higher dimensional space assuming the last dimension is collapsed */ 113665f567fSMatthew G. Knepley for (d = 0; d < dim; ++d) injpoint[d] = refpoint[d]; 114665f567fSMatthew G. Knepley refpoint = injpoint; 115665f567fSMatthew G. Knepley } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Reference spatial dimension %D != %D dual basis spatial dimension", fegeom->dim, dim); 116665f567fSMatthew G. Knepley } 117665f567fSMatthew G. Knepley CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, refpoint, x); 118c330f8ffSToby Isaac v0 = x; 1198c6c5593SMatthew G. Knepley } else { 120c330f8ffSToby Isaac v0 = &fegeom->v[tp*coordDim]; 1218c6c5593SMatthew G. Knepley } 122a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;} 123c330f8ffSToby Isaac ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr); 124c330f8ffSToby Isaac } 1254bee2e38SMatthew G. Knepley /* Transform point evaluations pointEval[q,c] */ 1262edcad52SToby Isaac ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr); 127c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 128ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 129ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 130c330f8ffSToby Isaac v += spDim; 13127f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 13227f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 13327f02ce8SMatthew G. Knepley } 134c330f8ffSToby Isaac } else { 135c330f8ffSToby Isaac for (d = 0; d < spDim; ++d, ++v) { 136c330f8ffSToby Isaac ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 137c330f8ffSToby Isaac } 138c330f8ffSToby Isaac } 139c330f8ffSToby Isaac } else { 14027f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 14127f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 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 165ef0bb6c7SMatthew G. Knepley . T - Input basis and derviatives 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 177a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private() 178a6e0b375SMatthew G. Knepley */ 179a6e0b375SMatthew 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, 180ef0bb6c7SMatthew G. Knepley PetscTabulation *T, PetscTabulation *TAux, 1818c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 1828c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1838c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 184191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 1858c6c5593SMatthew G. Knepley PetscScalar values[]) 1868c6c5593SMatthew G. Knepley { 1878c6c5593SMatthew G. Knepley PetscSection section, sectionAux = NULL; 1884bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 1898c6c5593SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 1908c6c5593SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 191191494d9SMatthew G. Knepley const PetscScalar *constants; 1928c6c5593SMatthew G. Knepley PetscReal *x; 193ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc; 1944bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 1954bee2e38SMatthew G. Knepley const PetscInt dE = cgeom->dimEmbed; 196ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0; 19727f02ce8SMatthew G. Knepley PetscBool isAffine, isHybrid, transform; 1988c6c5593SMatthew G. Knepley PetscErrorCode ierr; 1998c6c5593SMatthew G. Knepley 2008c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 201a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 202a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 20327f02ce8SMatthew G. Knepley ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr); 204a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr); 205a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsIn, &uOff);CHKERRQ(ierr); 206a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x);CHKERRQ(ierr); 207a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 208a6e0b375SMatthew G. Knepley ierr = PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 209a6e0b375SMatthew G. Knepley ierr = PetscDSGetConstants(dsIn, &numConstants, &constants);CHKERRQ(ierr); 210a6e0b375SMatthew G. Knepley ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr); 211a6e0b375SMatthew G. Knepley ierr = DMGetLocalSection(dmIn, §ion);CHKERRQ(ierr); 212a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp);CHKERRQ(ierr); 213a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr); 2148c6c5593SMatthew G. Knepley if (dmAux) { 21544171101SMatthew G. Knepley PetscInt subp; 2161ba34526SMatthew G. Knepley 217a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr); 218a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 21992fd8e1eSJed Brown ierr = DMGetLocalSection(dmAux, §ionAux);CHKERRQ(ierr); 220a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 221a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 222a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 2231ba34526SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 2248c6c5593SMatthew G. Knepley } 2258c6c5593SMatthew G. Knepley /* Get values for closure */ 2264bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 22727f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 22827f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 2294bee2e38SMatthew G. Knepley if (isAffine) { 2304bee2e38SMatthew G. Knepley fegeom.v = x; 2314bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2324bee2e38SMatthew G. Knepley fegeom.J = cgeom->J; 2334bee2e38SMatthew G. Knepley fegeom.invJ = cgeom->invJ; 2344bee2e38SMatthew G. Knepley fegeom.detJ = cgeom->detJ; 2354bee2e38SMatthew G. Knepley } 2368c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 237c330f8ffSToby Isaac PetscQuadrature allPoints; 238c330f8ffSToby Isaac PetscInt q, dim, numPoints; 239c330f8ffSToby Isaac const PetscReal *points; 240c330f8ffSToby Isaac PetscScalar *pointEval; 241c330f8ffSToby Isaac DM dm; 242c330f8ffSToby Isaac 2438c6c5593SMatthew G. Knepley if (!sp[f]) continue; 2448c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 245c330f8ffSToby Isaac if (!funcs[f]) { 246be1504a2SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 24727f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 24827f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 24927f02ce8SMatthew G. Knepley } 250c330f8ffSToby Isaac continue; 251c330f8ffSToby Isaac } 252c330f8ffSToby Isaac ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 253b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr); 254c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 255c330f8ffSToby Isaac ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 2560c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 257c330f8ffSToby Isaac if (isAffine) { 258ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, cgeom->dim, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x); 2591c531cf8SMatthew G. Knepley } else { 2604bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[tp*dE]; 2614bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[tp*dE*dE]; 2624bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[tp*dE*dE]; 2634bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[tp]; 2648c6c5593SMatthew G. Knepley } 265ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(dsIn, NfIn, 0, tp, T, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 266ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 267a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);} 268a6e0b375SMatthew 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]); 2691c531cf8SMatthew G. Knepley } 270c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 271c330f8ffSToby Isaac ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 272c330f8ffSToby Isaac v += spDim; 27327f02ce8SMatthew G. Knepley /* TODO: For now, set both sides equal, but this should use info from other support cell */ 27427f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) { 27527f02ce8SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = values[v - spDim]; 27627f02ce8SMatthew G. Knepley } 2778c6c5593SMatthew G. Knepley } 278a6e0b375SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr); 2798c6c5593SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 2808c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 2818c6c5593SMatthew G. Knepley } 2828c6c5593SMatthew G. Knepley 283a6e0b375SMatthew 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, 284ef0bb6c7SMatthew G. Knepley PetscTabulation *T, PetscTabulation *TAux, 285b18799e0SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 286b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 287b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 288b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 289b18799e0SMatthew G. Knepley PetscScalar values[]) 290b18799e0SMatthew G. Knepley { 291b18799e0SMatthew G. Knepley PetscSection section, sectionAux = NULL; 2924bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 293b18799e0SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 294b18799e0SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 295b18799e0SMatthew G. Knepley const PetscScalar *constants; 296b18799e0SMatthew G. Knepley PetscReal *x; 297ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc; 2989f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 2994bee2e38SMatthew G. Knepley const PetscInt dE = fgeom->dimEmbed; 300ef0bb6c7SMatthew G. Knepley PetscInt numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 301b18799e0SMatthew G. Knepley PetscBool isAffine; 302b18799e0SMatthew G. Knepley PetscErrorCode ierr; 303b18799e0SMatthew G. Knepley 304b18799e0SMatthew G. Knepley PetscFunctionBeginHot; 305a6e0b375SMatthew G. Knepley if (dm != dmIn) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM"); 306a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 307a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 308a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 309a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 310a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 311a6e0b375SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 312a6e0b375SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 31392fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 314a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 315b18799e0SMatthew G. Knepley if (dmAux) { 316a6e0b375SMatthew G. Knepley PetscInt subp; 317b18799e0SMatthew G. Knepley 318a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr); 319a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 32092fd8e1eSJed Brown ierr = DMGetLocalSection(dmAux, §ionAux);CHKERRQ(ierr); 321a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 322a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 323a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 324b18799e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 325b18799e0SMatthew G. Knepley } 326b18799e0SMatthew G. Knepley /* Get values for closure */ 3274bee2e38SMatthew G. Knepley isAffine = fgeom->isAffine; 328*ea78f98cSLisandro Dalcin fegeom.n = NULL; 329*ea78f98cSLisandro Dalcin fegeom.J = NULL; 330*ea78f98cSLisandro Dalcin fegeom.v = NULL; 331*ea78f98cSLisandro Dalcin fegeom.xi = NULL; 332a6e0b375SMatthew G. Knepley cgeom.dim = fgeom->dim; 333a6e0b375SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3344bee2e38SMatthew G. Knepley if (isAffine) { 3354bee2e38SMatthew G. Knepley fegeom.v = x; 3364bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3374bee2e38SMatthew G. Knepley fegeom.J = fgeom->J; 3384bee2e38SMatthew G. Knepley fegeom.invJ = fgeom->invJ; 3394bee2e38SMatthew G. Knepley fegeom.detJ = fgeom->detJ; 3404bee2e38SMatthew G. Knepley fegeom.n = fgeom->n; 3419f209ee4SMatthew G. Knepley 3429f209ee4SMatthew G. Knepley cgeom.J = fgeom->suppJ[0]; 3439f209ee4SMatthew G. Knepley cgeom.invJ = fgeom->suppInvJ[0]; 3449f209ee4SMatthew G. Knepley cgeom.detJ = fgeom->suppDetJ[0]; 3454bee2e38SMatthew G. Knepley } 346b18799e0SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 347b18799e0SMatthew G. Knepley PetscQuadrature allPoints; 348b18799e0SMatthew G. Knepley PetscInt q, dim, numPoints; 349b18799e0SMatthew G. Knepley const PetscReal *points; 350b18799e0SMatthew G. Knepley PetscScalar *pointEval; 351b18799e0SMatthew G. Knepley DM dm; 352b18799e0SMatthew G. Knepley 353b18799e0SMatthew G. Knepley if (!sp[f]) continue; 354b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 355b18799e0SMatthew G. Knepley if (!funcs[f]) { 356b18799e0SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 357b18799e0SMatthew G. Knepley continue; 358b18799e0SMatthew G. Knepley } 359b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 360b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f], &allPoints, NULL);CHKERRQ(ierr); 361b18799e0SMatthew G. Knepley ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 362b18799e0SMatthew G. Knepley ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 363b18799e0SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 364b18799e0SMatthew G. Knepley if (isAffine) { 365ef0bb6c7SMatthew G. Knepley CoordinatesRefToReal(dE, fgeom->dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x); 366b18799e0SMatthew G. Knepley } else { 3673fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[tp*dE]; 3689f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[tp*dE*dE]; 3699f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[tp*dE*dE]; 3704bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[tp]; 3714bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[tp*dE]; 3729f209ee4SMatthew G. Knepley 3739f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][tp*dE*dE]; 3749f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][tp*dE*dE]; 3759f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][tp]; 376b18799e0SMatthew G. Knepley } 377a6e0b375SMatthew 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 */ 378ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, tp, T, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 379ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, tp, TAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 3804bee2e38SMatthew 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]); 381b18799e0SMatthew G. Knepley } 382b18799e0SMatthew G. Knepley ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 383b18799e0SMatthew G. Knepley ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 384b18799e0SMatthew G. Knepley v += spDim; 385b18799e0SMatthew G. Knepley } 386a6e0b375SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 387b18799e0SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 388b18799e0SMatthew G. Knepley PetscFunctionReturn(0); 389b18799e0SMatthew G. Knepley } 390b18799e0SMatthew G. Knepley 391a6e0b375SMatthew 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[], 392ef0bb6c7SMatthew G. Knepley PetscDualSpace sp[], PetscInt p, PetscTabulation *T, PetscTabulation *TAux, 3938c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 3948c6c5593SMatthew G. Knepley { 3958c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom; 3968c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed; 3978c6c5593SMatthew G. Knepley PetscErrorCode ierr; 3988c6c5593SMatthew G. Knepley 3998c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 4008c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 4018c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 4028c6c5593SMatthew G. Knepley if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 4038c6c5593SMatthew G. Knepley switch (type) { 4048c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL: 4058c6c5593SMatthew G. Knepley case DM_BC_NATURAL: 406a6e0b375SMatthew 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; 4078c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 4088c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD: 409ef0bb6c7SMatthew G. Knepley ierr = DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, 4100c898c18SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 4110c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 4120c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 413191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 414b18799e0SMatthew G. Knepley case DM_BC_ESSENTIAL_BD_FIELD: 415ef0bb6c7SMatthew G. Knepley ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, T, TAux, 416b18799e0SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 417b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 418b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 419b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 4208c6c5593SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 4218c6c5593SMatthew G. Knepley } 4228c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 4238c6c5593SMatthew G. Knepley } 4248c6c5593SMatthew G. Knepley 425df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 426df5c1128SToby Isaac { 427df5c1128SToby Isaac PetscReal *points; 428df5c1128SToby Isaac PetscInt f, numPoints; 429df5c1128SToby Isaac PetscErrorCode ierr; 430df5c1128SToby Isaac 431df5c1128SToby Isaac PetscFunctionBegin; 432df5c1128SToby Isaac numPoints = 0; 433df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 434df5c1128SToby Isaac if (funcs[f]) { 435df5c1128SToby Isaac PetscQuadrature fAllPoints; 436df5c1128SToby Isaac PetscInt fNumPoints; 437df5c1128SToby Isaac 438b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr); 439df5c1128SToby Isaac ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr); 440df5c1128SToby Isaac numPoints += fNumPoints; 441df5c1128SToby Isaac } 442df5c1128SToby Isaac } 443df5c1128SToby Isaac ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 444df5c1128SToby Isaac numPoints = 0; 445df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 446df5c1128SToby Isaac if (funcs[f]) { 447df5c1128SToby Isaac PetscQuadrature fAllPoints; 44854ee0cceSMatthew G. Knepley PetscInt qdim, fNumPoints, q; 449df5c1128SToby Isaac const PetscReal *fPoints; 450df5c1128SToby Isaac 451b4457527SToby Isaac ierr = PetscDualSpaceGetAllData(sp[f],&fAllPoints, NULL);CHKERRQ(ierr); 45254ee0cceSMatthew G. Knepley ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr); 45354ee0cceSMatthew 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); 454df5c1128SToby Isaac for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q]; 455df5c1128SToby Isaac numPoints += fNumPoints; 456df5c1128SToby Isaac } 457df5c1128SToby Isaac } 458df5c1128SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 459df5c1128SToby Isaac ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 460df5c1128SToby Isaac PetscFunctionReturn(0); 461df5c1128SToby Isaac } 462df5c1128SToby Isaac 463a6e0b375SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *ds) 464e5e52638SMatthew G. Knepley { 465a6e0b375SMatthew G. Knepley DM plex; 466a6e0b375SMatthew G. Knepley DMEnclosureType enc; 467e5e52638SMatthew G. Knepley DMLabel depthLabel; 468e5e52638SMatthew G. Knepley PetscInt dim, cdepth, ls = -1, i; 469e5e52638SMatthew G. Knepley PetscErrorCode ierr; 470e5e52638SMatthew G. Knepley 471e5e52638SMatthew G. Knepley PetscFunctionBegin; 472e5e52638SMatthew G. Knepley if (lStart) *lStart = -1; 473e5e52638SMatthew G. Knepley if (!label) PetscFunctionReturn(0); 474a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dm, odm, &enc);CHKERRQ(ierr); 475e5e52638SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 476a6e0b375SMatthew G. Knepley ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 477a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr); 478e5e52638SMatthew G. Knepley cdepth = dim - height; 479e5e52638SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 480e5e52638SMatthew G. Knepley IS pointIS; 481e5e52638SMatthew G. Knepley const PetscInt *points; 482a6e0b375SMatthew G. Knepley PetscInt pdepth, point; 483e5e52638SMatthew G. Knepley 484e5e52638SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 485e5e52638SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 486e5e52638SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 487a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dm, odm, enc, points[0], &point);CHKERRQ(ierr); 488a6e0b375SMatthew G. Knepley ierr = DMLabelGetValue(depthLabel, point, &pdepth);CHKERRQ(ierr); 489e5e52638SMatthew G. Knepley if (pdepth == cdepth) { 490a6e0b375SMatthew G. Knepley ls = point; 491a6e0b375SMatthew G. Knepley if (ds) {ierr = DMGetCellDS(dm, ls, ds);CHKERRQ(ierr);} 492e5e52638SMatthew G. Knepley } 493e5e52638SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 494e5e52638SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 495e5e52638SMatthew G. Knepley if (ls >= 0) break; 496e5e52638SMatthew G. Knepley } 497e5e52638SMatthew G. Knepley if (lStart) *lStart = ls; 498a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 499e5e52638SMatthew G. Knepley PetscFunctionReturn(0); 500e5e52638SMatthew G. Knepley } 501e5e52638SMatthew G. Knepley 5020de51b56SMatthew G. Knepley /* 5030de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 5040de51b56SMatthew G. Knepley 5050de51b56SMatthew G. Knepley There are several different scenarios: 5060de51b56SMatthew G. Knepley 5070de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data 5080de51b56SMatthew G. Knepley 5090de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells. 5100de51b56SMatthew G. Knepley 5110de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data 5120de51b56SMatthew G. Knepley 5130de51b56SMatthew 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. 5140de51b56SMatthew G. Knepley 5150de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data 5160de51b56SMatthew G. Knepley 5170de51b56SMatthew 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. 5180de51b56SMatthew G. Knepley 5190de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest. 5200de51b56SMatthew G. Knepley 5210de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 5220de51b56SMatthew G. Knepley 5230de51b56SMatthew 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. 5240de51b56SMatthew 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 5255f790a90SMatthew 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 5260de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces. 5270de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 5280de51b56SMatthew G. Knepley 5290de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 5300de51b56SMatthew G. Knepley 5310de51b56SMatthew 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. 5320de51b56SMatthew G. Knepley */ 53346fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 5341c531cf8SMatthew G. Knepley PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], 5358c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 5368c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 5378c6c5593SMatthew G. Knepley { 538a6e0b375SMatthew G. Knepley DM plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm; 539a6e0b375SMatthew G. Knepley DMEnclosureType encIn, encAux; 540a6e0b375SMatthew G. Knepley PetscDS ds = NULL, dsIn = NULL, dsAux = NULL; 541ca3d3a14SMatthew G. Knepley Vec localA = NULL, tv; 542aa0eca99SMatthew G. Knepley IS fieldIS; 54347923291SMatthew G. Knepley PetscSection section; 5448c6c5593SMatthew G. Knepley PetscDualSpace *sp, *cellsp; 545ef0bb6c7SMatthew G. Knepley PetscTabulation *T = NULL, *TAux = NULL; 5468c6c5593SMatthew G. Knepley PetscInt *Nc; 547aa0eca99SMatthew G. Knepley PetscInt dim, dimEmbed, depth, minHeight, maxHeight, h, regionNum, Nf, NfIn, NfAux = 0, NfTot, f; 54827f02ce8SMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, isHybrid = PETSC_FALSE, transform; 549c330f8ffSToby Isaac DMField coordField; 550c330f8ffSToby Isaac DMLabel depthLabel; 551c330f8ffSToby Isaac PetscQuadrature allPoints = NULL; 55247923291SMatthew G. Knepley PetscErrorCode ierr; 55347923291SMatthew G. Knepley 55447923291SMatthew G. Knepley PetscFunctionBegin; 555a6e0b375SMatthew G. Knepley if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);} 556a6e0b375SMatthew G. Knepley else {dmIn = dm;} 55788aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 55888aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr); 559a6e0b375SMatthew G. Knepley ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 560a6e0b375SMatthew G. Knepley ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr); 561a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr); 562a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr); 5638c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 564a6e0b375SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr); 565ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr); 566ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr); 567ca3d3a14SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 5680de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */ 569a6e0b375SMatthew G. Knepley if (dmAux) { 570a6e0b375SMatthew G. Knepley ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr); 571a6e0b375SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 572a6e0b375SMatthew G. Knepley if (!localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector"); 573a6e0b375SMatthew G. Knepley if (!minHeight) { 57488aca1feSMatthew G. Knepley DMLabel spmap; 57588aca1feSMatthew G. Knepley 57688aca1feSMatthew G. Knepley /* If dmAux is a surface, then force the projection to take place over a surface */ 577a6e0b375SMatthew G. Knepley ierr = DMPlexGetSubpointMap(plexAux, &spmap);CHKERRQ(ierr); 578e39fcb4eSStefano Zampini if (spmap) { 579a6e0b375SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(plexAux, &minHeight);CHKERRQ(ierr); 580e39fcb4eSStefano Zampini auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE; 581e39fcb4eSStefano Zampini } 58288aca1feSMatthew G. Knepley } 5830de51b56SMatthew G. Knepley } 584a6e0b375SMatthew G. Knepley } 585a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 586a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr); 587a6e0b375SMatthew G. Knepley ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr); 5882af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight); 589e39fcb4eSStefano Zampini if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim); 590a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr); 591a6e0b375SMatthew G. Knepley if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);} 592a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr); 593a6e0b375SMatthew G. Knepley if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);} 594a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 595a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr); 596aa0eca99SMatthew G. Knepley ierr = DMGetNumFields(dm, &NfTot);CHKERRQ(ierr); 597aa0eca99SMatthew G. Knepley ierr = DMFindRegionNum(dm, ds, ®ionNum);CHKERRQ(ierr); 598aa0eca99SMatthew G. Knepley ierr = DMGetRegionNumDS(dm, regionNum, NULL, &fieldIS, NULL);CHKERRQ(ierr); 59927f02ce8SMatthew G. Knepley ierr = PetscDSGetHybrid(ds, &isHybrid);CHKERRQ(ierr); 6008c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 60192fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 6020c898c18SMatthew G. Knepley if (dmAux) { 603a6e0b375SMatthew G. Knepley ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr); 604a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 6050c898c18SMatthew G. Knepley } 606a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 607496733e2SMatthew G. Knepley ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr); 6088c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);} 6098c6c5593SMatthew G. Knepley else {cellsp = sp;} 610a6e0b375SMatthew G. Knepley if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 6118c6c5593SMatthew G. Knepley /* Get cell dual spaces */ 61247923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 613665f567fSMatthew G. Knepley PetscDiscType disctype; 61447923291SMatthew G. Knepley 615665f567fSMatthew G. Knepley ierr = PetscDSGetDiscType_Internal(ds, f, &disctype);CHKERRQ(ierr); 616665f567fSMatthew G. Knepley if (disctype == PETSC_DISC_FE) { 617665f567fSMatthew G. Knepley PetscFE fe; 61847923291SMatthew G. Knepley 61947923291SMatthew G. Knepley isFE[f] = PETSC_TRUE; 620665f567fSMatthew G. Knepley hasFE = PETSC_TRUE; 621665f567fSMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr); 62247923291SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 623665f567fSMatthew G. Knepley } else if (disctype == PETSC_DISC_FV) { 624665f567fSMatthew G. Knepley PetscFV fv; 62547923291SMatthew G. Knepley 62647923291SMatthew G. Knepley isFE[f] = PETSC_FALSE; 627665f567fSMatthew G. Knepley hasFV = PETSC_TRUE; 628665f567fSMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fv);CHKERRQ(ierr); 62947923291SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 630665f567fSMatthew G. Knepley } else { 631665f567fSMatthew G. Knepley isFE[f] = PETSC_FALSE; 632665f567fSMatthew G. Knepley cellsp[f] = NULL; 633665f567fSMatthew G. Knepley } 63447923291SMatthew G. Knepley } 635c330f8ffSToby Isaac ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr); 636ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 63774fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 63874fc349aSMatthew G. Knepley PetscFE fem, subfem; 639665f567fSMatthew G. Knepley PetscDiscType disctype; 6404a3e9fdbSToby Isaac const PetscReal *points; 6414a3e9fdbSToby Isaac PetscInt numPoints; 6420c898c18SMatthew G. Knepley 6432af58022SMatthew G. Knepley if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 64454ee0cceSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 64554ee0cceSMatthew G. Knepley if (!effectiveHeight) {sp[f] = cellsp[f];} 64654ee0cceSMatthew G. Knepley else {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);} 64754ee0cceSMatthew G. Knepley } 64854ee0cceSMatthew G. Knepley ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr); 649c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 650ef0bb6c7SMatthew G. Knepley ierr = PetscMalloc2(NfIn, &T, NfAux, &TAux);CHKERRQ(ierr); 651a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 652665f567fSMatthew G. Knepley ierr = PetscDSGetDiscType_Internal(dsIn, f, &disctype);CHKERRQ(ierr); 653665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 654a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr); 65574fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 65674fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 657ef0bb6c7SMatthew G. Knepley ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &T[f]);CHKERRQ(ierr); 6580c898c18SMatthew G. Knepley } 6590c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 660665f567fSMatthew G. Knepley ierr = PetscDSGetDiscType_Internal(dsAux, f, &disctype);CHKERRQ(ierr); 661665f567fSMatthew G. Knepley if (disctype != PETSC_DISC_FE) continue; 662a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 66374fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 66474fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 665ef0bb6c7SMatthew G. Knepley ierr = PetscFECreateTabulation(subfem, 1, numPoints, points, 1, &TAux[f]);CHKERRQ(ierr); 6660c898c18SMatthew G. Knepley } 6670c898c18SMatthew G. Knepley } 66847923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 6692af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) { 67088aca1feSMatthew G. Knepley PetscInt effectiveHeight = h - (auxBd ? 0 : minHeight); 671a6e0b375SMatthew G. Knepley PetscDS dsEff = ds; 6728c6c5593SMatthew G. Knepley PetscScalar *values; 673b7260050SToby Isaac PetscBool *fieldActive; 674b7260050SToby Isaac PetscInt maxDegree; 675e5e52638SMatthew G. Knepley PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 676c330f8ffSToby Isaac IS heightIS; 6778c6c5593SMatthew G. Knepley 678a6e0b375SMatthew G. Knepley /* Note we assume that dm and dmIn share the same topology */ 679412e9a14SMatthew G. Knepley ierr = DMPlexGetSimplexOrBoxCells(plex, h, &pStart, &pEnd);CHKERRQ(ierr); 680a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr); 681c330f8ffSToby Isaac ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr); 682c330f8ffSToby Isaac if (pEnd <= pStart) { 683c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 684c330f8ffSToby Isaac continue; 685c330f8ffSToby Isaac } 68647923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */ 68747923291SMatthew G. Knepley totDim = 0; 68847923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 689bec1f4aaSMatthew G. Knepley if (!effectiveHeight) { 69047923291SMatthew G. Knepley sp[f] = cellsp[f]; 69147923291SMatthew G. Knepley } else { 692bec1f4aaSMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr); 69347923291SMatthew G. Knepley } 694665f567fSMatthew G. Knepley if (!sp[f]) continue; 69547923291SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 6969c3cf19fSMatthew G. Knepley totDim += spDim; 69727f02ce8SMatthew G. Knepley if (isHybrid && (f < Nf-1)) totDim += spDim; 69847923291SMatthew G. Knepley } 699a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(plex, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr); 700083401c6SMatthew G. Knepley if (numValues != totDim) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point (%D) closure size %D != dual space dimension %D", lStart < 0 ? pStart : lStart, numValues, totDim); 701c330f8ffSToby Isaac if (!totDim) { 702c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 703c330f8ffSToby Isaac continue; 704c330f8ffSToby Isaac } 705a6e0b375SMatthew G. Knepley if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(ds, effectiveHeight, &dsEff);CHKERRQ(ierr);} 70647923291SMatthew G. Knepley /* Loop over points at this height */ 70769291d52SBarry Smith ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 708aa0eca99SMatthew G. Knepley ierr = DMGetWorkArray(dm, NfTot, MPI_INT, &fieldActive);CHKERRQ(ierr); 709aa0eca99SMatthew G. Knepley { 710aa0eca99SMatthew G. Knepley const PetscInt *fields; 711aa0eca99SMatthew G. Knepley 712aa0eca99SMatthew G. Knepley ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr); 713aa0eca99SMatthew G. Knepley for (f = 0; f < NfTot; ++f) {fieldActive[f] = PETSC_FALSE;} 714aa0eca99SMatthew G. Knepley for (f = 0; f < Nf; ++f) {fieldActive[fields[f]] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE;} 715aa0eca99SMatthew G. Knepley ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr); 716aa0eca99SMatthew G. Knepley } 7178c6c5593SMatthew G. Knepley if (label) { 7188c6c5593SMatthew G. Knepley PetscInt i; 71947923291SMatthew G. Knepley 72047923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 721c330f8ffSToby Isaac IS pointIS, isectIS; 72247923291SMatthew G. Knepley const PetscInt *points; 7238c6c5593SMatthew G. Knepley PetscInt n; 724c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 725c330f8ffSToby Isaac PetscQuadrature quad = NULL; 72647923291SMatthew G. Knepley 72747923291SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 72847923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 729c330f8ffSToby Isaac ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr); 730c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 731c330f8ffSToby Isaac if (!isectIS) continue; 732c330f8ffSToby Isaac ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr); 733c330f8ffSToby Isaac ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr); 734b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr); 735b7260050SToby Isaac if (maxDegree <= 1) { 736c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr); 737c330f8ffSToby Isaac } 738c330f8ffSToby Isaac if (!quad) { 739c330f8ffSToby Isaac if (!h && allPoints) { 740c330f8ffSToby Isaac quad = allPoints; 741c330f8ffSToby Isaac allPoints = NULL; 742c330f8ffSToby Isaac } else { 74327f02ce8SMatthew G. Knepley ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,isHybrid ? dim-h-1 : dim-h,funcs,&quad);CHKERRQ(ierr); 744c330f8ffSToby Isaac } 745c330f8ffSToby Isaac } 746a6e0b375SMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr); 74747923291SMatthew G. Knepley for (p = 0; p < n; ++p) { 74847923291SMatthew G. Knepley const PetscInt point = points[p]; 74947923291SMatthew G. Knepley 750580bdb30SBarry Smith ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr); 751c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 7521b32699bSMatthew G. Knepley ierr = DMPlexSetActivePoint(dm, point);CHKERRQ(ierr); 753ef0bb6c7SMatthew G. Knepley ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, T, TAux, type, funcs, ctxs, fieldActive, values); 75447923291SMatthew G. Knepley if (ierr) { 75547923291SMatthew G. Knepley PetscErrorCode ierr2; 75669291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 75769291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 75847923291SMatthew G. Knepley CHKERRQ(ierr); 75947923291SMatthew G. Knepley } 760a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 7615f790a90SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, label, ids[i], values, mode);CHKERRQ(ierr); 76247923291SMatthew G. Knepley } 763c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 764c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 765c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 766c330f8ffSToby Isaac ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr); 767c330f8ffSToby Isaac ierr = ISDestroy(&isectIS);CHKERRQ(ierr); 76847923291SMatthew G. Knepley } 7698c6c5593SMatthew G. Knepley } else { 770c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 771c330f8ffSToby Isaac PetscQuadrature quad = NULL; 772c330f8ffSToby Isaac IS pointIS; 773c330f8ffSToby Isaac 774c330f8ffSToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr); 775b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 776b7260050SToby Isaac if (maxDegree <= 1) { 777c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr); 778c330f8ffSToby Isaac } 779c330f8ffSToby Isaac if (!quad) { 780c330f8ffSToby Isaac if (!h && allPoints) { 781c330f8ffSToby Isaac quad = allPoints; 782c330f8ffSToby Isaac allPoints = NULL; 783c330f8ffSToby Isaac } else { 784c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 785c330f8ffSToby Isaac } 786c330f8ffSToby Isaac } 787a6e0b375SMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr); 7888c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 789580bdb30SBarry Smith ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr); 790c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr); 7911b32699bSMatthew G. Knepley ierr = DMPlexSetActivePoint(dm, p);CHKERRQ(ierr); 792ef0bb6c7SMatthew G. Knepley ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, T, TAux, type, funcs, ctxs, fieldActive, values); 7938c6c5593SMatthew G. Knepley if (ierr) { 7948c6c5593SMatthew G. Knepley PetscErrorCode ierr2; 79569291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 79669291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 7978c6c5593SMatthew G. Knepley CHKERRQ(ierr); 7988c6c5593SMatthew G. Knepley } 799a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 8005f790a90SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, NULL, -1, values, mode);CHKERRQ(ierr); 8018c6c5593SMatthew G. Knepley } 802c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr); 803c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 804c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 805c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 8068c6c5593SMatthew G. Knepley } 807c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 80869291d52SBarry Smith ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 80969291d52SBarry Smith ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 81047923291SMatthew G. Knepley } 8118c6c5593SMatthew G. Knepley /* Cleanup */ 812ece3a9fcSMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 81374fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 81474fc349aSMatthew G. Knepley PetscFE fem, subfem; 8150c898c18SMatthew G. Knepley 816a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 817a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr); 81874fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 81974fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 820ef0bb6c7SMatthew G. Knepley ierr = PetscTabulationDestroy(&T[f]);CHKERRQ(ierr); 8210c898c18SMatthew G. Knepley } 8220c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 823a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 82474fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 82574fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 826ef0bb6c7SMatthew G. Knepley ierr = PetscTabulationDestroy(&TAux[f]);CHKERRQ(ierr); 8270c898c18SMatthew G. Knepley } 828ef0bb6c7SMatthew G. Knepley ierr = PetscFree2(T, TAux);CHKERRQ(ierr); 8290c898c18SMatthew G. Knepley } 830c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr); 831496733e2SMatthew G. Knepley ierr = PetscFree2(isFE, sp);CHKERRQ(ierr); 8328c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);} 833a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 834a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plexIn);CHKERRQ(ierr); 835a6e0b375SMatthew G. Knepley if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);} 8368c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 83747923291SMatthew G. Knepley } 8388c6c5593SMatthew G. Knepley 8398c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 8408c6c5593SMatthew G. Knepley { 8418c6c5593SMatthew G. Knepley PetscErrorCode ierr; 8428c6c5593SMatthew G. Knepley 8438c6c5593SMatthew G. Knepley PetscFunctionBegin; 8441c531cf8SMatthew G. Knepley ierr = DMProjectLocal_Generic_Plex(dm, time, localX, 0, NULL, NULL, 0, NULL, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 8458c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 8468c6c5593SMatthew G. Knepley } 8478c6c5593SMatthew G. Knepley 8481c531cf8SMatthew 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) 8498c6c5593SMatthew G. Knepley { 8508c6c5593SMatthew G. Knepley PetscErrorCode ierr; 8518c6c5593SMatthew G. Knepley 8528c6c5593SMatthew G. Knepley PetscFunctionBegin; 8531c531cf8SMatthew G. Knepley ierr = DMProjectLocal_Generic_Plex(dm, time, localX, Ncc, comps, label, numIds, ids, DM_BC_ESSENTIAL, (void (**)(void)) funcs, ctxs, mode, localX);CHKERRQ(ierr); 85447923291SMatthew G. Knepley PetscFunctionReturn(0); 85547923291SMatthew G. Knepley } 85647923291SMatthew G. Knepley 8578c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 85847923291SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 85947923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 86047923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 861191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 86247923291SMatthew G. Knepley InsertMode mode, Vec localX) 86347923291SMatthew G. Knepley { 86447923291SMatthew G. Knepley PetscErrorCode ierr; 86547923291SMatthew G. Knepley 86647923291SMatthew G. Knepley PetscFunctionBegin; 8671c531cf8SMatthew 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); 8688c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 86947923291SMatthew G. Knepley } 87047923291SMatthew G. Knepley 8711c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 8728c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 8738c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8748c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 875191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8768c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 8778c6c5593SMatthew G. Knepley { 8788c6c5593SMatthew G. Knepley PetscErrorCode ierr; 87947923291SMatthew G. Knepley 8808c6c5593SMatthew G. Knepley PetscFunctionBegin; 8811c531cf8SMatthew 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); 88247923291SMatthew G. Knepley PetscFunctionReturn(0); 88347923291SMatthew G. Knepley } 884ece3a9fcSMatthew G. Knepley 885ece3a9fcSMatthew G. Knepley PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 886ece3a9fcSMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 887ece3a9fcSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 888ece3a9fcSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 889ece3a9fcSMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 890ece3a9fcSMatthew G. Knepley InsertMode mode, Vec localX) 891ece3a9fcSMatthew G. Knepley { 892ece3a9fcSMatthew G. Knepley PetscErrorCode ierr; 893ece3a9fcSMatthew G. Knepley 894ece3a9fcSMatthew G. Knepley PetscFunctionBegin; 895ece3a9fcSMatthew 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); 896ece3a9fcSMatthew G. Knepley PetscFunctionReturn(0); 897ece3a9fcSMatthew G. Knepley } 898