16dbbd306SMatthew G. Knepley #include <petscdmplex.h> /*I "petscdmplex.h" I*/ 2254c1ad2SMatthew G. Knepley #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 36dbbd306SMatthew G. Knepley #include <petscfv.h> 46dbbd306SMatthew G. Knepley 5*1475bf87SMatthew 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];} 6*1475bf87SMatthew 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;} 7*1475bf87SMatthew 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;} 8*1475bf87SMatthew 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);} 9254c1ad2SMatthew G. Knepley 10254c1ad2SMatthew G. Knepley typedef struct { 11254c1ad2SMatthew G. Knepley PetscBool setupGeom; /* Flag for geometry setup */ 12254c1ad2SMatthew G. Knepley PetscBool setupGrad; /* Flag for gradient calculation setup */ 13254c1ad2SMatthew G. Knepley Vec facegeom; /* FaceGeom struct for each face */ 14254c1ad2SMatthew G. Knepley Vec cellgeom; /* CellGeom struct for each cell */ 15254c1ad2SMatthew G. Knepley DM dmGrad; /* Layout for the gradient data */ 16254c1ad2SMatthew G. Knepley PetscReal minradius; /* Minimum distance from centroid to face */ 17254c1ad2SMatthew G. Knepley void (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *); 18254c1ad2SMatthew G. Knepley void *rhsfunctionlocalctx; 19254c1ad2SMatthew G. Knepley } DMTS_Plex; 206dbbd306SMatthew G. Knepley 216dbbd306SMatthew G. Knepley #undef __FUNCT__ 22254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMTSDestroy_Plex" 23254c1ad2SMatthew G. Knepley static PetscErrorCode DMTSDestroy_Plex(DMTS dmts) 246dbbd306SMatthew G. Knepley { 25254c1ad2SMatthew G. Knepley DMTS_Plex *dmplexts = (DMTS_Plex *) dmts->data; 266dbbd306SMatthew G. Knepley PetscErrorCode ierr; 276dbbd306SMatthew G. Knepley 28254c1ad2SMatthew G. Knepley PetscFunctionBegin; 29254c1ad2SMatthew G. Knepley ierr = DMDestroy(&dmplexts->dmGrad);CHKERRQ(ierr); 30254c1ad2SMatthew G. Knepley ierr = VecDestroy(&dmplexts->facegeom);CHKERRQ(ierr); 31254c1ad2SMatthew G. Knepley ierr = VecDestroy(&dmplexts->cellgeom);CHKERRQ(ierr); 32254c1ad2SMatthew G. Knepley ierr = PetscFree(dmts->data);CHKERRQ(ierr); 33254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 34254c1ad2SMatthew G. Knepley } 35254c1ad2SMatthew G. Knepley 36254c1ad2SMatthew G. Knepley #undef __FUNCT__ 37254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMTSDuplicate_Plex" 38254c1ad2SMatthew G. Knepley static PetscErrorCode DMTSDuplicate_Plex(DMTS olddmts, DMTS dmts) 39254c1ad2SMatthew G. Knepley { 40254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 41254c1ad2SMatthew G. Knepley 42254c1ad2SMatthew G. Knepley PetscFunctionBegin; 43254c1ad2SMatthew G. Knepley ierr = PetscNewLog(dmts, (DMTS_Plex **) &dmts->data);CHKERRQ(ierr); 44254c1ad2SMatthew G. Knepley if (olddmts->data) {ierr = PetscMemcpy(dmts->data, olddmts->data, sizeof(DMTS_Plex));CHKERRQ(ierr);} 45254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 46254c1ad2SMatthew G. Knepley } 47254c1ad2SMatthew G. Knepley 48254c1ad2SMatthew G. Knepley #undef __FUNCT__ 49254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetContext" 50254c1ad2SMatthew G. Knepley static PetscErrorCode DMPlexTSGetContext(DM dm, DMTS dmts, DMTS_Plex **dmplexts) 51254c1ad2SMatthew G. Knepley { 52254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 53254c1ad2SMatthew G. Knepley 54254c1ad2SMatthew G. Knepley PetscFunctionBegin; 55254c1ad2SMatthew G. Knepley *dmplexts = NULL; 56254c1ad2SMatthew G. Knepley if (!dmts->data) { 57254c1ad2SMatthew G. Knepley ierr = PetscNewLog(dm, (DMTS_Plex **) &dmts->data);CHKERRQ(ierr); 58254c1ad2SMatthew G. Knepley dmts->ops->destroy = DMTSDestroy_Plex; 59254c1ad2SMatthew G. Knepley dmts->ops->duplicate = DMTSDuplicate_Plex; 60254c1ad2SMatthew G. Knepley } 61254c1ad2SMatthew G. Knepley *dmplexts = (DMTS_Plex *) dmts->data; 62254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 63254c1ad2SMatthew G. Knepley } 64254c1ad2SMatthew G. Knepley 65254c1ad2SMatthew G. Knepley #undef __FUNCT__ 66a0ac79e7SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometry" 67a0ac79e7SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometry(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius) 68a0ac79e7SMatthew G. Knepley { 69a0ac79e7SMatthew G. Knepley DMTS dmts; 70a0ac79e7SMatthew G. Knepley DMTS_Plex *dmplexts; 71a0ac79e7SMatthew G. Knepley PetscErrorCode ierr; 72a0ac79e7SMatthew G. Knepley 73a0ac79e7SMatthew G. Knepley PetscFunctionBegin; 74a0ac79e7SMatthew G. Knepley ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr); 75a0ac79e7SMatthew G. Knepley ierr = DMPlexTSGetContext(dm, dmts, &dmplexts);CHKERRQ(ierr); 76a0ac79e7SMatthew G. Knepley if (facegeom) *facegeom = dmplexts->facegeom; 77a0ac79e7SMatthew G. Knepley if (cellgeom) *cellgeom = dmplexts->cellgeom; 78a0ac79e7SMatthew G. Knepley if (minRadius) *minRadius = dmplexts->minradius; 79a0ac79e7SMatthew G. Knepley PetscFunctionReturn(0); 80a0ac79e7SMatthew G. Knepley } 81a0ac79e7SMatthew G. Knepley 82a0ac79e7SMatthew G. Knepley #undef __FUNCT__ 83254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGeometry" 84254c1ad2SMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGeometry(DM dm, PetscFV fvm, DMTS_Plex *dmplexts) 85254c1ad2SMatthew G. Knepley { 86254c1ad2SMatthew G. Knepley DM dmFace, dmCell; 87254c1ad2SMatthew G. Knepley DMLabel ghostLabel; 88254c1ad2SMatthew G. Knepley PetscSection sectionFace, sectionCell; 89254c1ad2SMatthew G. Knepley PetscSection coordSection; 90254c1ad2SMatthew G. Knepley Vec coordinates; 91254c1ad2SMatthew G. Knepley PetscReal minradius; 92254c1ad2SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 93254c1ad2SMatthew G. Knepley PetscInt dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f; 94254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 95254c1ad2SMatthew G. Knepley 96254c1ad2SMatthew G. Knepley PetscFunctionBegin; 97254c1ad2SMatthew G. Knepley if (dmplexts->setupGeom) PetscFunctionReturn(0); 98254c1ad2SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 99254c1ad2SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 100254c1ad2SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 101254c1ad2SMatthew G. Knepley /* Make cell centroids and volumes */ 102254c1ad2SMatthew G. Knepley ierr = DMClone(dm, &dmCell);CHKERRQ(ierr); 103254c1ad2SMatthew G. Knepley ierr = DMSetCoordinateSection(dmCell, coordSection);CHKERRQ(ierr); 104254c1ad2SMatthew G. Knepley ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr); 105254c1ad2SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionCell);CHKERRQ(ierr); 106254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 107254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 108254c1ad2SMatthew G. Knepley ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr); 109254c1ad2SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, sizeof(CellGeom)/sizeof(PetscScalar));CHKERRQ(ierr);} 110254c1ad2SMatthew G. Knepley ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr); 111254c1ad2SMatthew G. Knepley ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr); 112254c1ad2SMatthew G. Knepley ierr = PetscSectionDestroy(§ionCell);CHKERRQ(ierr); 113254c1ad2SMatthew G. Knepley ierr = DMCreateLocalVector(dmCell, &dmplexts->cellgeom);CHKERRQ(ierr); 114254c1ad2SMatthew G. Knepley ierr = VecGetArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr); 115254c1ad2SMatthew G. Knepley for (c = cStart; c < cEndInterior; ++c) { 116254c1ad2SMatthew G. Knepley CellGeom *cg; 117254c1ad2SMatthew G. Knepley 118254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 119254c1ad2SMatthew G. Knepley ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr); 120254c1ad2SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr); 121254c1ad2SMatthew G. Knepley } 122254c1ad2SMatthew G. Knepley /* Compute face normals and minimum cell radius */ 123254c1ad2SMatthew G. Knepley ierr = DMClone(dm, &dmFace);CHKERRQ(ierr); 124254c1ad2SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionFace);CHKERRQ(ierr); 125254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 126254c1ad2SMatthew G. Knepley ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr); 127254c1ad2SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, sizeof(FaceGeom)/sizeof(PetscScalar));CHKERRQ(ierr);} 128254c1ad2SMatthew G. Knepley ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr); 129254c1ad2SMatthew G. Knepley ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr); 130254c1ad2SMatthew G. Knepley ierr = PetscSectionDestroy(§ionFace);CHKERRQ(ierr); 131254c1ad2SMatthew G. Knepley ierr = DMCreateLocalVector(dmFace, &dmplexts->facegeom);CHKERRQ(ierr); 132254c1ad2SMatthew G. Knepley ierr = VecGetArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr); 133254c1ad2SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 134254c1ad2SMatthew G. Knepley minradius = PETSC_MAX_REAL; 135254c1ad2SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 136254c1ad2SMatthew G. Knepley FaceGeom *fg; 137254c1ad2SMatthew G. Knepley PetscReal area; 138254c1ad2SMatthew G. Knepley PetscInt ghost, d; 139254c1ad2SMatthew G. Knepley 140254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr); 141254c1ad2SMatthew G. Knepley if (ghost >= 0) continue; 142254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr); 143254c1ad2SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr); 144254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] *= area; 145254c1ad2SMatthew G. Knepley /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */ 146254c1ad2SMatthew G. Knepley { 147254c1ad2SMatthew G. Knepley CellGeom *cL, *cR; 148254c1ad2SMatthew G. Knepley const PetscInt *cells; 149254c1ad2SMatthew G. Knepley PetscReal *lcentroid, *rcentroid; 150*1475bf87SMatthew G. Knepley PetscReal v[3]; 151254c1ad2SMatthew G. Knepley 152254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr); 153254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr); 154254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr); 155254c1ad2SMatthew G. Knepley lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid; 156254c1ad2SMatthew G. Knepley rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid; 157254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, lcentroid, rcentroid, v); 158*1475bf87SMatthew G. Knepley if (DotRealD(dim, fg->normal, v) < 0) { 159254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d]; 160254c1ad2SMatthew G. Knepley } 161*1475bf87SMatthew G. Knepley if (DotRealD(dim, fg->normal, v) <= 0) { 162254c1ad2SMatthew 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, fg->normal[0], fg->normal[1], v[0], v[1]); 163254c1ad2SMatthew 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, fg->normal[0], fg->normal[1], fg->normal[2], v[0], v[1], v[2]); 164254c1ad2SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f); 165254c1ad2SMatthew G. Knepley } 166254c1ad2SMatthew G. Knepley if (cells[0] < cEndInterior) { 167254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, fg->centroid, cL->centroid, v); 168254c1ad2SMatthew G. Knepley minradius = PetscMin(minradius, NormD(dim, v)); 169254c1ad2SMatthew G. Knepley } 170254c1ad2SMatthew G. Knepley if (cells[1] < cEndInterior) { 171254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, fg->centroid, cR->centroid, v); 172254c1ad2SMatthew G. Knepley minradius = PetscMin(minradius, NormD(dim, v)); 173254c1ad2SMatthew G. Knepley } 174254c1ad2SMatthew G. Knepley } 175254c1ad2SMatthew G. Knepley } 176254c1ad2SMatthew G. Knepley ierr = MPI_Allreduce(&minradius, &dmplexts->minradius, 1, MPIU_SCALAR, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 177254c1ad2SMatthew G. Knepley /* Compute centroids of ghost cells */ 178254c1ad2SMatthew G. Knepley for (c = cEndInterior; c < cEnd; ++c) { 179254c1ad2SMatthew G. Knepley FaceGeom *fg; 180254c1ad2SMatthew G. Knepley const PetscInt *cone, *support; 181254c1ad2SMatthew G. Knepley PetscInt coneSize, supportSize, s; 182254c1ad2SMatthew G. Knepley 183254c1ad2SMatthew G. Knepley ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr); 184254c1ad2SMatthew G. Knepley if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize); 185254c1ad2SMatthew G. Knepley ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr); 186254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr); 187254c1ad2SMatthew G. Knepley if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 1", cone[0], supportSize); 188254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr); 189254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr); 190254c1ad2SMatthew G. Knepley for (s = 0; s < 2; ++s) { 191254c1ad2SMatthew G. Knepley /* Reflect ghost centroid across plane of face */ 192254c1ad2SMatthew G. Knepley if (support[s] == c) { 193254c1ad2SMatthew G. Knepley const CellGeom *ci; 194254c1ad2SMatthew G. Knepley CellGeom *cg; 195*1475bf87SMatthew G. Knepley PetscReal c2f[3], a; 196254c1ad2SMatthew G. Knepley 197254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr); 198254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */ 199*1475bf87SMatthew G. Knepley a = DotRealD(dim, c2f, fg->normal)/DotRealD(dim, fg->normal, fg->normal); 200254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr); 201254c1ad2SMatthew G. Knepley WaxpyD(dim, 2*a, fg->normal, ci->centroid, cg->centroid); 202254c1ad2SMatthew G. Knepley cg->volume = ci->volume; 203254c1ad2SMatthew G. Knepley } 204254c1ad2SMatthew G. Knepley } 205254c1ad2SMatthew G. Knepley } 206254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr); 207254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr); 208254c1ad2SMatthew G. Knepley ierr = DMDestroy(&dmCell);CHKERRQ(ierr); 209254c1ad2SMatthew G. Knepley ierr = DMDestroy(&dmFace);CHKERRQ(ierr); 210254c1ad2SMatthew G. Knepley dmplexts->setupGeom = PETSC_TRUE; 211254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 212254c1ad2SMatthew G. Knepley } 213254c1ad2SMatthew G. Knepley 214254c1ad2SMatthew G. Knepley #undef __FUNCT__ 215c5148223SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction" 216c5148223SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 217254c1ad2SMatthew G. Knepley { 2188cd7fcbbSMatthew G. Knepley DMLabel ghostLabel; 219c5148223SMatthew G. Knepley PetscScalar *dx, *grad, **gref; 220c5148223SMatthew G. Knepley PetscInt dim, cStart, cEnd, c, cEndInterior, maxNumFaces; 221254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 222254c1ad2SMatthew G. Knepley 223254c1ad2SMatthew G. Knepley PetscFunctionBegin; 224254c1ad2SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 225254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 226c5148223SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 227254c1ad2SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr); 228c5148223SMatthew G. Knepley ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr); 229254c1ad2SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 230c5148223SMatthew G. Knepley ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr); 231254c1ad2SMatthew G. Knepley for (c = cStart; c < cEndInterior; c++) { 232254c1ad2SMatthew G. Knepley const PetscInt *faces; 233c5148223SMatthew G. Knepley PetscInt numFaces, usedFaces, f, d; 234254c1ad2SMatthew G. Knepley const CellGeom *cg; 2358cd7fcbbSMatthew G. Knepley PetscBool boundary; 2368cd7fcbbSMatthew G. Knepley PetscInt ghost; 2378cd7fcbbSMatthew G. Knepley 238254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr); 239c5148223SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr); 240c5148223SMatthew G. Knepley ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr); 241c5148223SMatthew 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); 242c5148223SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 243c5148223SMatthew G. Knepley const CellGeom *cg1; 244c5148223SMatthew G. Knepley FaceGeom *fg; 245254c1ad2SMatthew G. Knepley const PetscInt *fcells; 246254c1ad2SMatthew G. Knepley PetscInt ncell, side; 247254c1ad2SMatthew G. Knepley 248254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2498cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2508cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 251254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr); 252254c1ad2SMatthew G. Knepley side = (c != fcells[0]); /* c is on left=0 or right=1 of face */ 253254c1ad2SMatthew G. Knepley ncell = fcells[!side]; /* the neighbor */ 254254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr); 255254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr); 256c5148223SMatthew G. Knepley for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d]; 257254c1ad2SMatthew G. Knepley gref[usedFaces++] = fg->grad[side]; /* Gradient reconstruction term will go here */ 258254c1ad2SMatthew G. Knepley } 259254c1ad2SMatthew G. Knepley if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?"); 260c5148223SMatthew G. Knepley ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr); 261c5148223SMatthew G. Knepley for (f = 0, usedFaces = 0; f < numFaces; ++f) { 262254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr); 2638cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr); 2648cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 265c5148223SMatthew G. Knepley for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d]; 266c5148223SMatthew G. Knepley ++usedFaces; 267254c1ad2SMatthew G. Knepley } 268254c1ad2SMatthew G. Knepley } 269c5148223SMatthew G. Knepley ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr); 270254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 271254c1ad2SMatthew G. Knepley } 272254c1ad2SMatthew G. Knepley 273254c1ad2SMatthew G. Knepley #undef __FUNCT__ 274254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGradient" 275254c1ad2SMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGradient(DM dm, PetscFV fvm, DMTS_Plex *dmplexts) 276254c1ad2SMatthew G. Knepley { 277254c1ad2SMatthew G. Knepley DM dmFace, dmCell; 278254c1ad2SMatthew G. Knepley PetscScalar *fgeom, *cgeom; 279254c1ad2SMatthew G. Knepley PetscSection sectionGrad; 280254c1ad2SMatthew G. Knepley PetscInt dim, pdim, cStart, cEnd, cEndInterior, c; 281254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 282254c1ad2SMatthew G. Knepley 283254c1ad2SMatthew G. Knepley PetscFunctionBegin; 284254c1ad2SMatthew G. Knepley if (dmplexts->setupGrad) PetscFunctionReturn(0); 285254c1ad2SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 286254c1ad2SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 287254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 288254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 289254c1ad2SMatthew G. Knepley /* Construct the interpolant corresponding to each face from the leat-square solution over the cell neighborhood */ 290254c1ad2SMatthew G. Knepley ierr = VecGetDM(dmplexts->facegeom, &dmFace);CHKERRQ(ierr); 291254c1ad2SMatthew G. Knepley ierr = VecGetDM(dmplexts->cellgeom, &dmCell);CHKERRQ(ierr); 292254c1ad2SMatthew G. Knepley ierr = VecGetArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr); 293254c1ad2SMatthew G. Knepley ierr = VecGetArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr); 294c5148223SMatthew G. Knepley ierr = BuildGradientReconstruction(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr); 295254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr); 296254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr); 297254c1ad2SMatthew G. Knepley /* Create storage for gradients */ 298254c1ad2SMatthew G. Knepley ierr = DMClone(dm, &dmplexts->dmGrad);CHKERRQ(ierr); 299254c1ad2SMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ionGrad);CHKERRQ(ierr); 300254c1ad2SMatthew G. Knepley ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr); 301254c1ad2SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);} 302254c1ad2SMatthew G. Knepley ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr); 303254c1ad2SMatthew G. Knepley ierr = DMSetDefaultSection(dmplexts->dmGrad, sectionGrad);CHKERRQ(ierr); 304254c1ad2SMatthew G. Knepley ierr = PetscSectionDestroy(§ionGrad);CHKERRQ(ierr); 305254c1ad2SMatthew G. Knepley dmplexts->setupGrad = PETSC_TRUE; 306254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 307254c1ad2SMatthew G. Knepley } 308254c1ad2SMatthew G. Knepley 309254c1ad2SMatthew G. Knepley #undef __FUNCT__ 310254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM_Static" 311254c1ad2SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValuesFVM_Static(DM dm, PetscFV fvm, PetscReal time, Vec locX, Vec Grad, DMTS_Plex *dmplexts) 312254c1ad2SMatthew G. Knepley { 313254c1ad2SMatthew G. Knepley Vec faceGeometry = dmplexts->facegeom; 314254c1ad2SMatthew G. Knepley Vec cellGeometry = dmplexts->cellgeom; 315254c1ad2SMatthew G. Knepley DM dmFace, dmCell; 316254c1ad2SMatthew G. Knepley const PetscScalar *facegeom, *cellgeom, *grad; 317254c1ad2SMatthew G. Knepley PetscScalar *x, *fx; 318254c1ad2SMatthew G. Knepley PetscInt numBd, b, dim, pdim, fStart, fEnd; 319254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 320254c1ad2SMatthew G. Knepley 321254c1ad2SMatthew G. Knepley PetscFunctionBegin; 322254c1ad2SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 323254c1ad2SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 324254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 325254c1ad2SMatthew G. Knepley ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 326254c1ad2SMatthew G. Knepley if (Grad) { 327254c1ad2SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 328254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 329254c1ad2SMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 330254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 331254c1ad2SMatthew G. Knepley } 332254c1ad2SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 333254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 334254c1ad2SMatthew G. Knepley ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 335254c1ad2SMatthew G. Knepley for (b = 0; b < numBd; ++b) { 336*1475bf87SMatthew G. Knepley PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*); 3378cd7fcbbSMatthew G. Knepley DMLabel label; 3388cd7fcbbSMatthew G. Knepley const char *labelname; 339254c1ad2SMatthew G. Knepley const PetscInt *ids; 340254c1ad2SMatthew G. Knepley PetscInt numids, i; 341254c1ad2SMatthew G. Knepley void *ctx; 342254c1ad2SMatthew G. Knepley 3438cd7fcbbSMatthew G. Knepley ierr = DMPlexGetBoundary(dm, b, NULL, NULL, &labelname, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr); 3448cd7fcbbSMatthew G. Knepley ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 345254c1ad2SMatthew G. Knepley for (i = 0; i < numids; ++i) { 346254c1ad2SMatthew G. Knepley IS faceIS; 347254c1ad2SMatthew G. Knepley const PetscInt *faces; 348254c1ad2SMatthew G. Knepley PetscInt numFaces, f; 349254c1ad2SMatthew G. Knepley 350254c1ad2SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 351254c1ad2SMatthew G. Knepley if (!faceIS) continue; /* No points with that id on this process */ 352254c1ad2SMatthew G. Knepley ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 353254c1ad2SMatthew G. Knepley ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 354254c1ad2SMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 355254c1ad2SMatthew G. Knepley const PetscInt face = faces[f], *cells; 356254c1ad2SMatthew G. Knepley const FaceGeom *fg; 357254c1ad2SMatthew G. Knepley 358254c1ad2SMatthew G. Knepley if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 359254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 360254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 361254c1ad2SMatthew G. Knepley if (Grad) { 362254c1ad2SMatthew G. Knepley const CellGeom *cg; 363254c1ad2SMatthew G. Knepley const PetscScalar *cx, *cgrad; 364*1475bf87SMatthew G. Knepley PetscScalar *xG; 365*1475bf87SMatthew G. Knepley PetscReal dx[3]; 366254c1ad2SMatthew G. Knepley PetscInt d; 367254c1ad2SMatthew G. Knepley 368254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 369254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 370254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 371254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 372254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cg->centroid, fg->centroid, dx); 373254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DotD(dim, &cgrad[d*dim], dx); 374254c1ad2SMatthew G. Knepley ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr); 375254c1ad2SMatthew G. Knepley } else { 376254c1ad2SMatthew G. Knepley const PetscScalar *xI; 377254c1ad2SMatthew G. Knepley PetscScalar *xG; 378254c1ad2SMatthew G. Knepley 379254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 380254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 381254c1ad2SMatthew G. Knepley ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr); 382254c1ad2SMatthew G. Knepley } 383254c1ad2SMatthew G. Knepley } 384254c1ad2SMatthew G. Knepley ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 385254c1ad2SMatthew G. Knepley ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 386254c1ad2SMatthew G. Knepley } 387254c1ad2SMatthew G. Knepley } 388254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 389254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 390254c1ad2SMatthew G. Knepley if (Grad) { 391254c1ad2SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 392254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 393254c1ad2SMatthew G. Knepley } 394254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 395254c1ad2SMatthew G. Knepley } 396254c1ad2SMatthew G. Knepley 397254c1ad2SMatthew G. Knepley #undef __FUNCT__ 398254c1ad2SMatthew G. Knepley #define __FUNCT__ "TSComputeRHSFunction_DMPlex" 399254c1ad2SMatthew G. Knepley static PetscErrorCode TSComputeRHSFunction_DMPlex(TS ts, PetscReal time, Vec X, Vec F, void *ctx) 400254c1ad2SMatthew G. Knepley { 401254c1ad2SMatthew G. Knepley DM dm; 402254c1ad2SMatthew G. Knepley DMTS_Plex *dmplexts = (DMTS_Plex *) ctx; 403254c1ad2SMatthew G. Knepley void (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *) = dmplexts->riemann; 404254c1ad2SMatthew G. Knepley PetscFV fvm; 405254c1ad2SMatthew G. Knepley PetscLimiter lim; 406254c1ad2SMatthew G. Knepley Vec faceGeometry = dmplexts->facegeom; 407254c1ad2SMatthew G. Knepley Vec cellGeometry = dmplexts->cellgeom; 408254c1ad2SMatthew G. Knepley Vec Grad = NULL, locGrad, locX; 409254c1ad2SMatthew G. Knepley DM dmFace, dmCell; 4108cd7fcbbSMatthew G. Knepley DMLabel ghostLabel; 411254c1ad2SMatthew G. Knepley PetscCellGeometry fgeom, cgeom; 412254c1ad2SMatthew G. Knepley const PetscScalar *facegeom, *cellgeom, *x, *lgrad; 413254c1ad2SMatthew G. Knepley PetscScalar *grad, *f, *uL, *uR, *fluxL, *fluxR; 414254c1ad2SMatthew G. Knepley PetscReal *centroid, *normal, *vol, *cellPhi; 415254c1ad2SMatthew G. Knepley PetscBool computeGradients; 416254c1ad2SMatthew G. Knepley PetscInt Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface, cell, cStart, cEnd, cEndInterior; 417254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 418254c1ad2SMatthew G. Knepley 419254c1ad2SMatthew G. Knepley PetscFunctionBegin; 420254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(ts,TS_CLASSID,1); 421254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(X,VEC_CLASSID,3); 422254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(F,VEC_CLASSID,5); 423254c1ad2SMatthew G. Knepley ierr = TSGetDM(ts, &dm);CHKERRQ(ierr); 424254c1ad2SMatthew G. Knepley ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr); 425254c1ad2SMatthew G. Knepley ierr = VecZeroEntries(locX);CHKERRQ(ierr); 426254c1ad2SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 427254c1ad2SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr); 428254c1ad2SMatthew G. Knepley ierr = VecZeroEntries(F);CHKERRQ(ierr); 4296dbbd306SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 4306dbbd306SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 4316dbbd306SMatthew G. Knepley ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 432254c1ad2SMatthew G. Knepley ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr); 4336dbbd306SMatthew G. Knepley ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr); 434254c1ad2SMatthew G. Knepley ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); 435254c1ad2SMatthew G. Knepley if (computeGradients) { 436254c1ad2SMatthew G. Knepley ierr = DMGetGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr); 437254c1ad2SMatthew G. Knepley ierr = VecZeroEntries(Grad);CHKERRQ(ierr); 438254c1ad2SMatthew G. Knepley ierr = VecGetArray(Grad, &grad);CHKERRQ(ierr); 439254c1ad2SMatthew G. Knepley } 4406dbbd306SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr); 4416dbbd306SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 4426dbbd306SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 4436dbbd306SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 4446dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 4456dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 4466dbbd306SMatthew G. Knepley ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr); 447254c1ad2SMatthew G. Knepley /* Count faces and reconstruct gradients */ 4486dbbd306SMatthew G. Knepley for (face = fStart; face < fEnd; ++face) { 449254c1ad2SMatthew G. Knepley const PetscInt *cells; 450254c1ad2SMatthew G. Knepley const FaceGeom *fg; 451254c1ad2SMatthew G. Knepley const PetscScalar *cx[2]; 452254c1ad2SMatthew G. Knepley PetscScalar *cgrad[2]; 4538cd7fcbbSMatthew G. Knepley PetscBool boundary; 4548cd7fcbbSMatthew G. Knepley PetscInt ghost, c, pd, d; 4556dbbd306SMatthew G. Knepley 4566dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 4576dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 4586dbbd306SMatthew G. Knepley ++numFaces; 459254c1ad2SMatthew G. Knepley if (!computeGradients) continue; 4608cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 4618cd7fcbbSMatthew G. Knepley if (boundary) continue; 462254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 463254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 464254c1ad2SMatthew G. Knepley for (c = 0; c < 2; ++c) { 465254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr); 466254c1ad2SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cells[c], grad, &cgrad[c]);CHKERRQ(ierr); 467254c1ad2SMatthew G. Knepley } 468254c1ad2SMatthew G. Knepley for (pd = 0; pd < pdim; ++pd) { 469254c1ad2SMatthew G. Knepley PetscScalar delta = cx[1][pd] - cx[0][pd]; 470254c1ad2SMatthew G. Knepley 471254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) { 472254c1ad2SMatthew G. Knepley if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta; 473254c1ad2SMatthew G. Knepley if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta; 474254c1ad2SMatthew G. Knepley } 475254c1ad2SMatthew G. Knepley } 476254c1ad2SMatthew G. Knepley } 477254c1ad2SMatthew G. Knepley /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */ 478254c1ad2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 479254c1ad2SMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 480254c1ad2SMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 481254c1ad2SMatthew G. Knepley for (cell = computeGradients && lim ? cStart : cEnd; cell < cEndInterior; ++cell) { 482254c1ad2SMatthew G. Knepley const PetscInt *faces; 483254c1ad2SMatthew G. Knepley const PetscScalar *cx; 484254c1ad2SMatthew G. Knepley const CellGeom *cg; 485254c1ad2SMatthew G. Knepley PetscScalar *cgrad; 486254c1ad2SMatthew G. Knepley PetscInt coneSize, f, pd, d; 487254c1ad2SMatthew G. Knepley 488254c1ad2SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); 489254c1ad2SMatthew G. Knepley ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr); 490254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr); 491254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr); 492254c1ad2SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cell, grad, &cgrad);CHKERRQ(ierr); 493254c1ad2SMatthew G. Knepley if (!cgrad) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Supposedly ghost cell %d, but this should be impossible", cell); 494254c1ad2SMatthew G. Knepley /* Limiter will be minimum value over all neighbors */ 495254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) cellPhi[d] = PETSC_MAX_REAL; 496254c1ad2SMatthew G. Knepley for (f = 0; f < coneSize; ++f) { 497254c1ad2SMatthew G. Knepley const PetscScalar *ncx; 498254c1ad2SMatthew G. Knepley const CellGeom *ncg; 499254c1ad2SMatthew G. Knepley const PetscInt *fcells; 5008cd7fcbbSMatthew G. Knepley PetscInt face = faces[f], ncell, ghost; 501*1475bf87SMatthew G. Knepley PetscReal v[3]; 5028cd7fcbbSMatthew G. Knepley PetscBool boundary; 503254c1ad2SMatthew G. Knepley 504254c1ad2SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 5058cd7fcbbSMatthew G. Knepley ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr); 5068cd7fcbbSMatthew G. Knepley if ((ghost >= 0) || boundary) continue; 507254c1ad2SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr); 508254c1ad2SMatthew G. Knepley ncell = cell == fcells[0] ? fcells[1] : fcells[0]; 509254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr); 510254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr); 511254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cg->centroid, ncg->centroid, v); 512254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 513254c1ad2SMatthew G. Knepley /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */ 514*1475bf87SMatthew G. Knepley PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DotD(dim, &cgrad[d*dim], v); 515254c1ad2SMatthew G. Knepley 516254c1ad2SMatthew G. Knepley ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr); 517254c1ad2SMatthew G. Knepley cellPhi[d] = PetscMin(cellPhi[d], phi); 518254c1ad2SMatthew G. Knepley } 519254c1ad2SMatthew G. Knepley } 520254c1ad2SMatthew G. Knepley /* Apply limiter to gradient */ 521254c1ad2SMatthew G. Knepley for (pd = 0; pd < pdim; ++pd) 522254c1ad2SMatthew G. Knepley /* Scalar limiter applied to each component separately */ 523254c1ad2SMatthew G. Knepley for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd]; 524254c1ad2SMatthew G. Knepley } 525254c1ad2SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr); 526254c1ad2SMatthew G. Knepley ierr = DMPlexInsertBoundaryValuesFVM_Static(dm, fvm, time, locX, Grad, dmplexts);CHKERRQ(ierr); 527254c1ad2SMatthew G. Knepley if (computeGradients) { 528254c1ad2SMatthew G. Knepley ierr = VecRestoreArray(Grad, &grad);CHKERRQ(ierr); 529254c1ad2SMatthew G. Knepley ierr = DMGetLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr); 530254c1ad2SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 531254c1ad2SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 532254c1ad2SMatthew G. Knepley ierr = DMRestoreGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr); 533254c1ad2SMatthew G. Knepley ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 5346dbbd306SMatthew G. Knepley } 5356dbbd306SMatthew 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); 536254c1ad2SMatthew G. Knepley /* Read out values */ 5376dbbd306SMatthew G. Knepley for (face = fStart, iface = 0; face < fEnd; ++face) { 5386dbbd306SMatthew G. Knepley const PetscInt *cells; 5396dbbd306SMatthew G. Knepley const FaceGeom *fg; 5406dbbd306SMatthew G. Knepley const CellGeom *cgL, *cgR; 541254c1ad2SMatthew G. Knepley const PetscScalar *xL, *xR, *gL, *gR; 5426dbbd306SMatthew G. Knepley PetscInt ghost, d; 5436dbbd306SMatthew G. Knepley 5446dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 5456dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 5466dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 5476dbbd306SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 5486dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr); 5496dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr); 5506dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr); 5516dbbd306SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr); 552254c1ad2SMatthew G. Knepley if (computeGradients) { 553*1475bf87SMatthew G. Knepley PetscReal dxL[3], dxR[3]; 554254c1ad2SMatthew G. Knepley 555254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr); 556254c1ad2SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr); 557254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cgL->centroid, fg->centroid, dxL); 558254c1ad2SMatthew G. Knepley WaxpyD(dim, -1, cgR->centroid, fg->centroid, dxR); 559254c1ad2SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 560254c1ad2SMatthew G. Knepley uL[iface*pdim+d] = xL[d] + DotD(dim, &gL[d*dim], dxL); 561254c1ad2SMatthew G. Knepley uR[iface*pdim+d] = xR[d] + DotD(dim, &gR[d*dim], dxR); 562254c1ad2SMatthew G. Knepley } 563254c1ad2SMatthew G. Knepley } else { 5646dbbd306SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 5656dbbd306SMatthew G. Knepley uL[iface*pdim+d] = xL[d]; 5666dbbd306SMatthew G. Knepley uR[iface*pdim+d] = xR[d]; 5676dbbd306SMatthew G. Knepley } 568254c1ad2SMatthew G. Knepley } 5696dbbd306SMatthew G. Knepley for (d = 0; d < dim; ++d) { 5706dbbd306SMatthew G. Knepley centroid[iface*dim+d] = fg->centroid[d]; 5716dbbd306SMatthew G. Knepley normal[iface*dim+d] = fg->normal[d]; 5726dbbd306SMatthew G. Knepley } 5736dbbd306SMatthew G. Knepley vol[iface*2+0] = cgL->volume; 5746dbbd306SMatthew G. Knepley vol[iface*2+1] = cgR->volume; 5756dbbd306SMatthew G. Knepley ++iface; 5766dbbd306SMatthew G. Knepley } 577254c1ad2SMatthew G. Knepley if (computeGradients) { 578254c1ad2SMatthew G. Knepley ierr = VecRestoreArrayRead(locGrad,&lgrad);CHKERRQ(ierr); 579254c1ad2SMatthew G. Knepley ierr = DMRestoreLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr); 580254c1ad2SMatthew G. Knepley } 5816dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr); 5826dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 5836dbbd306SMatthew G. Knepley ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 5846dbbd306SMatthew G. Knepley fgeom.v0 = centroid; 5856dbbd306SMatthew G. Knepley fgeom.n = normal; 5866dbbd306SMatthew G. Knepley cgeom.vol = vol; 587254c1ad2SMatthew G. Knepley /* Riemann solve */ 588254c1ad2SMatthew G. Knepley ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, cgeom, uL, uR, riemann, fluxL, fluxR, dmplexts->rhsfunctionlocalctx);CHKERRQ(ierr); 589254c1ad2SMatthew G. Knepley /* Insert fluxes */ 5906dbbd306SMatthew G. Knepley ierr = VecGetArray(F, &f);CHKERRQ(ierr); 5916dbbd306SMatthew G. Knepley for (face = fStart, iface = 0; face < fEnd; ++face) { 5926dbbd306SMatthew G. Knepley const PetscInt *cells; 5936dbbd306SMatthew G. Knepley PetscScalar *fL, *fR; 5946dbbd306SMatthew G. Knepley PetscInt ghost, d; 5956dbbd306SMatthew G. Knepley 5966dbbd306SMatthew G. Knepley ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr); 5976dbbd306SMatthew G. Knepley if (ghost >= 0) continue; 5986dbbd306SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 5996dbbd306SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr); 6006dbbd306SMatthew G. Knepley ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr); 6016dbbd306SMatthew G. Knepley for (d = 0; d < pdim; ++d) { 6026dbbd306SMatthew G. Knepley if (fL) fL[d] -= fluxL[iface*pdim+d]; 6036dbbd306SMatthew G. Knepley if (fR) fR[d] += fluxR[iface*pdim+d]; 6046dbbd306SMatthew G. Knepley } 6056dbbd306SMatthew G. Knepley ++iface; 6066dbbd306SMatthew G. Knepley } 6076dbbd306SMatthew G. Knepley ierr = VecRestoreArray(F, &f);CHKERRQ(ierr); 6086dbbd306SMatthew G. Knepley ierr = PetscFree7(centroid,normal,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr); 609254c1ad2SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr); 610254c1ad2SMatthew G. Knepley PetscFunctionReturn(0); 611254c1ad2SMatthew G. Knepley } 612254c1ad2SMatthew G. Knepley 613254c1ad2SMatthew G. Knepley #undef __FUNCT__ 614254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetRHSFunctionLocal" 615254c1ad2SMatthew G. Knepley /*@C 616254c1ad2SMatthew G. Knepley DMPlexTSSetRHSFunctionLocal - set a local residual evaluation function 617254c1ad2SMatthew G. Knepley 618254c1ad2SMatthew G. Knepley Logically Collective 619254c1ad2SMatthew G. Knepley 620254c1ad2SMatthew G. Knepley Input Arguments: 621254c1ad2SMatthew G. Knepley + dm - DM to associate callback with 622254c1ad2SMatthew G. Knepley . riemann - Riemann solver 623254c1ad2SMatthew G. Knepley - ctx - optional context for Riemann solve 624254c1ad2SMatthew G. Knepley 625254c1ad2SMatthew G. Knepley Calling sequence for riemann: 626254c1ad2SMatthew G. Knepley 627254c1ad2SMatthew G. Knepley $ riemann(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx) 628254c1ad2SMatthew G. Knepley 629254c1ad2SMatthew G. Knepley + x - The coordinates at a point on the interface 630254c1ad2SMatthew G. Knepley . n - The normal vector to the interface 631254c1ad2SMatthew G. Knepley . uL - The state vector to the left of the interface 632254c1ad2SMatthew G. Knepley . uR - The state vector to the right of the interface 633254c1ad2SMatthew G. Knepley . flux - output array of flux through the interface 634254c1ad2SMatthew G. Knepley - ctx - optional user context 635254c1ad2SMatthew G. Knepley 636254c1ad2SMatthew G. Knepley Level: beginner 637254c1ad2SMatthew G. Knepley 638254c1ad2SMatthew G. Knepley .seealso: DMTSSetRHSFunctionLocal() 639254c1ad2SMatthew G. Knepley @*/ 640254c1ad2SMatthew G. Knepley PetscErrorCode DMPlexTSSetRHSFunctionLocal(DM dm, void (*riemann)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx), void *ctx) 641254c1ad2SMatthew G. Knepley { 642254c1ad2SMatthew G. Knepley DMTS dmts; 643254c1ad2SMatthew G. Knepley DMTS_Plex *dmplexts; 644254c1ad2SMatthew G. Knepley PetscFV fvm; 645254c1ad2SMatthew G. Knepley PetscInt Nf; 646254c1ad2SMatthew G. Knepley PetscBool computeGradients; 647254c1ad2SMatthew G. Knepley PetscErrorCode ierr; 648254c1ad2SMatthew G. Knepley 649254c1ad2SMatthew G. Knepley PetscFunctionBegin; 650254c1ad2SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 651254c1ad2SMatthew G. Knepley ierr = DMGetDMTSWrite(dm, &dmts);CHKERRQ(ierr); 652254c1ad2SMatthew G. Knepley ierr = DMPlexTSGetContext(dm, dmts, &dmplexts);CHKERRQ(ierr); 653254c1ad2SMatthew G. Knepley dmplexts->riemann = riemann; 654254c1ad2SMatthew G. Knepley dmplexts->rhsfunctionlocalctx = ctx; 655254c1ad2SMatthew G. Knepley ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 656254c1ad2SMatthew G. Knepley ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr); 657254c1ad2SMatthew G. Knepley ierr = DMPlexTSSetupGeometry(dm, fvm, dmplexts);CHKERRQ(ierr); 658254c1ad2SMatthew G. Knepley ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr); 659254c1ad2SMatthew G. Knepley if (computeGradients) {ierr = DMPlexTSSetupGradient(dm, fvm, dmplexts);CHKERRQ(ierr);} 660254c1ad2SMatthew G. Knepley ierr = DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMPlex, dmplexts);CHKERRQ(ierr); 6616dbbd306SMatthew G. Knepley PetscFunctionReturn(0); 6626dbbd306SMatthew G. Knepley } 663