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 5*a6e0b375SMatthew G. Knepley /* 6*a6e0b375SMatthew G. Knepley DMProjectPoint_Func_Private - Interpolate the given function in the output basis on the given point 7*a6e0b375SMatthew G. Knepley 8*a6e0b375SMatthew G. Knepley Input Parameters: 9*a6e0b375SMatthew G. Knepley + dm - The output DM 10*a6e0b375SMatthew G. Knepley . ds - The output DS 11*a6e0b375SMatthew G. Knepley . dmIn - The input DM 12*a6e0b375SMatthew G. Knepley . dsIn - The input DS 13*a6e0b375SMatthew G. Knepley . time - The time for this evaluation 14*a6e0b375SMatthew G. Knepley . fegeom - The FE geometry for this point 15*a6e0b375SMatthew G. Knepley . fvgeom - The FV geometry for this point 16*a6e0b375SMatthew G. Knepley . isFE - Flag indicating whether each output field has an FE discretization 17*a6e0b375SMatthew G. Knepley . sp - The output PetscDualSpace for each field 18*a6e0b375SMatthew G. Knepley . funcs - The evaluation function for each field 19*a6e0b375SMatthew G. Knepley - ctxs - The user context for each field 20*a6e0b375SMatthew G. Knepley 21*a6e0b375SMatthew G. Knepley Output Parameter: 22*a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space 23*a6e0b375SMatthew G. Knepley 24*a6e0b375SMatthew G. Knepley Level: developer 25*a6e0b375SMatthew G. Knepley 26*a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private() 27*a6e0b375SMatthew G. Knepley */ 28*a6e0b375SMatthew 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[], 298c6c5593SMatthew G. Knepley PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, 308c6c5593SMatthew G. Knepley PetscScalar values[]) 3147923291SMatthew G. Knepley { 32*a6e0b375SMatthew G. Knepley PetscInt coordDim, Nf, *Nc, f, spDim, d, v, tp; 33ca3d3a14SMatthew G. Knepley PetscBool isAffine, transform; 348c6c5593SMatthew G. Knepley PetscErrorCode ierr; 358c6c5593SMatthew G. Knepley 368c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 37*a6e0b375SMatthew G. Knepley ierr = DMGetCoordinateDim(dmIn, &coordDim);CHKERRQ(ierr); 38*a6e0b375SMatthew G. Knepley ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr); 39*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 40*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 418c6c5593SMatthew G. Knepley /* Get values for closure */ 42c330f8ffSToby Isaac isAffine = fegeom->isAffine; 43c330f8ffSToby Isaac for (f = 0, v = 0, tp = 0; f < Nf; ++f) { 448c6c5593SMatthew G. Knepley void * const ctx = ctxs ? ctxs[f] : NULL; 458c6c5593SMatthew G. Knepley 468c6c5593SMatthew G. Knepley if (!sp[f]) continue; 478c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 488c6c5593SMatthew G. Knepley if (funcs[f]) { 49c330f8ffSToby Isaac if (isFE[f]) { 50c330f8ffSToby Isaac PetscQuadrature allPoints; 51c330f8ffSToby Isaac PetscInt q, dim, numPoints; 52c330f8ffSToby Isaac const PetscReal *points; 53c330f8ffSToby Isaac PetscScalar *pointEval; 54c330f8ffSToby Isaac PetscReal *x; 55ca3d3a14SMatthew G. Knepley DM rdm; 56c330f8ffSToby Isaac 57ca3d3a14SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&rdm);CHKERRQ(ierr); 58c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 59c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 60ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 61ca3d3a14SMatthew G. Knepley ierr = DMGetWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 62c330f8ffSToby Isaac for (q = 0; q < numPoints; q++, tp++) { 63c330f8ffSToby Isaac const PetscReal *v0; 64c330f8ffSToby Isaac 65c330f8ffSToby Isaac if (isAffine) { 66c330f8ffSToby Isaac CoordinatesRefToReal(coordDim, fegeom->dim, fegeom->xi, fegeom->v, fegeom->J, &points[q*dim], x); 67c330f8ffSToby Isaac v0 = x; 688c6c5593SMatthew G. Knepley } else { 69c330f8ffSToby Isaac v0 = &fegeom->v[tp*coordDim]; 708c6c5593SMatthew G. Knepley } 71*a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, v0, PETSC_TRUE, coordDim, v0, x, dm->transformCtx);CHKERRQ(ierr); v0 = x;} 72c330f8ffSToby Isaac ierr = (*funcs[f])(coordDim, time, v0, Nc[f], &pointEval[Nc[f]*q], ctx);CHKERRQ(ierr); 73c330f8ffSToby Isaac } 744bee2e38SMatthew G. Knepley /* Transform point evaluations pointEval[q,c] */ 754bee2e38SMatthew G. Knepley ierr = PetscDualSpacePullback(sp[f], fegeom, numPoints, Nc[f], pointEval);CHKERRQ(ierr); 76c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 77ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,coordDim,MPIU_REAL,&x);CHKERRQ(ierr); 78ca3d3a14SMatthew G. Knepley ierr = DMRestoreWorkArray(rdm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 79c330f8ffSToby Isaac v += spDim; 80c330f8ffSToby Isaac } else { 81c330f8ffSToby Isaac for (d = 0; d < spDim; ++d, ++v) { 82c330f8ffSToby Isaac ierr = PetscDualSpaceApplyFVM(sp[f], d, time, fvgeom, Nc[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 83c330f8ffSToby Isaac } 84c330f8ffSToby Isaac } 85c330f8ffSToby Isaac } else { 86c330f8ffSToby Isaac for (d = 0; d < spDim; d++, v++) {values[v] = 0.0;} 878c6c5593SMatthew G. Knepley } 889c3cf19fSMatthew G. Knepley } 898c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 908c6c5593SMatthew G. Knepley } 918c6c5593SMatthew G. Knepley 92*a6e0b375SMatthew G. Knepley /* 93*a6e0b375SMatthew G. Knepley DMProjectPoint_Field_Private - Interpolate a function of the given field, in the input basis, using the output basis on the given point 94*a6e0b375SMatthew G. Knepley 95*a6e0b375SMatthew G. Knepley Input Parameters: 96*a6e0b375SMatthew G. Knepley + dm - The output DM 97*a6e0b375SMatthew G. Knepley . ds - The output DS 98*a6e0b375SMatthew G. Knepley . dmIn - The input DM 99*a6e0b375SMatthew G. Knepley . dsIn - The input DS 100*a6e0b375SMatthew G. Knepley . dmAux - The auxiliary DM, which is always for the input space 101*a6e0b375SMatthew G. Knepley . dsAux - The auxiliary DS, which is always for the input space 102*a6e0b375SMatthew G. Knepley . time - The time for this evaluation 103*a6e0b375SMatthew G. Knepley . localU - The local solution 104*a6e0b375SMatthew G. Knepley . localA - The local auziliary fields 105*a6e0b375SMatthew G. Knepley . cgeom - The FE geometry for this point 106*a6e0b375SMatthew G. Knepley . sp - The output PetscDualSpace for each field 107*a6e0b375SMatthew G. Knepley . p - The point in the output DM 108*a6e0b375SMatthew G. Knepley . basisTab - Input basis for each field tabulated on the quadrature points 109*a6e0b375SMatthew G. Knepley . basisDerTab - Input basis derivatives for each field tabulated on the quadrature points 110*a6e0b375SMatthew G. Knepley . basisTabAux - Auxiliary basis for each aux field tabulated on the quadrature points 111*a6e0b375SMatthew G. Knepley . basisDerTabAux - Auxiliary basis for each aux field tabulated on the quadrature points 112*a6e0b375SMatthew G. Knepley . funcs - The evaluation function for each field 113*a6e0b375SMatthew G. Knepley - ctxs - The user context for each field 114*a6e0b375SMatthew G. Knepley 115*a6e0b375SMatthew G. Knepley Output Parameter: 116*a6e0b375SMatthew G. Knepley . values - The value for each dual basis vector in the output dual space 117*a6e0b375SMatthew G. Knepley 118*a6e0b375SMatthew G. Knepley Note: Not supported for FV 119*a6e0b375SMatthew G. Knepley 120*a6e0b375SMatthew G. Knepley Level: developer 121*a6e0b375SMatthew G. Knepley 122*a6e0b375SMatthew G. Knepley .seealso: DMProjectPoint_Field_Private() 123*a6e0b375SMatthew G. Knepley */ 124*a6e0b375SMatthew 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, 1250c898c18SMatthew G. Knepley PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 1268c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 1278c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1288c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 129191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 1308c6c5593SMatthew G. Knepley PetscScalar values[]) 1318c6c5593SMatthew G. Knepley { 1328c6c5593SMatthew G. Knepley PetscSection section, sectionAux = NULL; 1334bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 1348c6c5593SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 1358c6c5593SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 136191494d9SMatthew G. Knepley const PetscScalar *constants; 1378c6c5593SMatthew G. Knepley PetscReal *x; 138*a6e0b375SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nc, *NbIn, *NcIn, *NbAux = NULL, *NcAux = NULL; 1394bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 1404bee2e38SMatthew G. Knepley const PetscInt dE = cgeom->dimEmbed; 141*a6e0b375SMatthew G. Knepley PetscInt dimIn, dimAux = 0, numConstants, Nf, NfIn, NfAux = 0, f, spDim, d, v, inp, tp = 0; 142ca3d3a14SMatthew G. Knepley PetscBool isAffine, transform; 1438c6c5593SMatthew G. Knepley PetscErrorCode ierr; 1448c6c5593SMatthew G. Knepley 1458c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 146*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 147*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 148*a6e0b375SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsIn, &dimIn);CHKERRQ(ierr); 149*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr); 150*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDimensions(dsIn, &NbIn);CHKERRQ(ierr); 151*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(dsIn, &NcIn);CHKERRQ(ierr); 152*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsIn, &uOff);CHKERRQ(ierr); 153*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsIn, &uOff_x);CHKERRQ(ierr); 154*a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsIn, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 155*a6e0b375SMatthew G. Knepley ierr = PetscDSGetWorkspace(dsIn, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 156*a6e0b375SMatthew G. Knepley ierr = PetscDSGetConstants(dsIn, &numConstants, &constants);CHKERRQ(ierr); 157*a6e0b375SMatthew G. Knepley ierr = DMHasBasisTransform(dmIn, &transform);CHKERRQ(ierr); 158*a6e0b375SMatthew G. Knepley ierr = DMGetLocalSection(dmIn, §ion);CHKERRQ(ierr); 159*a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmIn, dm, encIn, p, &inp);CHKERRQ(ierr); 160*a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr); 1618c6c5593SMatthew G. Knepley if (dmAux) { 16244171101SMatthew G. Knepley PetscInt subp; 1631ba34526SMatthew G. Knepley 164*a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr); 165*a6e0b375SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 166*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 167*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 168*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 16992fd8e1eSJed Brown ierr = DMGetLocalSection(dmAux, §ionAux);CHKERRQ(ierr); 170*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 171*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 172*a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 1731ba34526SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 1748c6c5593SMatthew G. Knepley } 1758c6c5593SMatthew G. Knepley /* Get values for closure */ 1764bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 1774bee2e38SMatthew G. Knepley if (isAffine) { 1784bee2e38SMatthew G. Knepley fegeom.v = x; 1794bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 1804bee2e38SMatthew G. Knepley fegeom.J = cgeom->J; 1814bee2e38SMatthew G. Knepley fegeom.invJ = cgeom->invJ; 1824bee2e38SMatthew G. Knepley fegeom.detJ = cgeom->detJ; 1834bee2e38SMatthew G. Knepley } 1848c6c5593SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 185c330f8ffSToby Isaac PetscQuadrature allPoints; 186c330f8ffSToby Isaac PetscInt q, dim, numPoints; 187c330f8ffSToby Isaac const PetscReal *points; 188c330f8ffSToby Isaac PetscScalar *pointEval; 189c330f8ffSToby Isaac DM dm; 190c330f8ffSToby Isaac 1918c6c5593SMatthew G. Knepley if (!sp[f]) continue; 1928c6c5593SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 193c330f8ffSToby Isaac if (!funcs[f]) { 194be1504a2SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 195c330f8ffSToby Isaac continue; 196c330f8ffSToby Isaac } 197c330f8ffSToby Isaac ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 198c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 199c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 200c330f8ffSToby Isaac ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 2010c898c18SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 202c330f8ffSToby Isaac if (isAffine) { 203*a6e0b375SMatthew G. Knepley CoordinatesRefToReal(dE, dimIn, fegeom.xi, cgeom->v, fegeom.J, &points[q*dim], x); 2041c531cf8SMatthew G. Knepley } else { 2054bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[tp*dE]; 2064bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[tp*dE*dE]; 2074bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[tp*dE*dE]; 2084bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[tp]; 2098c6c5593SMatthew G. Knepley } 210*a6e0b375SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(dsIn, dimIn, NfIn, NbIn, NcIn, tp, basisTab, basisDerTab, &fegeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 211*a6e0b375SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &fegeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 212*a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformApplyReal_Internal(dmIn, fegeom.v, PETSC_TRUE, dE, fegeom.v, fegeom.v, dm->transformCtx);CHKERRQ(ierr);} 213*a6e0b375SMatthew 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]); 2141c531cf8SMatthew G. Knepley } 215c330f8ffSToby Isaac ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 216c330f8ffSToby Isaac ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 217c330f8ffSToby Isaac v += spDim; 2188c6c5593SMatthew G. Knepley } 219*a6e0b375SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmIn, section, localU, inp, NULL, &coefficients);CHKERRQ(ierr); 2208c6c5593SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 2218c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 2228c6c5593SMatthew G. Knepley } 2238c6c5593SMatthew G. Knepley 224*a6e0b375SMatthew 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, 225b18799e0SMatthew G. Knepley PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 226b18799e0SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 227b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 228b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 229b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void **ctxs, 230b18799e0SMatthew G. Knepley PetscScalar values[]) 231b18799e0SMatthew G. Knepley { 232b18799e0SMatthew G. Knepley PetscSection section, sectionAux = NULL; 2334bee2e38SMatthew G. Knepley PetscScalar *u, *u_t = NULL, *u_x, *a = NULL, *a_t = NULL, *a_x = NULL, *bc; 234b18799e0SMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 235b18799e0SMatthew G. Knepley PetscScalar *coefficients_t = NULL, *coefficientsAux_t = NULL; 236b18799e0SMatthew G. Knepley const PetscScalar *constants; 237b18799e0SMatthew G. Knepley PetscReal *x; 238b18799e0SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 2399f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 2404bee2e38SMatthew G. Knepley const PetscInt dE = fgeom->dimEmbed; 241b18799e0SMatthew G. Knepley PetscInt dimAux = 0, numConstants, Nf, NfAux = 0, f, spDim, d, v, tp = 0; 242b18799e0SMatthew G. Knepley PetscBool isAffine; 243b18799e0SMatthew G. Knepley PetscErrorCode ierr; 244b18799e0SMatthew G. Knepley 245b18799e0SMatthew G. Knepley PetscFunctionBeginHot; 246*a6e0b375SMatthew G. Knepley if (dm != dmIn) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Not yet upgraded to use different input DM"); 247*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 248*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 249*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 250*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 251*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 252*a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, &bc /*&u_t*/, &u_x);CHKERRQ(ierr); 253*a6e0b375SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 254*a6e0b375SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 25592fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 256*a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 257b18799e0SMatthew G. Knepley if (dmAux) { 258*a6e0b375SMatthew G. Knepley PetscInt subp; 259b18799e0SMatthew G. Knepley 260*a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dmAux, dm, encAux, p, &subp);CHKERRQ(ierr); 261*a6e0b375SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 262*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 263*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 264*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 26592fd8e1eSJed Brown ierr = DMGetLocalSection(dmAux, §ionAux);CHKERRQ(ierr); 266*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 267*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 268*a6e0b375SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL /*&a_t*/, &a_x);CHKERRQ(ierr); 269b18799e0SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, localA, subp, NULL, &coefficientsAux);CHKERRQ(ierr); 270b18799e0SMatthew G. Knepley } 271b18799e0SMatthew G. Knepley /* Get values for closure */ 2724bee2e38SMatthew G. Knepley isAffine = fgeom->isAffine; 2734dcf62a8SSatish Balay fegeom.n = 0; 2744dcf62a8SSatish Balay fegeom.J = 0; 2754dcf62a8SSatish Balay fegeom.v = 0; 2764dcf62a8SSatish Balay fegeom.xi = 0; 277*a6e0b375SMatthew G. Knepley cgeom.dim = fgeom->dim; 278*a6e0b375SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 2794bee2e38SMatthew G. Knepley if (isAffine) { 2804bee2e38SMatthew G. Knepley fegeom.v = x; 2814bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 2824bee2e38SMatthew G. Knepley fegeom.J = fgeom->J; 2834bee2e38SMatthew G. Knepley fegeom.invJ = fgeom->invJ; 2844bee2e38SMatthew G. Knepley fegeom.detJ = fgeom->detJ; 2854bee2e38SMatthew G. Knepley fegeom.n = fgeom->n; 2869f209ee4SMatthew G. Knepley 2879f209ee4SMatthew G. Knepley cgeom.J = fgeom->suppJ[0]; 2889f209ee4SMatthew G. Knepley cgeom.invJ = fgeom->suppInvJ[0]; 2899f209ee4SMatthew G. Knepley cgeom.detJ = fgeom->suppDetJ[0]; 2904bee2e38SMatthew G. Knepley } 291b18799e0SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 292b18799e0SMatthew G. Knepley PetscQuadrature allPoints; 293b18799e0SMatthew G. Knepley PetscInt q, dim, numPoints; 294b18799e0SMatthew G. Knepley const PetscReal *points; 295b18799e0SMatthew G. Knepley PetscScalar *pointEval; 296b18799e0SMatthew G. Knepley DM dm; 297b18799e0SMatthew G. Knepley 298b18799e0SMatthew G. Knepley if (!sp[f]) continue; 299b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 300b18799e0SMatthew G. Knepley if (!funcs[f]) { 301b18799e0SMatthew G. Knepley for (d = 0; d < spDim; d++, v++) values[v] = 0.; 302b18799e0SMatthew G. Knepley continue; 303b18799e0SMatthew G. Knepley } 304b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp[f],&dm);CHKERRQ(ierr); 305b18799e0SMatthew G. Knepley ierr = PetscDualSpaceGetAllPoints(sp[f], &allPoints);CHKERRQ(ierr); 306b18799e0SMatthew G. Knepley ierr = PetscQuadratureGetData(allPoints,&dim,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 307b18799e0SMatthew G. Knepley ierr = DMGetWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 308b18799e0SMatthew G. Knepley for (q = 0; q < numPoints; ++q, ++tp) { 309b18799e0SMatthew G. Knepley if (isAffine) { 3104bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, fgeom->v, fegeom.J, &points[q*dim], x); 311b18799e0SMatthew G. Knepley } else { 3123fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[tp*dE]; 3139f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[tp*dE*dE]; 3149f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[tp*dE*dE]; 3154bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[tp]; 3164bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[tp*dE]; 3179f209ee4SMatthew G. Knepley 3189f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][tp*dE*dE]; 3199f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][tp*dE*dE]; 3209f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][tp]; 321b18799e0SMatthew G. Knepley } 322*a6e0b375SMatthew 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 */ 323*a6e0b375SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, tp, basisTab, basisDerTab, &cgeom, coefficients, coefficients_t, u, u_x, u_t);CHKERRQ(ierr); 324*a6e0b375SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dimAux, NfAux, NbAux, NcAux, tp, basisTabAux, basisDerTabAux, &cgeom, coefficientsAux, coefficientsAux_t, a, a_x, a_t);CHKERRQ(ierr);} 3254bee2e38SMatthew 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]); 326b18799e0SMatthew G. Knepley } 327b18799e0SMatthew G. Knepley ierr = PetscDualSpaceApplyAll(sp[f], pointEval, &values[v]);CHKERRQ(ierr); 328b18799e0SMatthew G. Knepley ierr = DMRestoreWorkArray(dm,numPoints*Nc[f],MPIU_SCALAR,&pointEval);CHKERRQ(ierr); 329b18799e0SMatthew G. Knepley v += spDim; 330b18799e0SMatthew G. Knepley } 331*a6e0b375SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmIn, section, localU, p, NULL, &coefficients);CHKERRQ(ierr); 332b18799e0SMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, localA, p, NULL, &coefficientsAux);CHKERRQ(ierr);} 333b18799e0SMatthew G. Knepley PetscFunctionReturn(0); 334b18799e0SMatthew G. Knepley } 335b18799e0SMatthew G. Knepley 336*a6e0b375SMatthew 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[], 337*a6e0b375SMatthew G. Knepley PetscDualSpace sp[], PetscInt p, 3381c531cf8SMatthew G. Knepley PetscReal **basisTab, PetscReal **basisDerTab, PetscReal **basisTabAux, PetscReal **basisDerTabAux, 3398c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, PetscBool fieldActive[], PetscScalar values[]) 3408c6c5593SMatthew G. Knepley { 3418c6c5593SMatthew G. Knepley PetscFVCellGeom fvgeom; 3428c6c5593SMatthew G. Knepley PetscInt dim, dimEmbed; 3438c6c5593SMatthew G. Knepley PetscErrorCode ierr; 3448c6c5593SMatthew G. Knepley 3458c6c5593SMatthew G. Knepley PetscFunctionBeginHot; 3468c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 3478c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 3488c6c5593SMatthew G. Knepley if (hasFV) {ierr = DMPlexComputeCellGeometryFVM(dm, p, &fvgeom.volume, fvgeom.centroid, NULL);CHKERRQ(ierr);} 3498c6c5593SMatthew G. Knepley switch (type) { 3508c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL: 3518c6c5593SMatthew G. Knepley case DM_BC_NATURAL: 352*a6e0b375SMatthew 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; 3538c6c5593SMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 3548c6c5593SMatthew G. Knepley case DM_BC_NATURAL_FIELD: 355*a6e0b375SMatthew G. Knepley ierr = DMProjectPoint_Field_Private(dm, ds, dmIn, encIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, 3560c898c18SMatthew G. Knepley basisTab, basisDerTab, basisTabAux, basisDerTabAux, 3570c898c18SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 3580c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 3590c898c18SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 360191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 361b18799e0SMatthew G. Knepley case DM_BC_ESSENTIAL_BD_FIELD: 362*a6e0b375SMatthew G. Knepley ierr = DMProjectPoint_BdField_Private(dm, ds, dmIn, dsIn, dmAux, encAux, dsAux, time, localU, localA, fegeom, sp, p, 363b18799e0SMatthew G. Knepley basisTab, basisDerTab, basisTabAux, basisDerTabAux, 364b18799e0SMatthew G. Knepley (void (**)(PetscInt, PetscInt, PetscInt, 365b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 366b18799e0SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 367b18799e0SMatthew G. Knepley PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])) funcs, ctxs, values);CHKERRQ(ierr);break; 3688c6c5593SMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown boundary condition type: %d", (int) type); 3698c6c5593SMatthew G. Knepley } 3708c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 3718c6c5593SMatthew G. Knepley } 3728c6c5593SMatthew G. Knepley 373df5c1128SToby Isaac static PetscErrorCode PetscDualSpaceGetAllPointsUnion(PetscInt Nf, PetscDualSpace *sp, PetscInt dim, void (**funcs)(void), PetscQuadrature *allPoints) 374df5c1128SToby Isaac { 375df5c1128SToby Isaac PetscReal *points; 376df5c1128SToby Isaac PetscInt f, numPoints; 377df5c1128SToby Isaac PetscErrorCode ierr; 378df5c1128SToby Isaac 379df5c1128SToby Isaac PetscFunctionBegin; 380df5c1128SToby Isaac numPoints = 0; 381df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 382df5c1128SToby Isaac if (funcs[f]) { 383df5c1128SToby Isaac PetscQuadrature fAllPoints; 384df5c1128SToby Isaac PetscInt fNumPoints; 385df5c1128SToby Isaac 386df5c1128SToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 387df5c1128SToby Isaac ierr = PetscQuadratureGetData(fAllPoints, NULL, NULL, &fNumPoints, NULL, NULL);CHKERRQ(ierr); 388df5c1128SToby Isaac numPoints += fNumPoints; 389df5c1128SToby Isaac } 390df5c1128SToby Isaac } 391df5c1128SToby Isaac ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 392df5c1128SToby Isaac numPoints = 0; 393df5c1128SToby Isaac for (f = 0; f < Nf; ++f) { 394df5c1128SToby Isaac PetscInt spDim; 395df5c1128SToby Isaac 396df5c1128SToby Isaac ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 397df5c1128SToby Isaac if (funcs[f]) { 398df5c1128SToby Isaac PetscQuadrature fAllPoints; 39954ee0cceSMatthew G. Knepley PetscInt qdim, fNumPoints, q; 400df5c1128SToby Isaac const PetscReal *fPoints; 401df5c1128SToby Isaac 402df5c1128SToby Isaac ierr = PetscDualSpaceGetAllPoints(sp[f],&fAllPoints);CHKERRQ(ierr); 40354ee0cceSMatthew G. Knepley ierr = PetscQuadratureGetData(fAllPoints, &qdim, NULL, &fNumPoints, &fPoints, NULL);CHKERRQ(ierr); 40454ee0cceSMatthew 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); 405df5c1128SToby Isaac for (q = 0; q < fNumPoints*dim; ++q) points[numPoints*dim+q] = fPoints[q]; 406df5c1128SToby Isaac numPoints += fNumPoints; 407df5c1128SToby Isaac } 408df5c1128SToby Isaac } 409df5c1128SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 410df5c1128SToby Isaac ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 411df5c1128SToby Isaac PetscFunctionReturn(0); 412df5c1128SToby Isaac } 413df5c1128SToby Isaac 414*a6e0b375SMatthew G. Knepley static PetscErrorCode DMGetFirstLabelEntry_Private(DM dm, DM odm, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt height, PetscInt *lStart, PetscDS *ds) 415e5e52638SMatthew G. Knepley { 416*a6e0b375SMatthew G. Knepley DM plex; 417*a6e0b375SMatthew G. Knepley DMEnclosureType enc; 418e5e52638SMatthew G. Knepley DMLabel depthLabel; 419e5e52638SMatthew G. Knepley PetscInt dim, cdepth, ls = -1, i; 420e5e52638SMatthew G. Knepley PetscErrorCode ierr; 421e5e52638SMatthew G. Knepley 422e5e52638SMatthew G. Knepley PetscFunctionBegin; 423e5e52638SMatthew G. Knepley if (lStart) *lStart = -1; 424e5e52638SMatthew G. Knepley if (!label) PetscFunctionReturn(0); 425*a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dm, odm, &enc);CHKERRQ(ierr); 426e5e52638SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 427*a6e0b375SMatthew G. Knepley ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 428*a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr); 429e5e52638SMatthew G. Knepley cdepth = dim - height; 430e5e52638SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 431e5e52638SMatthew G. Knepley IS pointIS; 432e5e52638SMatthew G. Knepley const PetscInt *points; 433*a6e0b375SMatthew G. Knepley PetscInt pdepth, point; 434e5e52638SMatthew G. Knepley 435e5e52638SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 436e5e52638SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 437e5e52638SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 438*a6e0b375SMatthew G. Knepley ierr = DMGetEnclosurePoint(dm, odm, enc, points[0], &point);CHKERRQ(ierr); 439*a6e0b375SMatthew G. Knepley ierr = DMLabelGetValue(depthLabel, point, &pdepth);CHKERRQ(ierr); 440e5e52638SMatthew G. Knepley if (pdepth == cdepth) { 441*a6e0b375SMatthew G. Knepley ls = point; 442*a6e0b375SMatthew G. Knepley if (ds) {ierr = DMGetCellDS(dm, ls, ds);CHKERRQ(ierr);} 443e5e52638SMatthew G. Knepley } 444e5e52638SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 445e5e52638SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 446e5e52638SMatthew G. Knepley if (ls >= 0) break; 447e5e52638SMatthew G. Knepley } 448e5e52638SMatthew G. Knepley if (lStart) *lStart = ls; 449*a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 450e5e52638SMatthew G. Knepley PetscFunctionReturn(0); 451e5e52638SMatthew G. Knepley } 452e5e52638SMatthew G. Knepley 4530de51b56SMatthew G. Knepley /* 4540de51b56SMatthew G. Knepley This function iterates over a manifold, and interpolates the input function/field using the basis provided by the DS in our DM 4550de51b56SMatthew G. Knepley 4560de51b56SMatthew G. Knepley There are several different scenarios: 4570de51b56SMatthew G. Knepley 4580de51b56SMatthew G. Knepley 1) Volumetric mesh with volumetric auxiliary data 4590de51b56SMatthew G. Knepley 4600de51b56SMatthew G. Knepley Here minHeight=0 since we loop over cells. 4610de51b56SMatthew G. Knepley 4620de51b56SMatthew G. Knepley 2) Boundary mesh with boundary auxiliary data 4630de51b56SMatthew G. Knepley 4640de51b56SMatthew 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. 4650de51b56SMatthew G. Knepley 4660de51b56SMatthew G. Knepley 3) Volumetric mesh with boundary auxiliary data 4670de51b56SMatthew G. Knepley 4680de51b56SMatthew 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. 4690de51b56SMatthew G. Knepley 4700de51b56SMatthew G. Knepley The maxHeight is used to support enforcement of constraints in DMForest. 4710de51b56SMatthew G. Knepley 4720de51b56SMatthew G. Knepley If localU is given and not equal to localX, we call DMPlexInsertBoundaryValues() to complete it. 4730de51b56SMatthew G. Knepley 4740de51b56SMatthew 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. 4750de51b56SMatthew 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 4760de51b56SMatthew G. Knepley is for the volumetric mesh, but we are iterating over a surface, then the effective height is nonzero. When th effective height is nonzero, we need to extract 4770de51b56SMatthew G. Knepley dual spaces for the boundary from our input spaces. 4780de51b56SMatthew G. Knepley - After extracting all quadrature points, we tabulate the input fields and auxiliary fields on them. 4790de51b56SMatthew G. Knepley 4800de51b56SMatthew G. Knepley We check that the #dof(closure(p)) == #dual basis functionals(p) for a representative p in the iteration 4810de51b56SMatthew G. Knepley 4820de51b56SMatthew 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. 4830de51b56SMatthew G. Knepley */ 48446fa42a0SMatthew G. Knepley static PetscErrorCode DMProjectLocal_Generic_Plex(DM dm, PetscReal time, Vec localU, 4851c531cf8SMatthew G. Knepley PetscInt Ncc, const PetscInt comps[], DMLabel label, PetscInt numIds, const PetscInt ids[], 4868c6c5593SMatthew G. Knepley DMBoundaryConditionType type, void (**funcs)(void), void **ctxs, 4878c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 4888c6c5593SMatthew G. Knepley { 489*a6e0b375SMatthew G. Knepley DM plex, dmIn, plexIn, dmAux = NULL, plexAux = NULL, tdm; 490*a6e0b375SMatthew G. Knepley DMEnclosureType encIn, encAux; 491*a6e0b375SMatthew G. Knepley PetscDS ds = NULL, dsIn = NULL, dsAux = NULL; 492ca3d3a14SMatthew G. Knepley Vec localA = NULL, tv; 49347923291SMatthew G. Knepley PetscSection section; 4948c6c5593SMatthew G. Knepley PetscDualSpace *sp, *cellsp; 4950c898c18SMatthew G. Knepley PetscReal **basisTab = NULL, **basisDerTab = NULL, **basisTabAux = NULL, **basisDerTabAux = NULL; 4968c6c5593SMatthew G. Knepley PetscInt *Nc; 497*a6e0b375SMatthew G. Knepley PetscInt dim, dimEmbed, depth, minHeight, maxHeight, h, Nf, NfIn, NfAux = 0, f; 498ca3d3a14SMatthew G. Knepley PetscBool *isFE, hasFE = PETSC_FALSE, hasFV = PETSC_FALSE, auxBd = PETSC_FALSE, transform; 499c330f8ffSToby Isaac DMField coordField; 500c330f8ffSToby Isaac DMLabel depthLabel; 501c330f8ffSToby Isaac PetscQuadrature allPoints = NULL; 50247923291SMatthew G. Knepley PetscErrorCode ierr; 50347923291SMatthew G. Knepley 50447923291SMatthew G. Knepley PetscFunctionBegin; 505*a6e0b375SMatthew G. Knepley if (localU) {ierr = VecGetDM(localU, &dmIn);CHKERRQ(ierr);} 506*a6e0b375SMatthew G. Knepley else {dmIn = dm;} 50788aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 50888aca1feSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &localA);CHKERRQ(ierr); 509*a6e0b375SMatthew G. Knepley ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 510*a6e0b375SMatthew G. Knepley ierr = DMConvert(dmIn, DMPLEX, &plexIn);CHKERRQ(ierr); 511*a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dmIn, dm, &encIn);CHKERRQ(ierr); 512*a6e0b375SMatthew G. Knepley ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr); 5138c6c5593SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 514*a6e0b375SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(plex, &minHeight);CHKERRQ(ierr); 515ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformDM_Internal(dm, &tdm);CHKERRQ(ierr); 516ca3d3a14SMatthew G. Knepley ierr = DMGetBasisTransformVec_Internal(dm, &tv);CHKERRQ(ierr); 517ca3d3a14SMatthew G. Knepley ierr = DMHasBasisTransform(dm, &transform);CHKERRQ(ierr); 5180de51b56SMatthew G. Knepley /* Auxiliary information can only be used with interpolation of field functions */ 519*a6e0b375SMatthew G. Knepley if (dmAux) { 520*a6e0b375SMatthew G. Knepley ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr); 521*a6e0b375SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_ESSENTIAL_BD_FIELD || type == DM_BC_NATURAL_FIELD) { 522*a6e0b375SMatthew G. Knepley if (!localA) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Missing localA vector"); 523*a6e0b375SMatthew G. Knepley if (!minHeight) { 52488aca1feSMatthew G. Knepley DMLabel spmap; 52588aca1feSMatthew G. Knepley 52688aca1feSMatthew G. Knepley /* If dmAux is a surface, then force the projection to take place over a surface */ 527*a6e0b375SMatthew G. Knepley ierr = DMPlexGetSubpointMap(plexAux, &spmap);CHKERRQ(ierr); 528e39fcb4eSStefano Zampini if (spmap) { 529*a6e0b375SMatthew G. Knepley ierr = DMPlexGetVTKCellHeight(plexAux, &minHeight);CHKERRQ(ierr); 530e39fcb4eSStefano Zampini auxBd = minHeight ? PETSC_TRUE : PETSC_FALSE; 531e39fcb4eSStefano Zampini } 53288aca1feSMatthew G. Knepley } 5330de51b56SMatthew G. Knepley } 534*a6e0b375SMatthew G. Knepley } 535*a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr); 536*a6e0b375SMatthew G. Knepley ierr = DMPlexGetDepthLabel(plex, &depthLabel);CHKERRQ(ierr); 537*a6e0b375SMatthew G. Knepley ierr = DMPlexGetMaxProjectionHeight(plex, &maxHeight);CHKERRQ(ierr); 5382af58022SMatthew G. Knepley maxHeight = PetscMax(maxHeight, minHeight); 539e39fcb4eSStefano Zampini if (maxHeight < 0 || maxHeight > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum projection height %D not in [0, %D)", maxHeight, dim); 540*a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, 0, NULL, &ds);CHKERRQ(ierr); 541*a6e0b375SMatthew G. Knepley if (!ds) {ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);} 542*a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dmIn, dm, label, numIds, ids, 0, NULL, &dsIn);CHKERRQ(ierr); 543*a6e0b375SMatthew G. Knepley if (!dsIn) {ierr = DMGetDS(dmIn, &dsIn);CHKERRQ(ierr);} 544*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 545*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsIn, &NfIn);CHKERRQ(ierr); 5468c6c5593SMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 54792fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 5480c898c18SMatthew G. Knepley if (dmAux) { 549*a6e0b375SMatthew G. Knepley ierr = DMGetDS(dmAux, &dsAux);CHKERRQ(ierr); 550*a6e0b375SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 5510c898c18SMatthew G. Knepley } 552*a6e0b375SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 553496733e2SMatthew G. Knepley ierr = PetscMalloc2(Nf, &isFE, Nf, &sp);CHKERRQ(ierr); 5548c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscMalloc1(Nf, &cellsp);CHKERRQ(ierr);} 5558c6c5593SMatthew G. Knepley else {cellsp = sp;} 556*a6e0b375SMatthew G. Knepley if (localU && localU != localX) {ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, localU, time, NULL, NULL, NULL);CHKERRQ(ierr);} 5578c6c5593SMatthew G. Knepley /* Get cell dual spaces */ 55847923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 55947923291SMatthew G. Knepley PetscObject obj; 56047923291SMatthew G. Knepley PetscClassId id; 56147923291SMatthew G. Knepley 562*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, f, &obj);CHKERRQ(ierr); 56347923291SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 56447923291SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 56547923291SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 56647923291SMatthew G. Knepley 56747923291SMatthew G. Knepley hasFE = PETSC_TRUE; 56847923291SMatthew G. Knepley isFE[f] = PETSC_TRUE; 56947923291SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &cellsp[f]);CHKERRQ(ierr); 57047923291SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 57147923291SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 57247923291SMatthew G. Knepley 57347923291SMatthew G. Knepley hasFV = PETSC_TRUE; 57447923291SMatthew G. Knepley isFE[f] = PETSC_FALSE; 57547923291SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &cellsp[f]);CHKERRQ(ierr); 57647923291SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 57747923291SMatthew G. Knepley } 578c330f8ffSToby Isaac ierr = DMGetCoordinateField(dm,&coordField);CHKERRQ(ierr); 5790c898c18SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 58074fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 58174fc349aSMatthew G. Knepley PetscFE fem, subfem; 582*a6e0b375SMatthew G. Knepley PetscBool isfe; 5834a3e9fdbSToby Isaac const PetscReal *points; 5844a3e9fdbSToby Isaac PetscInt numPoints; 5850c898c18SMatthew G. Knepley 5862af58022SMatthew G. Knepley if (maxHeight > minHeight) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Field projection not supported for face interpolation"); 58754ee0cceSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 58854ee0cceSMatthew G. Knepley if (!effectiveHeight) {sp[f] = cellsp[f];} 58954ee0cceSMatthew G. Knepley else {ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr);} 59054ee0cceSMatthew G. Knepley } 59154ee0cceSMatthew G. Knepley ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-effectiveHeight,funcs,&allPoints);CHKERRQ(ierr); 592c330f8ffSToby Isaac ierr = PetscQuadratureGetData(allPoints,NULL,NULL,&numPoints,&points,NULL);CHKERRQ(ierr); 593*a6e0b375SMatthew G. Knepley ierr = PetscMalloc4(NfIn, &basisTab, NfIn, &basisDerTab, NfAux, &basisTabAux, NfAux, &basisDerTabAux);CHKERRQ(ierr); 594*a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 595*a6e0b375SMatthew G. Knepley ierr = PetscDSIsFE_Internal(dsIn, f, &isfe);CHKERRQ(ierr); 596*a6e0b375SMatthew G. Knepley if (!isfe) continue; 597*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr); 59874fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 59974fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 60074fc349aSMatthew G. Knepley ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 6010c898c18SMatthew G. Knepley } 6020c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 603*a6e0b375SMatthew G. Knepley ierr = PetscDSIsFE_Internal(dsAux, f, &isfe);CHKERRQ(ierr); 604*a6e0b375SMatthew G. Knepley if (!isfe) continue; 605*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 60674fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 60774fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 60874fc349aSMatthew G. Knepley ierr = PetscFEGetTabulation(subfem, numPoints, points, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 6090c898c18SMatthew G. Knepley } 6100c898c18SMatthew G. Knepley } 61147923291SMatthew G. Knepley /* Note: We make no attempt to optimize for height. Higher height things just overwrite the lower height results. */ 6122af58022SMatthew G. Knepley for (h = minHeight; h <= maxHeight; h++) { 61388aca1feSMatthew G. Knepley PetscInt effectiveHeight = h - (auxBd ? 0 : minHeight); 614*a6e0b375SMatthew G. Knepley PetscDS dsEff = ds; 6158c6c5593SMatthew G. Knepley PetscScalar *values; 616b7260050SToby Isaac PetscBool *fieldActive; 617b7260050SToby Isaac PetscInt maxDegree; 618e5e52638SMatthew G. Knepley PetscInt pStart, pEnd, p, lStart, spDim, totDim, numValues; 619c330f8ffSToby Isaac IS heightIS; 6208c6c5593SMatthew G. Knepley 621*a6e0b375SMatthew G. Knepley /* Note we assume that dm and dmIn share the same topology */ 622*a6e0b375SMatthew G. Knepley ierr = DMPlexGetHeightStratum(plex, h, &pStart, &pEnd);CHKERRQ(ierr); 623*a6e0b375SMatthew G. Knepley ierr = DMGetFirstLabelEntry_Private(dm, dm, label, numIds, ids, h, &lStart, NULL);CHKERRQ(ierr); 624c330f8ffSToby Isaac ierr = DMLabelGetStratumIS(depthLabel, depth - h, &heightIS);CHKERRQ(ierr); 6258c6c5593SMatthew G. Knepley if (!h) { 6268c6c5593SMatthew G. Knepley PetscInt cEndInterior; 6278c6c5593SMatthew G. Knepley 628485ad865SMatthew G. Knepley ierr = DMPlexGetInteriorCellStratum(dm, NULL, &cEndInterior);CHKERRQ(ierr); 6298c6c5593SMatthew G. Knepley pEnd = cEndInterior < 0 ? pEnd : cEndInterior; 6308c6c5593SMatthew G. Knepley } 631c330f8ffSToby Isaac if (pEnd <= pStart) { 632c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 633c330f8ffSToby Isaac continue; 634c330f8ffSToby Isaac } 63547923291SMatthew G. Knepley /* Compute totDim, the number of dofs in the closure of a point at this height */ 63647923291SMatthew G. Knepley totDim = 0; 63747923291SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 638bec1f4aaSMatthew G. Knepley if (!effectiveHeight) { 63947923291SMatthew G. Knepley sp[f] = cellsp[f]; 64047923291SMatthew G. Knepley } else { 641bec1f4aaSMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(cellsp[f], effectiveHeight, &sp[f]);CHKERRQ(ierr); 64247923291SMatthew G. Knepley if (!sp[f]) continue; 64347923291SMatthew G. Knepley } 64447923291SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 6459c3cf19fSMatthew G. Knepley totDim += spDim; 64647923291SMatthew G. Knepley } 647*a6e0b375SMatthew G. Knepley ierr = DMPlexVecGetClosure(plex, section, localX, lStart < 0 ? pStart : lStart, &numValues, NULL);CHKERRQ(ierr); 6488c6c5593SMatthew G. Knepley if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The section point closure size %d != dual space dimension %d", numValues, totDim); 649c330f8ffSToby Isaac if (!totDim) { 650c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 651c330f8ffSToby Isaac continue; 652c330f8ffSToby Isaac } 653*a6e0b375SMatthew G. Knepley if (effectiveHeight) {ierr = PetscDSGetHeightSubspace(ds, effectiveHeight, &dsEff);CHKERRQ(ierr);} 65447923291SMatthew G. Knepley /* Loop over points at this height */ 65569291d52SBarry Smith ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 65669291d52SBarry Smith ierr = DMGetWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 6578c6c5593SMatthew G. Knepley for (f = 0; f < Nf; ++f) fieldActive[f] = (funcs[f] && sp[f]) ? PETSC_TRUE : PETSC_FALSE; 6588c6c5593SMatthew G. Knepley if (label) { 6598c6c5593SMatthew G. Knepley PetscInt i; 66047923291SMatthew G. Knepley 66147923291SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 662c330f8ffSToby Isaac IS pointIS, isectIS; 66347923291SMatthew G. Knepley const PetscInt *points; 6648c6c5593SMatthew G. Knepley PetscInt n; 665c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 666c330f8ffSToby Isaac PetscQuadrature quad = NULL; 66747923291SMatthew G. Knepley 66847923291SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 66947923291SMatthew G. Knepley if (!pointIS) continue; /* No points with that id on this process */ 670c330f8ffSToby Isaac ierr = ISIntersect(pointIS,heightIS,&isectIS);CHKERRQ(ierr); 671c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 672c330f8ffSToby Isaac if (!isectIS) continue; 673c330f8ffSToby Isaac ierr = ISGetLocalSize(isectIS, &n);CHKERRQ(ierr); 674c330f8ffSToby Isaac ierr = ISGetIndices(isectIS, &points);CHKERRQ(ierr); 675b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,isectIS,NULL,&maxDegree);CHKERRQ(ierr); 676b7260050SToby Isaac if (maxDegree <= 1) { 677c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,isectIS,&quad);CHKERRQ(ierr); 678c330f8ffSToby Isaac } 679c330f8ffSToby Isaac if (!quad) { 680c330f8ffSToby Isaac if (!h && allPoints) { 681c330f8ffSToby Isaac quad = allPoints; 682c330f8ffSToby Isaac allPoints = NULL; 683c330f8ffSToby Isaac } else { 684c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 685c330f8ffSToby Isaac } 686c330f8ffSToby Isaac } 687*a6e0b375SMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField,isectIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr); 68847923291SMatthew G. Knepley for (p = 0; p < n; ++p) { 68947923291SMatthew G. Knepley const PetscInt point = points[p]; 69047923291SMatthew G. Knepley 691580bdb30SBarry Smith ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr); 692c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 693*a6e0b375SMatthew G. Knepley ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, point, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 69447923291SMatthew G. Knepley if (ierr) { 69547923291SMatthew G. Knepley PetscErrorCode ierr2; 69669291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 69769291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 69847923291SMatthew G. Knepley CHKERRQ(ierr); 69947923291SMatthew G. Knepley } 700*a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, point, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 701*a6e0b375SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, point, Ncc, comps, values, mode);CHKERRQ(ierr); 70247923291SMatthew G. Knepley } 703c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p,p+1,&chunkgeom);CHKERRQ(ierr); 704c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 705c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 706c330f8ffSToby Isaac ierr = ISRestoreIndices(isectIS, &points);CHKERRQ(ierr); 707c330f8ffSToby Isaac ierr = ISDestroy(&isectIS);CHKERRQ(ierr); 70847923291SMatthew G. Knepley } 7098c6c5593SMatthew G. Knepley } else { 710c330f8ffSToby Isaac PetscFEGeom *fegeom = NULL, *chunkgeom = NULL; 711c330f8ffSToby Isaac PetscQuadrature quad = NULL; 712c330f8ffSToby Isaac IS pointIS; 713c330f8ffSToby Isaac 714c330f8ffSToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF,pEnd-pStart,pStart,1,&pointIS);CHKERRQ(ierr); 715b7260050SToby Isaac ierr = DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);CHKERRQ(ierr); 716b7260050SToby Isaac if (maxDegree <= 1) { 717c330f8ffSToby Isaac ierr = DMFieldCreateDefaultQuadrature(coordField,pointIS,&quad);CHKERRQ(ierr); 718c330f8ffSToby Isaac } 719c330f8ffSToby Isaac if (!quad) { 720c330f8ffSToby Isaac if (!h && allPoints) { 721c330f8ffSToby Isaac quad = allPoints; 722c330f8ffSToby Isaac allPoints = NULL; 723c330f8ffSToby Isaac } else { 724c330f8ffSToby Isaac ierr = PetscDualSpaceGetAllPointsUnion(Nf,sp,dim-h,funcs,&quad);CHKERRQ(ierr); 725c330f8ffSToby Isaac } 726c330f8ffSToby Isaac } 727*a6e0b375SMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,(effectiveHeight && h == minHeight)?PETSC_TRUE:PETSC_FALSE,&fegeom);CHKERRQ(ierr); 7288c6c5593SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 729580bdb30SBarry Smith ierr = PetscArrayzero(values, numValues);CHKERRQ(ierr); 730c330f8ffSToby Isaac ierr = PetscFEGeomGetChunk(fegeom,p-pStart,p-pStart+1,&chunkgeom);CHKERRQ(ierr); 731*a6e0b375SMatthew G. Knepley ierr = DMProjectPoint_Private(dm, dsEff, plexIn, encIn, dsIn, plexAux, encAux, dsAux, chunkgeom, effectiveHeight, time, localU, localA, hasFE, hasFV, isFE, sp, p, basisTab, basisDerTab, basisTabAux, basisDerTabAux, type, funcs, ctxs, fieldActive, values); 7328c6c5593SMatthew G. Knepley if (ierr) { 7338c6c5593SMatthew G. Knepley PetscErrorCode ierr2; 73469291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr2); 73569291d52SBarry Smith ierr2 = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr2); 7368c6c5593SMatthew G. Knepley CHKERRQ(ierr); 7378c6c5593SMatthew G. Knepley } 738*a6e0b375SMatthew G. Knepley if (transform) {ierr = DMPlexBasisTransformPoint_Internal(plex, tdm, tv, p, fieldActive, PETSC_FALSE, values);CHKERRQ(ierr);} 739*a6e0b375SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(plex, section, localX, fieldActive, p, Ncc, comps, values, mode);CHKERRQ(ierr); 7408c6c5593SMatthew G. Knepley } 741c330f8ffSToby Isaac ierr = PetscFEGeomRestoreChunk(fegeom,p-pStart,pStart-p+1,&chunkgeom);CHKERRQ(ierr); 742c330f8ffSToby Isaac ierr = PetscFEGeomDestroy(&fegeom);CHKERRQ(ierr); 743c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr); 744c330f8ffSToby Isaac ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 7458c6c5593SMatthew G. Knepley } 746c330f8ffSToby Isaac ierr = ISDestroy(&heightIS);CHKERRQ(ierr); 74769291d52SBarry Smith ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr); 74869291d52SBarry Smith ierr = DMRestoreWorkArray(dm, Nf, MPI_INT, &fieldActive);CHKERRQ(ierr); 74947923291SMatthew G. Knepley } 7508c6c5593SMatthew G. Knepley /* Cleanup */ 7510c898c18SMatthew G. Knepley if (type == DM_BC_ESSENTIAL_FIELD || type == DM_BC_NATURAL_FIELD) { 75274fc349aSMatthew G. Knepley PetscInt effectiveHeight = auxBd ? minHeight : 0; 75374fc349aSMatthew G. Knepley PetscFE fem, subfem; 754*a6e0b375SMatthew G. Knepley PetscBool isfe; 7550c898c18SMatthew G. Knepley 756*a6e0b375SMatthew G. Knepley for (f = 0; f < NfIn; ++f) { 757*a6e0b375SMatthew G. Knepley ierr = PetscDSIsFE_Internal(dsIn, f, &isfe);CHKERRQ(ierr); 758*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsIn, f, (PetscObject *) &fem);CHKERRQ(ierr); 75974fc349aSMatthew G. Knepley if (!effectiveHeight) {subfem = fem;} 76074fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 76174fc349aSMatthew G. Knepley ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTab[f], &basisDerTab[f], NULL);CHKERRQ(ierr); 7620c898c18SMatthew G. Knepley } 7630c898c18SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 764*a6e0b375SMatthew G. Knepley ierr = PetscDSIsFE_Internal(dsAux, f, &isfe);CHKERRQ(ierr); 765*a6e0b375SMatthew G. Knepley ierr = PetscDSGetDiscretization(dsAux, f, (PetscObject *) &fem);CHKERRQ(ierr); 76674fc349aSMatthew G. Knepley if (!effectiveHeight || auxBd) {subfem = fem;} 76774fc349aSMatthew G. Knepley else {ierr = PetscFEGetHeightSubspace(fem, effectiveHeight, &subfem);CHKERRQ(ierr);} 76874fc349aSMatthew G. Knepley ierr = PetscFERestoreTabulation(subfem, 0, NULL, &basisTabAux[f], &basisDerTabAux[f], NULL);CHKERRQ(ierr); 7690c898c18SMatthew G. Knepley } 7700c898c18SMatthew G. Knepley ierr = PetscFree4(basisTab, basisDerTab, basisTabAux, basisDerTabAux);CHKERRQ(ierr); 7710c898c18SMatthew G. Knepley } 772c330f8ffSToby Isaac ierr = PetscQuadratureDestroy(&allPoints);CHKERRQ(ierr); 773496733e2SMatthew G. Knepley ierr = PetscFree2(isFE, sp);CHKERRQ(ierr); 7748c6c5593SMatthew G. Knepley if (maxHeight > 0) {ierr = PetscFree(cellsp);CHKERRQ(ierr);} 775*a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plex);CHKERRQ(ierr); 776*a6e0b375SMatthew G. Knepley ierr = DMDestroy(&plexIn);CHKERRQ(ierr); 777*a6e0b375SMatthew G. Knepley if (dmAux) {ierr = DMDestroy(&plexAux);CHKERRQ(ierr);} 7788c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 77947923291SMatthew G. Knepley } 7808c6c5593SMatthew G. Knepley 7818c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFunctionLocal_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 7828c6c5593SMatthew G. Knepley { 7838c6c5593SMatthew G. Knepley PetscErrorCode ierr; 7848c6c5593SMatthew G. Knepley 7858c6c5593SMatthew G. Knepley PetscFunctionBegin; 7861c531cf8SMatthew 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); 7878c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 7888c6c5593SMatthew G. Knepley } 7898c6c5593SMatthew G. Knepley 7901c531cf8SMatthew 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) 7918c6c5593SMatthew G. Knepley { 7928c6c5593SMatthew G. Knepley PetscErrorCode ierr; 7938c6c5593SMatthew G. Knepley 7948c6c5593SMatthew G. Knepley PetscFunctionBegin; 7951c531cf8SMatthew 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); 79647923291SMatthew G. Knepley PetscFunctionReturn(0); 79747923291SMatthew G. Knepley } 79847923291SMatthew G. Knepley 7998c6c5593SMatthew G. Knepley PetscErrorCode DMProjectFieldLocal_Plex(DM dm, PetscReal time, Vec localU, 80047923291SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 80147923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 80247923291SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 803191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 80447923291SMatthew G. Knepley InsertMode mode, Vec localX) 80547923291SMatthew G. Knepley { 80647923291SMatthew G. Knepley PetscErrorCode ierr; 80747923291SMatthew G. Knepley 80847923291SMatthew G. Knepley PetscFunctionBegin; 8091c531cf8SMatthew 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); 8108c6c5593SMatthew G. Knepley PetscFunctionReturn(0); 81147923291SMatthew G. Knepley } 81247923291SMatthew G. Knepley 8131c531cf8SMatthew G. Knepley PetscErrorCode DMProjectFieldLabelLocal_Plex(DM dm, PetscReal time, DMLabel label, PetscInt numIds, const PetscInt ids[], PetscInt Ncc, const PetscInt comps[], Vec localU, 8148c6c5593SMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 8158c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 8168c6c5593SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 817191494d9SMatthew G. Knepley PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 8188c6c5593SMatthew G. Knepley InsertMode mode, Vec localX) 8198c6c5593SMatthew G. Knepley { 8208c6c5593SMatthew G. Knepley PetscErrorCode ierr; 82147923291SMatthew G. Knepley 8228c6c5593SMatthew G. Knepley PetscFunctionBegin; 8231c531cf8SMatthew 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); 82447923291SMatthew G. Knepley PetscFunctionReturn(0); 82547923291SMatthew G. Knepley } 826