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