16dbbd306SMatthew G. Knepley #include <petscdmplex.h> /*I "petscdmplex.h" I*/ 2254c1ad2SMatthew G. Knepley #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 3*924a1b8fSMatthew G. Knepley #include <petscds.h> 46dbbd306SMatthew G. Knepley #include <petscfv.h> 56dbbd306SMatthew G. Knepley 61475bf87SMatthew G. Knepley PETSC_STATIC_INLINE void WaxpyD(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) {PetscInt d; for (d = 0; d < dim; ++d) w[d] = a*x[d] + y[d];} 71475bf87SMatthew G. Knepley PETSC_STATIC_INLINE PetscReal DotD(PetscInt dim, const PetscScalar *x, const PetscReal *y) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d])*y[d]; return sum;} 81475bf87SMatthew G. Knepley PETSC_STATIC_INLINE PetscReal DotRealD(PetscInt dim, const PetscReal *x, const PetscReal *y) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += x[d]*y[d]; return sum;} 91475bf87SMatthew G. Knepley PETSC_STATIC_INLINE PetscReal NormD(PetscInt dim, const PetscReal *x) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += x[d]*x[d]; return PetscSqrtReal(sum);} 10254c1ad2SMatthew G. Knepley 11254c1ad2SMatthew G. Knepley #undef __FUNCT__ 12a0ac79e7SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometry" 13c510411aSMatthew G. Knepley /*@ 14c510411aSMatthew G. Knepley DMPlexTSGetGeometry - Return precomputed geometric data 15c510411aSMatthew G. Knepley 16c510411aSMatthew G. Knepley Input Parameter: 17c510411aSMatthew G. Knepley . dm - The DM 18c510411aSMatthew G. Knepley 19c510411aSMatthew G. Knepley Output Parameters: 20c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry 21c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry 22c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell 23c510411aSMatthew G. Knepley 24c510411aSMatthew G. Knepley Level: developer 25c510411aSMatthew G. Knepley 26c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal() 27c510411aSMatthew G. Knepley @*/ 28a0ac79e7SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometry(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius) 29a0ac79e7SMatthew G. Knepley { 30a0ac79e7SMatthew G. Knepley DMTS dmts; 31a0ac79e7SMatthew G. Knepley PetscErrorCode ierr; 32a0ac79e7SMatthew G. Knepley 33a0ac79e7SMatthew G. Knepley PetscFunctionBegin; 34*924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 35a0ac79e7SMatthew G. Knepley ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr); 36*924a1b8fSMatthew G. Knepley if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom", (PetscObject *) facegeom);CHKERRQ(ierr);} 37*924a1b8fSMatthew G. Knepley if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom", (PetscObject *) cellgeom);CHKERRQ(ierr);} 38*924a1b8fSMatthew G. Knepley if (minRadius) {ierr = DMTSGetMinRadius(dm, minRadius);CHKERRQ(ierr);} 39*924a1b8fSMatthew G. Knepley PetscFunctionReturn(0); 40*924a1b8fSMatthew G. Knepley } 41*924a1b8fSMatthew G. Knepley 42*924a1b8fSMatthew G. Knepley #undef __FUNCT__ 43*924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM" 44*924a1b8fSMatthew G. Knepley /*@ 45*924a1b8fSMatthew G. Knepley DMPlexTSGetGradientDM - Return gradient data layout 46*924a1b8fSMatthew G. Knepley 47*924a1b8fSMatthew G. Knepley Input Parameter: 48*924a1b8fSMatthew G. Knepley . dm - The DM 49*924a1b8fSMatthew G. Knepley 50*924a1b8fSMatthew G. Knepley Output Parameter: 51*924a1b8fSMatthew G. Knepley . dmGrad - The layout for gradient values 52*924a1b8fSMatthew G. Knepley 53*924a1b8fSMatthew G. Knepley Level: developer 54*924a1b8fSMatthew G. Knepley 55*924a1b8fSMatthew G. Knepley .seealso: DMPlexTSGetGeometry(), DMPlexTSSetRHSFunctionLocal() 56*924a1b8fSMatthew G. Knepley @*/ 57*924a1b8fSMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, DM *dmGrad) 58*924a1b8fSMatthew G. Knepley { 59*924a1b8fSMatthew G. Knepley DMTS dmts; 60*924a1b8fSMatthew G. Knepley PetscErrorCode ierr; 61*924a1b8fSMatthew G. Knepley 62*924a1b8fSMatthew G. Knepley PetscFunctionBegin; 63*924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 64*924a1b8fSMatthew G. Knepley PetscValidPointer(dmGrad,2); 65*924a1b8fSMatthew G. Knepley ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr); 66*924a1b8fSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad", (PetscObject *) dmGrad);CHKERRQ(ierr); 67a0ac79e7SMatthew G. Knepley PetscFunctionReturn(0); 68a0ac79e7SMatthew G. Knepley } 69a0ac79e7SMatthew G. Knepley 70a0ac79e7SMatthew G. Knepley #undef __FUNCT__ 71254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGeometry" 72*924a1b8fSMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGeometry(DM dm, PetscFV fvm, DMTS dmts) 73254c1ad2SMatthew G. Knepley { 74254c1ad2SMatthew G. Knepley DM dmFace, dmCell; 75254c1ad2SMatthew G. Knepley DMLabel ghostLabel; 76254c1ad2SMatthew G. Knepley PetscSection sectionFace, sectionCell; 77254c1ad2SMatthew G. Knepley PetscSection coordSection; 78*924a1b8fSMatthew G. Knepley Vec coordinates, cellgeom, facegeom; 79254c1ad2SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 80*924a1b8fSMatthew G. Knepley PetscReal minradius, gminradius; 81254c1ad2SMatthew G. Knepley PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 82254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 83254c1ad2SMatthew G. Knepley 84254c1ad2SMatthew G. Knepley PetscFunctionBegin; 85c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 86254c1ad2SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 87254c1ad2SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 88254c1ad2SMatthew G. Knepley /* Make cell centroids and volumes */ 89254c1ad2SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 9046e270d4SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, PETSC_DETERMINE, coordSection);CHKERRQ(ierr); 91254c1ad2SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 92254c1ad2SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 93254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 94254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 95254c1ad2SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 96254c1ad2SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, sizeof(CellGeom)/sizeof(PetscScalar));CHKERRQ(ierr);} 97254c1ad2SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 98254c1ad2SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 99254c1ad2SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 100*924a1b8fSMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, &cellgeom);CHKERRQ(ierr); 101*924a1b8fSMatthew G. Knepley ierr = VecGetArray(cellgeom, &cgeom);CHKERRQ(ierr); 102254c1ad2SMatthew G. Knepley for (c = cStart; c < cEndInterior; ++c) { 103254c1ad2SMatthew G. Knepley CellGeom *cg; 104254c1ad2SMatthew G. Knepley 105254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 106254c1ad2SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 107254c1ad2SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 108254c1ad2SMatthew G. Knepley } 109254c1ad2SMatthew G. Knepley /* Compute face normals and minimum cell radius */ 110254c1ad2SMatthew G. Knepley ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 111254c1ad2SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 112254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 113254c1ad2SMatthew G. Knepley ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 114254c1ad2SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, sizeof(FaceGeom)/sizeof(PetscScalar));CHKERRQ(ierr);} 115254c1ad2SMatthew G. Knepley ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 116254c1ad2SMatthew G. Knepley ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 117254c1ad2SMatthew G. Knepley ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 118*924a1b8fSMatthew G. Knepley ierr = DMCreateLocalVector(dmFace, &facegeom);CHKERRQ(ierr); 119*924a1b8fSMatthew G. Knepley ierr = VecGetArray(facegeom, &fgeom);CHKERRQ(ierr); 120254c1ad2SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 121254c1ad2SMatthew G. Knepley minradius = PETSC_MAX_REAL; 122254c1ad2SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 123254c1ad2SMatthew G. Knepley FaceGeom *fg; 124254c1ad2SMatthew G. Knepley PetscReal area; 125254c1ad2SMatthew G. Knepley PetscInt ghost, d; 126254c1ad2SMatthew G. Knepley 127254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr); 128254c1ad2SMatthew G. Knepley if (ghost >= 0) continue; 129254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 130254c1ad2SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 131254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] *= area; 132254c1ad2SMatthew G. Knepley /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 133254c1ad2SMatthew G. Knepley { 134254c1ad2SMatthew G. Knepley CellGeom *cL, *cR; 135254c1ad2SMatthew G. Knepley const PetscInt *cells; 136254c1ad2SMatthew G. Knepley PetscReal *lcentroid, *rcentroid; 1371475bf87SMatthew G. Knepley PetscReal v[3]; 138254c1ad2SMatthew G. Knepley 139254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 140254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 141254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 142254c1ad2SMatthew G. Knepley lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 143254c1ad2SMatthew G. Knepley rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 144254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, lcentroid, rcentroid, v); 1451475bf87SMatthew G. Knepley if (DotRealD(dim, fg->normal, v) < 0) { 146254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 147254c1ad2SMatthew G. Knepley } 1481475bf87SMatthew G. Knepley if (DotRealD(dim, fg->normal, v) <= 0) { 149809582abSMatthew G. Knepley if (dim == 2) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g) v (%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) v[0], (double) v[1]); 150809582abSMatthew G. Knepley if (dim == 3) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) fg->normal[2], (double) v[0], (double) v[1], (double) v[2]); 151254c1ad2SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 152254c1ad2SMatthew G. Knepley } 153254c1ad2SMatthew G. Knepley if (cells[0] < cEndInterior) { 154254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, fg->centroid, cL->centroid, v); 155254c1ad2SMatthew G. Knepley minradius = PetscMin(minradius, NormD(dim, v)); 156254c1ad2SMatthew G. Knepley } 157254c1ad2SMatthew G. Knepley if (cells[1] < cEndInterior) { 158254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, fg->centroid, cR->centroid, v); 159254c1ad2SMatthew G. Knepley minradius = PetscMin(minradius, NormD(dim, v)); 160254c1ad2SMatthew G. Knepley } 161254c1ad2SMatthew G. Knepley } 162254c1ad2SMatthew G. Knepley } 163*924a1b8fSMatthew G. Knepley ierr = MPI_Allreduce(&minradius, &gminradius, 1, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 164*924a1b8fSMatthew G. Knepley ierr = DMTSSetMinRadius(dm, gminradius);CHKERRQ(ierr); 165254c1ad2SMatthew G. Knepley /* Compute centroids of ghost cells */ 166254c1ad2SMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) { 167254c1ad2SMatthew G. Knepley FaceGeom *fg; 168254c1ad2SMatthew G. Knepley const PetscInt *cone, *support; 169254c1ad2SMatthew G. Knepley PetscInt coneSize, supportSize, s; 170254c1ad2SMatthew G. Knepley 171254c1ad2SMatthew G. Knepley ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 172254c1ad2SMatthew G. Knepley if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 173254c1ad2SMatthew G. Knepley ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 174254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 175254c1ad2SMatthew G. Knepley if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 1", cone[0], supportSize); 176254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 177254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 178254c1ad2SMatthew G. Knepley for (s = 0; s < 2; ++s) { 179254c1ad2SMatthew G. Knepley /* Reflect ghost centroid across plane of face */ 180254c1ad2SMatthew G. Knepley if (support[s] == c) { 181254c1ad2SMatthew G. Knepley const CellGeom *ci; 182254c1ad2SMatthew G. Knepley CellGeom *cg; 1831475bf87SMatthew G. Knepley PetscReal c2f[3], a; 184254c1ad2SMatthew G. Knepley 185254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 186254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 1871475bf87SMatthew G. Knepley a = DotRealD(dim, c2f, fg->normal)/DotRealD(dim, fg->normal, fg->normal); 188254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 189254c1ad2SMatthew G. Knepley WaxpyD(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 190254c1ad2SMatthew G. Knepley cg->volume = ci->volume; 191254c1ad2SMatthew G. Knepley } 192254c1ad2SMatthew G. Knepley } 193254c1ad2SMatthew G. Knepley } 194*924a1b8fSMatthew G. Knepley ierr = VecRestoreArray(facegeom, &fgeom);CHKERRQ(ierr); 195*924a1b8fSMatthew G. Knepley ierr = VecRestoreArray(cellgeom, &cgeom);CHKERRQ(ierr); 196*924a1b8fSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom", (PetscObject) facegeom);CHKERRQ(ierr); 197*924a1b8fSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom", (PetscObject) cellgeom);CHKERRQ(ierr); 198*924a1b8fSMatthew G. Knepley ierr = VecDestroy(&facegeom);CHKERRQ(ierr); 199*924a1b8fSMatthew G. Knepley ierr = VecDestroy(&cellgeom);CHKERRQ(ierr); 200254c1ad2SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 201254c1ad2SMatthew G. Knepley ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 202254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 203254c1ad2SMatthew G. Knepley } 204254c1ad2SMatthew G. Knepley 205254c1ad2SMatthew G. Knepley #undef __FUNCT__ 206c5148223SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction" 207c5148223SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 208254c1ad2SMatthew G. Knepley { 2098cd7fcbbSMatthew G. Knepley DMLabel ghostLabel; 210c5148223SMatthew G. Knepley PetscScalar *dx, *grad, **gref; 211c5148223SMatthew G. Knepley PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 212254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 213254c1ad2SMatthew G. Knepley 214254c1ad2SMatthew G. Knepley PetscFunctionBegin; 215c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 216254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 217c5148223SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 218254c1ad2SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 219c5148223SMatthew G. Knepley ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 220254c1ad2SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 221c5148223SMatthew G. Knepley ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 222254c1ad2SMatthew G. Knepley for (c = cStart; c < cEndInterior; c++) { 223254c1ad2SMatthew G. Knepley const PetscInt *faces; 224c5148223SMatthew G. Knepley PetscInt numFaces, usedFaces, f, d; 225254c1ad2SMatthew G. Knepley const CellGeom *cg; 2268cd7fcbbSMatthew G. Knepley PetscBool boundary; 2278cd7fcbbSMatthew G. Knepley PetscInt ghost; 2288cd7fcbbSMatthew G. Knepley 229254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 230c5148223SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 231c5148223SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 232c5148223SMatthew G. Knepley if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction", c, numFaces); 233c5148223SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 234c5148223SMatthew G. Knepley const CellGeom *cg1; 235c5148223SMatthew G. Knepley FaceGeom *fg; 236254c1ad2SMatthew G. Knepley const PetscInt *fcells; 237254c1ad2SMatthew G. Knepley PetscInt ncell, side; 238254c1ad2SMatthew G. Knepley 239254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2408cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2418cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 242254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 243254c1ad2SMatthew G. Knepley side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 244254c1ad2SMatthew G. Knepley ncell = fcells[!side]; /* the neighbor */ 245254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 246254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 247c5148223SMatthew G. Knepley for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 248254c1ad2SMatthew G. Knepley gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 249254c1ad2SMatthew G. Knepley } 250254c1ad2SMatthew G. Knepley if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 251c5148223SMatthew G. Knepley ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 252c5148223SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 253254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2548cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2558cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 256c5148223SMatthew G. Knepley for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 257c5148223SMatthew G. Knepley ++usedFaces; 258254c1ad2SMatthew G. Knepley } 259254c1ad2SMatthew G. Knepley } 260c5148223SMatthew G. Knepley ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 261254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 262254c1ad2SMatthew G. Knepley } 263254c1ad2SMatthew G. Knepley 264254c1ad2SMatthew G. Knepley #undef __FUNCT__ 265254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGradient" 266*924a1b8fSMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGradient(DM dm, PetscFV fvm, DMTS dmts) 267254c1ad2SMatthew G. Knepley { 268*924a1b8fSMatthew G. Knepley DM dmFace, dmCell, dmGrad; 269*924a1b8fSMatthew G. Knepley Vec facegeom, cellgeom; 270254c1ad2SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 271254c1ad2SMatthew G. Knepley PetscSection sectionGrad; 272254c1ad2SMatthew G. Knepley PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 273254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 274254c1ad2SMatthew G. Knepley 275254c1ad2SMatthew G. Knepley PetscFunctionBegin; 276c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 277254c1ad2SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 278254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 279254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 280254c1ad2SMatthew G. Knepley /* Construct the interpolant corresponding to each face from the leat-square solution over the cell neighborhood */ 281*924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGeometry(dm, &facegeom, &cellgeom, NULL);CHKERRQ(ierr); 282*924a1b8fSMatthew G. Knepley ierr = VecGetDM(facegeom, &dmFace);CHKERRQ(ierr); 283*924a1b8fSMatthew G. Knepley ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr); 284*924a1b8fSMatthew G. Knepley ierr = VecGetArray(facegeom, &fgeom);CHKERRQ(ierr); 285*924a1b8fSMatthew G. Knepley ierr = VecGetArray(cellgeom, &cgeom);CHKERRQ(ierr); 286c5148223SMatthew G. Knepley ierr = BuildGradientReconstruction(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 287*924a1b8fSMatthew G. Knepley ierr = VecRestoreArray(facegeom, &fgeom);CHKERRQ(ierr); 288*924a1b8fSMatthew G. Knepley ierr = VecRestoreArray(cellgeom, &cgeom);CHKERRQ(ierr); 289254c1ad2SMatthew G. Knepley /* Create storage for gradients */ 290*924a1b8fSMatthew G. Knepley ierr = DMClone(dm, &dmGrad);CHKERRQ(ierr); 291254c1ad2SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 292254c1ad2SMatthew G. Knepley ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 293254c1ad2SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 294254c1ad2SMatthew G. Knepley ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 295*924a1b8fSMatthew G. Knepley ierr = DMSetDefaultSection(dmGrad, sectionGrad);CHKERRQ(ierr); 296254c1ad2SMatthew G. Knepley ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 297*924a1b8fSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad", (PetscObject) dmGrad);CHKERRQ(ierr); 298*924a1b8fSMatthew G. Knepley ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 299254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 300254c1ad2SMatthew G. Knepley } 301254c1ad2SMatthew G. Knepley 302254c1ad2SMatthew G. Knepley #undef __FUNCT__ 303254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM_Static" 304*924a1b8fSMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValuesFVM_Static(DM dm, PetscFV fvm, PetscReal time, Vec locX, Vec Grad) 305254c1ad2SMatthew G. Knepley { 306*924a1b8fSMatthew G. Knepley Vec faceGeometry, cellGeometry; 307*924a1b8fSMatthew G. Knepley DM dmFace, dmCell, dmGrad; 308254c1ad2SMatthew G. Knepley const PetscScalar *facegeom, *cellgeom, *grad; 309254c1ad2SMatthew G. Knepley PetscScalar *x, *fx; 310254c1ad2SMatthew G. Knepley PetscInt numBd, b, dim, pdim, fStart, fEnd; 311254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 312254c1ad2SMatthew G. Knepley 313254c1ad2SMatthew G. Knepley PetscFunctionBegin; 314c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 315*924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGeometry(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr); 316*924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGradientDM(dm, &dmGrad);CHKERRQ(ierr); 317254c1ad2SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 318254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 319254c1ad2SMatthew G. Knepley ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 320254c1ad2SMatthew G. Knepley if (Grad) { 321254c1ad2SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 322254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 323254c1ad2SMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 324254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 325254c1ad2SMatthew G. Knepley } 326254c1ad2SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 327254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 328254c1ad2SMatthew G. Knepley ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 329254c1ad2SMatthew G. Knepley for (b = 0; b < numBd; ++b) { 3301475bf87SMatthew G. Knepley PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*); 3318cd7fcbbSMatthew G. Knepley DMLabel label; 3328cd7fcbbSMatthew G. Knepley const char *labelname; 333254c1ad2SMatthew G. Knepley const PetscInt *ids; 334254c1ad2SMatthew G. Knepley PetscInt numids, i; 335254c1ad2SMatthew G. Knepley void *ctx; 336254c1ad2SMatthew G. Knepley 3378cd7fcbbSMatthew G. Knepley ierr = DMPlexGetBoundary(dm, b, NULL, NULL, &labelname, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr); 3388cd7fcbbSMatthew G. Knepley ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 339254c1ad2SMatthew G. Knepley for (i = 0; i < numids; ++i) { 340254c1ad2SMatthew G. Knepley IS faceIS; 341254c1ad2SMatthew G. Knepley const PetscInt *faces; 342254c1ad2SMatthew G. Knepley PetscInt numFaces, f; 343254c1ad2SMatthew G. Knepley 344254c1ad2SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 345254c1ad2SMatthew G. Knepley if (!faceIS) continue; /* No points with that id on this process */ 346254c1ad2SMatthew G. Knepley ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 347254c1ad2SMatthew G. Knepley ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 348254c1ad2SMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 349254c1ad2SMatthew G. Knepley const PetscInt face = faces[f], *cells; 350254c1ad2SMatthew G. Knepley const FaceGeom *fg; 351254c1ad2SMatthew G. Knepley 352254c1ad2SMatthew G. Knepley if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 353254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 354254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 355254c1ad2SMatthew G. Knepley if (Grad) { 356254c1ad2SMatthew G. Knepley const CellGeom *cg; 357254c1ad2SMatthew G. Knepley const PetscScalar *cx, *cgrad; 3581475bf87SMatthew G. Knepley PetscScalar *xG; 3591475bf87SMatthew G. Knepley PetscReal dx[3]; 360254c1ad2SMatthew G. Knepley PetscInt d; 361254c1ad2SMatthew G. Knepley 362254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 363254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 364*924a1b8fSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 365254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 366254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cg->centroid, fg->centroid, dx); 367254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DotD(dim, &cgrad[d*dim], dx); 368254c1ad2SMatthew G. Knepley ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr); 369254c1ad2SMatthew G. Knepley } else { 370254c1ad2SMatthew G. Knepley const PetscScalar *xI; 371254c1ad2SMatthew G. Knepley PetscScalar *xG; 372254c1ad2SMatthew G. Knepley 373254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 374254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 375254c1ad2SMatthew G. Knepley ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr); 376254c1ad2SMatthew G. Knepley } 377254c1ad2SMatthew G. Knepley } 378254c1ad2SMatthew G. Knepley ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 379254c1ad2SMatthew G. Knepley ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 380254c1ad2SMatthew G. Knepley } 381254c1ad2SMatthew G. Knepley } 382254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 383254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 384254c1ad2SMatthew G. Knepley if (Grad) { 385254c1ad2SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 386254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 387254c1ad2SMatthew G. Knepley } 388254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 389254c1ad2SMatthew G. Knepley } 390254c1ad2SMatthew G. Knepley 391254c1ad2SMatthew G. Knepley #undef __FUNCT__ 392*924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM" 393*924a1b8fSMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *ctx) 394254c1ad2SMatthew G. Knepley { 395*924a1b8fSMatthew G. Knepley PetscDS prob; 396*924a1b8fSMatthew G. Knepley void (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *); 397*924a1b8fSMatthew G. Knepley void *rctx; 398254c1ad2SMatthew G. Knepley PetscFV fvm; 399254c1ad2SMatthew G. Knepley PetscLimiter lim; 400*924a1b8fSMatthew G. Knepley Vec faceGeometry, cellGeometry; 401*924a1b8fSMatthew G. Knepley Vec Grad = NULL, locGrad; 402*924a1b8fSMatthew G. Knepley DM dmFace, dmCell, dmGrad; 4038cd7fcbbSMatthew G. Knepley DMLabel ghostLabel; 404254c1ad2SMatthew G. Knepley PetscCellGeometry fgeom, cgeom; 405254c1ad2SMatthew G. Knepley const PetscScalar *facegeom, *cellgeom, *x, *lgrad; 406254c1ad2SMatthew G. Knepley PetscScalar *grad, *f, *uL, *uR, *fluxL, *fluxR; 407254c1ad2SMatthew G. Knepley PetscReal *centroid, *normal, *vol, *cellPhi; 408254c1ad2SMatthew G. Knepley PetscBool computeGradients; 409254c1ad2SMatthew G. Knepley PetscInt Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface, cell, cStart, cEnd, cEndInterior; 410254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 411254c1ad2SMatthew G. Knepley 412254c1ad2SMatthew G. Knepley PetscFunctionBegin; 413*924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 414*924a1b8fSMatthew G. Knepley PetscValidHeaderSpecific(locX,VEC_CLASSID,3); 415254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(F,VEC_CLASSID,5); 416*924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGeometry(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr); 417*924a1b8fSMatthew G. Knepley ierr = DMPlexTSGetGradientDM(dm, &dmGrad);CHKERRQ(ierr); 418*924a1b8fSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 419*924a1b8fSMatthew G. Knepley ierr = PetscDSGetRiemannSolver(prob, 0, &riemann);CHKERRQ(ierr); 420*924a1b8fSMatthew G. Knepley ierr = PetscDSGetContext(prob, 0, &rctx);CHKERRQ(ierr); 421c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 4226dbbd306SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 4236dbbd306SMatthew G. Knepley ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 424254c1ad2SMatthew G. Knepley ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr); 4256dbbd306SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 426254c1ad2SMatthew G. Knepley ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); 427254c1ad2SMatthew G. Knepley if (computeGradients) { 428*924a1b8fSMatthew G. Knepley ierr = DMGetGlobalVector(dmGrad, &Grad);CHKERRQ(ierr); 429254c1ad2SMatthew G. Knepley ierr = VecZeroEntries(Grad);CHKERRQ(ierr); 430254c1ad2SMatthew G. Knepley ierr = VecGetArray(Grad, &grad);CHKERRQ(ierr); 431254c1ad2SMatthew G. Knepley } 4326dbbd306SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 4336dbbd306SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 4346dbbd306SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 4356dbbd306SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 4366dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 4376dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 4386dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); 439254c1ad2SMatthew G. Knepley /* Count faces and reconstruct gradients */ 4406dbbd306SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 441254c1ad2SMatthew G. Knepley const PetscInt *cells; 442254c1ad2SMatthew G. Knepley const FaceGeom *fg; 443254c1ad2SMatthew G. Knepley const PetscScalar *cx[2]; 444254c1ad2SMatthew G. Knepley PetscScalar *cgrad[2]; 4458cd7fcbbSMatthew G. Knepley PetscBool boundary; 4468cd7fcbbSMatthew G. Knepley PetscInt ghost, c, pd, d; 4476dbbd306SMatthew G. Knepley 4486dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 4496dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 4506dbbd306SMatthew G. Knepley ++numFaces; 451254c1ad2SMatthew G. Knepley if (!computeGradients) continue; 4528cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 4538cd7fcbbSMatthew G. Knepley if (boundary) continue; 454254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 455254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 456254c1ad2SMatthew G. Knepley for (c = 0; c < 2; ++c) { 457254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr); 458*924a1b8fSMatthew G. Knepley ierr = DMPlexPointGlobalRef(dmGrad, cells[c], grad, &cgrad[c]);CHKERRQ(ierr); 459254c1ad2SMatthew G. Knepley } 460254c1ad2SMatthew G. Knepley for (pd = 0; pd < pdim; ++pd) { 461254c1ad2SMatthew G. Knepley PetscScalar delta = cx[1][pd] - cx[0][pd]; 462254c1ad2SMatthew G. Knepley 463254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) { 464254c1ad2SMatthew G. Knepley if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta; 465254c1ad2SMatthew G. Knepley if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta; 466254c1ad2SMatthew G. Knepley } 467254c1ad2SMatthew G. Knepley } 468254c1ad2SMatthew G. Knepley } 469254c1ad2SMatthew G. Knepley /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */ 470254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 471254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 472254c1ad2SMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 473254c1ad2SMatthew G. Knepley for (cell = computeGradients && lim ? cStart : cEnd; cell < cEndInterior; ++cell) { 474254c1ad2SMatthew G. Knepley const PetscInt *faces; 475254c1ad2SMatthew G. Knepley const PetscScalar *cx; 476254c1ad2SMatthew G. Knepley const CellGeom *cg; 477254c1ad2SMatthew G. Knepley PetscScalar *cgrad; 478254c1ad2SMatthew G. Knepley PetscInt coneSize, f, pd, d; 479254c1ad2SMatthew G. Knepley 480254c1ad2SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 481254c1ad2SMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 482254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr); 483254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr); 484*924a1b8fSMatthew G. Knepley ierr = DMPlexPointGlobalRef(dmGrad, cell, grad, &cgrad);CHKERRQ(ierr); 48571ea8452SMatthew G. Knepley if (!cgrad) continue; /* Unowned overlap cell, we do not compute */ 486254c1ad2SMatthew G. Knepley /* Limiter will be minimum value over all neighbors */ 487254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) cellPhi[d] = PETSC_MAX_REAL; 488254c1ad2SMatthew G. Knepley for (f = 0; f < coneSize; ++f) { 489254c1ad2SMatthew G. Knepley const PetscScalar *ncx; 490254c1ad2SMatthew G. Knepley const CellGeom *ncg; 491254c1ad2SMatthew G. Knepley const PetscInt *fcells; 4928cd7fcbbSMatthew G. Knepley PetscInt face = faces[f], ncell, ghost; 4931475bf87SMatthew G. Knepley PetscReal v[3]; 4948cd7fcbbSMatthew G. Knepley PetscBool boundary; 495254c1ad2SMatthew G. Knepley 496254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 4978cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 4988cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 499254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr); 500254c1ad2SMatthew G. Knepley ncell = cell == fcells[0] ? fcells[1] : fcells[0]; 501254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr); 502254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr); 503254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cg->centroid, ncg->centroid, v); 504254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 505254c1ad2SMatthew G. Knepley /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */ 5061475bf87SMatthew G. Knepley PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DotD(dim, &cgrad[d*dim], v); 507254c1ad2SMatthew G. Knepley 508254c1ad2SMatthew G. Knepley ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr); 509254c1ad2SMatthew G. Knepley cellPhi[d] = PetscMin(cellPhi[d], phi); 510254c1ad2SMatthew G. Knepley } 511254c1ad2SMatthew G. Knepley } 512254c1ad2SMatthew G. Knepley /* Apply limiter to gradient */ 513254c1ad2SMatthew G. Knepley for (pd = 0; pd < pdim; ++pd) 514254c1ad2SMatthew G. Knepley /* Scalar limiter applied to each component separately */ 515254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd]; 516254c1ad2SMatthew G. Knepley } 517254c1ad2SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 518*924a1b8fSMatthew G. Knepley ierr = DMPlexInsertBoundaryValuesFVM_Static(dm, fvm, time, locX, Grad);CHKERRQ(ierr); 519254c1ad2SMatthew G. Knepley if (computeGradients) { 520254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(Grad, &grad);CHKERRQ(ierr); 521*924a1b8fSMatthew G. Knepley ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 522*924a1b8fSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 523*924a1b8fSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 524*924a1b8fSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmGrad, &Grad);CHKERRQ(ierr); 525254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 5266dbbd306SMatthew G. Knepley } 5276dbbd306SMatthew G. Knepley ierr = PetscMalloc7(numFaces*dim,¢roid,numFaces*dim,&normal,numFaces*2,&vol,numFaces*pdim,&uL,numFaces*pdim,&uR,numFaces*pdim,&fluxL,numFaces*pdim,&fluxR);CHKERRQ(ierr); 528254c1ad2SMatthew G. Knepley /* Read out values */ 5296dbbd306SMatthew G. Knepley for (face = fStart, iface = 0; face < fEnd; ++face) { 5306dbbd306SMatthew G. Knepley const PetscInt *cells; 5316dbbd306SMatthew G. Knepley const FaceGeom *fg; 5326dbbd306SMatthew G. Knepley const CellGeom *cgL, *cgR; 533254c1ad2SMatthew G. Knepley const PetscScalar *xL, *xR, *gL, *gR; 5346dbbd306SMatthew G. Knepley PetscInt ghost, d; 5356dbbd306SMatthew G. Knepley 5366dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 5376dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 5386dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 5396dbbd306SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 5406dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); 5416dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); 5426dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr); 5436dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr); 544254c1ad2SMatthew G. Knepley if (computeGradients) { 5451475bf87SMatthew G. Knepley PetscReal dxL[3], dxR[3]; 546254c1ad2SMatthew G. Knepley 547*924a1b8fSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr); 548*924a1b8fSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr); 549254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cgL->centroid, fg->centroid, dxL); 550254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cgR->centroid, fg->centroid, dxR); 551254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 552254c1ad2SMatthew G. Knepley uL[iface*pdim+d] = xL[d] + DotD(dim, &gL[d*dim], dxL); 553254c1ad2SMatthew G. Knepley uR[iface*pdim+d] = xR[d] + DotD(dim, &gR[d*dim], dxR); 554254c1ad2SMatthew G. Knepley } 555254c1ad2SMatthew G. Knepley } else { 5566dbbd306SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 5576dbbd306SMatthew G. Knepley uL[iface*pdim+d] = xL[d]; 5586dbbd306SMatthew G. Knepley uR[iface*pdim+d] = xR[d]; 5596dbbd306SMatthew G. Knepley } 560254c1ad2SMatthew G. Knepley } 5616dbbd306SMatthew G. Knepley for (d = 0; d < dim; ++d) { 5626dbbd306SMatthew G. Knepley centroid[iface*dim+d] = fg->centroid[d]; 5636dbbd306SMatthew G. Knepley normal[iface*dim+d] = fg->normal[d]; 5646dbbd306SMatthew G. Knepley } 5656dbbd306SMatthew G. Knepley vol[iface*2+0] = cgL->volume; 5666dbbd306SMatthew G. Knepley vol[iface*2+1] = cgR->volume; 5676dbbd306SMatthew G. Knepley ++iface; 5686dbbd306SMatthew G. Knepley } 569254c1ad2SMatthew G. Knepley if (computeGradients) { 570254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(locGrad,&lgrad);CHKERRQ(ierr); 571*924a1b8fSMatthew G. Knepley ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 572254c1ad2SMatthew G. Knepley } 5736dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); 5746dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 5756dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 5766dbbd306SMatthew G. Knepley fgeom.v0 = centroid; 5776dbbd306SMatthew G. Knepley fgeom.n = normal; 5786dbbd306SMatthew G. Knepley cgeom.vol = vol; 579254c1ad2SMatthew G. Knepley /* Riemann solve */ 580*924a1b8fSMatthew G. Knepley ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, cgeom, uL, uR, riemann, fluxL, fluxR, rctx);CHKERRQ(ierr); 581254c1ad2SMatthew G. Knepley /* Insert fluxes */ 5826dbbd306SMatthew G. Knepley ierr = VecGetArray(F, &f);CHKERRQ(ierr); 5836dbbd306SMatthew G. Knepley for (face = fStart, iface = 0; face < fEnd; ++face) { 5846dbbd306SMatthew G. Knepley const PetscInt *cells; 5856dbbd306SMatthew G. Knepley PetscScalar *fL, *fR; 5866dbbd306SMatthew G. Knepley PetscInt ghost, d; 5876dbbd306SMatthew G. Knepley 5886dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 5896dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 5906dbbd306SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 5916dbbd306SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr); 5926dbbd306SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr); 5936dbbd306SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 5946dbbd306SMatthew G. Knepley if (fL) fL[d] -= fluxL[iface*pdim+d]; 5956dbbd306SMatthew G. Knepley if (fR) fR[d] += fluxR[iface*pdim+d]; 5966dbbd306SMatthew G. Knepley } 5976dbbd306SMatthew G. Knepley ++iface; 5986dbbd306SMatthew G. Knepley } 5996dbbd306SMatthew G. Knepley ierr = VecRestoreArray(F, &f);CHKERRQ(ierr); 6006dbbd306SMatthew G. Knepley ierr = PetscFree7(centroid,normal,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr); 601254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 602254c1ad2SMatthew G. Knepley } 603254c1ad2SMatthew G. Knepley 604254c1ad2SMatthew G. Knepley #undef __FUNCT__ 605254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetRHSFunctionLocal" 606254c1ad2SMatthew G. Knepley /*@C 607254c1ad2SMatthew G. Knepley DMPlexTSSetRHSFunctionLocal - set a local residual evaluation function 608254c1ad2SMatthew G. Knepley 609254c1ad2SMatthew G. Knepley Logically Collective 610254c1ad2SMatthew G. Knepley 611254c1ad2SMatthew G. Knepley Input Arguments: 612254c1ad2SMatthew G. Knepley + dm - DM to associate callback with 613*924a1b8fSMatthew G. Knepley . func - the RHS evaluation function 614*924a1b8fSMatthew G. Knepley - ctx - an optional user context 615254c1ad2SMatthew G. Knepley 616254c1ad2SMatthew G. Knepley Level: beginner 617254c1ad2SMatthew G. Knepley 618254c1ad2SMatthew G. Knepley .seealso: DMTSSetRHSFunctionLocal() 619254c1ad2SMatthew G. Knepley @*/ 620*924a1b8fSMatthew G. Knepley PetscErrorCode DMPlexTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM dm, PetscReal time, Vec X, Vec F, void *ctx), void *ctx) 621254c1ad2SMatthew G. Knepley { 622254c1ad2SMatthew G. Knepley DMTS dmts; 623254c1ad2SMatthew G. Knepley PetscFV fvm; 624254c1ad2SMatthew G. Knepley PetscInt Nf; 625254c1ad2SMatthew G. Knepley PetscBool computeGradients; 626254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 627254c1ad2SMatthew G. Knepley 628254c1ad2SMatthew G. Knepley PetscFunctionBegin; 629254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 630254c1ad2SMatthew G. Knepley ierr = DMGetDMTSWrite(dm, &dmts);CHKERRQ(ierr); 631254c1ad2SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 632254c1ad2SMatthew G. Knepley ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 633*924a1b8fSMatthew G. Knepley ierr = DMPlexTSSetupGeometry(dm, fvm, dmts);CHKERRQ(ierr); 634254c1ad2SMatthew G. Knepley ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); 635*924a1b8fSMatthew G. Knepley if (computeGradients) {ierr = DMPlexTSSetupGradient(dm, fvm, dmts);CHKERRQ(ierr);} 636*924a1b8fSMatthew G. Knepley ierr = DMTSSetRHSFunctionLocal(dm, func, ctx);CHKERRQ(ierr); 6376dbbd306SMatthew G. Knepley PetscFunctionReturn(0); 6386dbbd306SMatthew G. Knepley } 639