xref: /petsc/src/ts/utils/dmplexts.c (revision 71ea845219e042c33e7cf8c77706e102bc3bc773)
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 
51475bf87SMatthew 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];}
61475bf87SMatthew 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;}
71475bf87SMatthew 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;}
81475bf87SMatthew 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"
67c510411aSMatthew G. Knepley /*@
68c510411aSMatthew G. Knepley   DMPlexTSGetGeometry - Return precomputed geometric data
69c510411aSMatthew G. Knepley 
70c510411aSMatthew G. Knepley   Input Parameter:
71c510411aSMatthew G. Knepley . dm - The DM
72c510411aSMatthew G. Knepley 
73c510411aSMatthew G. Knepley   Output Parameters:
74c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry
75c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry
76c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
77c510411aSMatthew G. Knepley 
78c510411aSMatthew G. Knepley   Level: developer
79c510411aSMatthew G. Knepley 
80c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal()
81c510411aSMatthew G. Knepley @*/
82a0ac79e7SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometry(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
83a0ac79e7SMatthew G. Knepley {
84a0ac79e7SMatthew G. Knepley   DMTS           dmts;
85a0ac79e7SMatthew G. Knepley   DMTS_Plex     *dmplexts;
86a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
87a0ac79e7SMatthew G. Knepley 
88a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
89a0ac79e7SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
90a0ac79e7SMatthew G. Knepley   ierr = DMPlexTSGetContext(dm, dmts, &dmplexts);CHKERRQ(ierr);
91a0ac79e7SMatthew G. Knepley   if (facegeom)  *facegeom  = dmplexts->facegeom;
92a0ac79e7SMatthew G. Knepley   if (cellgeom)  *cellgeom  = dmplexts->cellgeom;
93a0ac79e7SMatthew G. Knepley   if (minRadius) *minRadius = dmplexts->minradius;
94a0ac79e7SMatthew G. Knepley   PetscFunctionReturn(0);
95a0ac79e7SMatthew G. Knepley }
96a0ac79e7SMatthew G. Knepley 
97a0ac79e7SMatthew G. Knepley #undef __FUNCT__
98254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGeometry"
99254c1ad2SMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGeometry(DM dm, PetscFV fvm, DMTS_Plex *dmplexts)
100254c1ad2SMatthew G. Knepley {
101254c1ad2SMatthew G. Knepley   DM             dmFace, dmCell;
102254c1ad2SMatthew G. Knepley   DMLabel        ghostLabel;
103254c1ad2SMatthew G. Knepley   PetscSection   sectionFace, sectionCell;
104254c1ad2SMatthew G. Knepley   PetscSection   coordSection;
105254c1ad2SMatthew G. Knepley   Vec            coordinates;
106254c1ad2SMatthew G. Knepley   PetscReal      minradius;
107254c1ad2SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
108254c1ad2SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
109254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
110254c1ad2SMatthew G. Knepley 
111254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
112254c1ad2SMatthew G. Knepley   if (dmplexts->setupGeom) PetscFunctionReturn(0);
113254c1ad2SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
114254c1ad2SMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
115254c1ad2SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
116254c1ad2SMatthew G. Knepley   /* Make cell centroids and volumes */
117254c1ad2SMatthew G. Knepley   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
118254c1ad2SMatthew G. Knepley   ierr = DMSetCoordinateSection(dmCell, coordSection);CHKERRQ(ierr);
119254c1ad2SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
120254c1ad2SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
121254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
122254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
123254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
124254c1ad2SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, sizeof(CellGeom)/sizeof(PetscScalar));CHKERRQ(ierr);}
125254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
126254c1ad2SMatthew G. Knepley   ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr);
127254c1ad2SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
128254c1ad2SMatthew G. Knepley   ierr = DMCreateLocalVector(dmCell, &dmplexts->cellgeom);CHKERRQ(ierr);
129254c1ad2SMatthew G. Knepley   ierr = VecGetArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
130254c1ad2SMatthew G. Knepley   for (c = cStart; c < cEndInterior; ++c) {
131254c1ad2SMatthew G. Knepley     CellGeom *cg;
132254c1ad2SMatthew G. Knepley 
133254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
134254c1ad2SMatthew G. Knepley     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
135254c1ad2SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr);
136254c1ad2SMatthew G. Knepley   }
137254c1ad2SMatthew G. Knepley   /* Compute face normals and minimum cell radius */
138254c1ad2SMatthew G. Knepley   ierr = DMClone(dm, &dmFace);CHKERRQ(ierr);
139254c1ad2SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionFace);CHKERRQ(ierr);
140254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
141254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr);
142254c1ad2SMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, sizeof(FaceGeom)/sizeof(PetscScalar));CHKERRQ(ierr);}
143254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr);
144254c1ad2SMatthew G. Knepley   ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr);
145254c1ad2SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionFace);CHKERRQ(ierr);
146254c1ad2SMatthew G. Knepley   ierr = DMCreateLocalVector(dmFace, &dmplexts->facegeom);CHKERRQ(ierr);
147254c1ad2SMatthew G. Knepley   ierr = VecGetArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
148254c1ad2SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
149254c1ad2SMatthew G. Knepley   minradius = PETSC_MAX_REAL;
150254c1ad2SMatthew G. Knepley   for (f = fStart; f < fEnd; ++f) {
151254c1ad2SMatthew G. Knepley     FaceGeom *fg;
152254c1ad2SMatthew G. Knepley     PetscReal area;
153254c1ad2SMatthew G. Knepley     PetscInt  ghost, d;
154254c1ad2SMatthew G. Knepley 
155254c1ad2SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);
156254c1ad2SMatthew G. Knepley     if (ghost >= 0) continue;
157254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr);
158254c1ad2SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr);
159254c1ad2SMatthew G. Knepley     for (d = 0; d < dim; ++d) fg->normal[d] *= area;
160254c1ad2SMatthew G. Knepley     /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
161254c1ad2SMatthew G. Knepley     {
162254c1ad2SMatthew G. Knepley       CellGeom       *cL, *cR;
163254c1ad2SMatthew G. Knepley       const PetscInt *cells;
164254c1ad2SMatthew G. Knepley       PetscReal      *lcentroid, *rcentroid;
1651475bf87SMatthew G. Knepley       PetscReal       v[3];
166254c1ad2SMatthew G. Knepley 
167254c1ad2SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr);
168254c1ad2SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr);
169254c1ad2SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr);
170254c1ad2SMatthew G. Knepley       lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
171254c1ad2SMatthew G. Knepley       rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
172254c1ad2SMatthew G. Knepley       WaxpyD(dim, -1, lcentroid, rcentroid, v);
1731475bf87SMatthew G. Knepley       if (DotRealD(dim, fg->normal, v) < 0) {
174254c1ad2SMatthew G. Knepley         for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
175254c1ad2SMatthew G. Knepley       }
1761475bf87SMatthew G. Knepley       if (DotRealD(dim, fg->normal, v) <= 0) {
177809582abSMatthew G. Knepley         if (dim == 2) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g) v (%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) v[0], (double) v[1]);
178809582abSMatthew G. Knepley         if (dim == 3) SETERRQ7(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed, normal (%g,%g,%g) v (%g,%g,%g)", f, (double) fg->normal[0], (double) fg->normal[1], (double) fg->normal[2], (double) v[0], (double) v[1], (double) v[2]);
179254c1ad2SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
180254c1ad2SMatthew G. Knepley       }
181254c1ad2SMatthew G. Knepley       if (cells[0] < cEndInterior) {
182254c1ad2SMatthew G. Knepley         WaxpyD(dim, -1, fg->centroid, cL->centroid, v);
183254c1ad2SMatthew G. Knepley         minradius = PetscMin(minradius, NormD(dim, v));
184254c1ad2SMatthew G. Knepley       }
185254c1ad2SMatthew G. Knepley       if (cells[1] < cEndInterior) {
186254c1ad2SMatthew G. Knepley         WaxpyD(dim, -1, fg->centroid, cR->centroid, v);
187254c1ad2SMatthew G. Knepley         minradius = PetscMin(minradius, NormD(dim, v));
188254c1ad2SMatthew G. Knepley       }
189254c1ad2SMatthew G. Knepley     }
190254c1ad2SMatthew G. Knepley   }
191809582abSMatthew G. Knepley   ierr = MPI_Allreduce(&minradius, &dmplexts->minradius, 1, MPIU_REAL, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
192254c1ad2SMatthew G. Knepley   /* Compute centroids of ghost cells */
193254c1ad2SMatthew G. Knepley   for (c = cEndInterior; c < cEnd; ++c) {
194254c1ad2SMatthew G. Knepley     FaceGeom       *fg;
195254c1ad2SMatthew G. Knepley     const PetscInt *cone,    *support;
196254c1ad2SMatthew G. Knepley     PetscInt        coneSize, supportSize, s;
197254c1ad2SMatthew G. Knepley 
198254c1ad2SMatthew G. Knepley     ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr);
199254c1ad2SMatthew G. Knepley     if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
200254c1ad2SMatthew G. Knepley     ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr);
201254c1ad2SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr);
202254c1ad2SMatthew G. Knepley     if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 1", cone[0], supportSize);
203254c1ad2SMatthew G. Knepley     ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr);
204254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr);
205254c1ad2SMatthew G. Knepley     for (s = 0; s < 2; ++s) {
206254c1ad2SMatthew G. Knepley       /* Reflect ghost centroid across plane of face */
207254c1ad2SMatthew G. Knepley       if (support[s] == c) {
208254c1ad2SMatthew G. Knepley         const CellGeom *ci;
209254c1ad2SMatthew G. Knepley         CellGeom       *cg;
2101475bf87SMatthew G. Knepley         PetscReal       c2f[3], a;
211254c1ad2SMatthew G. Knepley 
212254c1ad2SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr);
213254c1ad2SMatthew G. Knepley         WaxpyD(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
2141475bf87SMatthew G. Knepley         a    = DotRealD(dim, c2f, fg->normal)/DotRealD(dim, fg->normal, fg->normal);
215254c1ad2SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr);
216254c1ad2SMatthew G. Knepley         WaxpyD(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
217254c1ad2SMatthew G. Knepley         cg->volume = ci->volume;
218254c1ad2SMatthew G. Knepley       }
219254c1ad2SMatthew G. Knepley     }
220254c1ad2SMatthew G. Knepley   }
221254c1ad2SMatthew G. Knepley   ierr = VecRestoreArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
222254c1ad2SMatthew G. Knepley   ierr = VecRestoreArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
223254c1ad2SMatthew G. Knepley   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
224254c1ad2SMatthew G. Knepley   ierr = DMDestroy(&dmFace);CHKERRQ(ierr);
225254c1ad2SMatthew G. Knepley   dmplexts->setupGeom = PETSC_TRUE;
226254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
227254c1ad2SMatthew G. Knepley }
228254c1ad2SMatthew G. Knepley 
229254c1ad2SMatthew G. Knepley #undef __FUNCT__
230c5148223SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction"
231c5148223SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
232254c1ad2SMatthew G. Knepley {
2338cd7fcbbSMatthew G. Knepley   DMLabel        ghostLabel;
234c5148223SMatthew G. Knepley   PetscScalar   *dx, *grad, **gref;
235c5148223SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
236254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
237254c1ad2SMatthew G. Knepley 
238254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
239254c1ad2SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
240254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
241c5148223SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
242254c1ad2SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr);
243c5148223SMatthew G. Knepley   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
244254c1ad2SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
245c5148223SMatthew G. Knepley   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
246254c1ad2SMatthew G. Knepley   for (c = cStart; c < cEndInterior; c++) {
247254c1ad2SMatthew G. Knepley     const PetscInt *faces;
248c5148223SMatthew G. Knepley     PetscInt        numFaces, usedFaces, f, d;
249254c1ad2SMatthew G. Knepley     const CellGeom *cg;
2508cd7fcbbSMatthew G. Knepley     PetscBool       boundary;
2518cd7fcbbSMatthew G. Knepley     PetscInt        ghost;
2528cd7fcbbSMatthew G. Knepley 
253254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
254c5148223SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr);
255c5148223SMatthew G. Knepley     ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr);
256c5148223SMatthew 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);
257c5148223SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
258c5148223SMatthew G. Knepley       const CellGeom *cg1;
259c5148223SMatthew G. Knepley       FaceGeom       *fg;
260254c1ad2SMatthew G. Knepley       const PetscInt *fcells;
261254c1ad2SMatthew G. Knepley       PetscInt        ncell, side;
262254c1ad2SMatthew G. Knepley 
263254c1ad2SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2648cd7fcbbSMatthew G. Knepley       ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2658cd7fcbbSMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
266254c1ad2SMatthew G. Knepley       ierr  = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
267254c1ad2SMatthew G. Knepley       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
268254c1ad2SMatthew G. Knepley       ncell = fcells[!side];    /* the neighbor */
269254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
270254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
271c5148223SMatthew G. Knepley       for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
272254c1ad2SMatthew G. Knepley       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
273254c1ad2SMatthew G. Knepley     }
274254c1ad2SMatthew G. Knepley     if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
275c5148223SMatthew G. Knepley     ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr);
276c5148223SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
277254c1ad2SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
2788cd7fcbbSMatthew G. Knepley       ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
2798cd7fcbbSMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
280c5148223SMatthew G. Knepley       for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
281c5148223SMatthew G. Knepley       ++usedFaces;
282254c1ad2SMatthew G. Knepley     }
283254c1ad2SMatthew G. Knepley   }
284c5148223SMatthew G. Knepley   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
285254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
286254c1ad2SMatthew G. Knepley }
287254c1ad2SMatthew G. Knepley 
288254c1ad2SMatthew G. Knepley #undef __FUNCT__
289254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGradient"
290254c1ad2SMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGradient(DM dm, PetscFV fvm, DMTS_Plex *dmplexts)
291254c1ad2SMatthew G. Knepley {
292254c1ad2SMatthew G. Knepley   DM             dmFace, dmCell;
293254c1ad2SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
294254c1ad2SMatthew G. Knepley   PetscSection   sectionGrad;
295254c1ad2SMatthew G. Knepley   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
296254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
297254c1ad2SMatthew G. Knepley 
298254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
299254c1ad2SMatthew G. Knepley   if (dmplexts->setupGrad) PetscFunctionReturn(0);
300254c1ad2SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
301254c1ad2SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
302254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
303254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
304254c1ad2SMatthew G. Knepley   /* Construct the interpolant corresponding to each face from the leat-square solution over the cell neighborhood */
305254c1ad2SMatthew G. Knepley   ierr = VecGetDM(dmplexts->facegeom, &dmFace);CHKERRQ(ierr);
306254c1ad2SMatthew G. Knepley   ierr = VecGetDM(dmplexts->cellgeom, &dmCell);CHKERRQ(ierr);
307254c1ad2SMatthew G. Knepley   ierr = VecGetArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
308254c1ad2SMatthew G. Knepley   ierr = VecGetArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
309c5148223SMatthew G. Knepley   ierr = BuildGradientReconstruction(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
310254c1ad2SMatthew G. Knepley   ierr = VecRestoreArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
311254c1ad2SMatthew G. Knepley   ierr = VecRestoreArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
312254c1ad2SMatthew G. Knepley   /* Create storage for gradients */
313254c1ad2SMatthew G. Knepley   ierr = DMClone(dm, &dmplexts->dmGrad);CHKERRQ(ierr);
314254c1ad2SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
315254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
316254c1ad2SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
317254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
318254c1ad2SMatthew G. Knepley   ierr = DMSetDefaultSection(dmplexts->dmGrad, sectionGrad);CHKERRQ(ierr);
319254c1ad2SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
320254c1ad2SMatthew G. Knepley   dmplexts->setupGrad = PETSC_TRUE;
321254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
322254c1ad2SMatthew G. Knepley }
323254c1ad2SMatthew G. Knepley 
324254c1ad2SMatthew G. Knepley #undef __FUNCT__
325254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM_Static"
326254c1ad2SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValuesFVM_Static(DM dm, PetscFV fvm, PetscReal time, Vec locX, Vec Grad, DMTS_Plex *dmplexts)
327254c1ad2SMatthew G. Knepley {
328254c1ad2SMatthew G. Knepley   Vec                faceGeometry = dmplexts->facegeom;
329254c1ad2SMatthew G. Knepley   Vec                cellGeometry = dmplexts->cellgeom;
330254c1ad2SMatthew G. Knepley   DM                 dmFace, dmCell;
331254c1ad2SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
332254c1ad2SMatthew G. Knepley   PetscScalar       *x, *fx;
333254c1ad2SMatthew G. Knepley   PetscInt           numBd, b, dim, pdim, fStart, fEnd;
334254c1ad2SMatthew G. Knepley   PetscErrorCode     ierr;
335254c1ad2SMatthew G. Knepley 
336254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
337254c1ad2SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
338254c1ad2SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
339254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
340254c1ad2SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
341254c1ad2SMatthew G. Knepley   if (Grad) {
342254c1ad2SMatthew G. Knepley     ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
343254c1ad2SMatthew G. Knepley     ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
344254c1ad2SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
345254c1ad2SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
346254c1ad2SMatthew G. Knepley   }
347254c1ad2SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
348254c1ad2SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
349254c1ad2SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
350254c1ad2SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
3511475bf87SMatthew G. Knepley     PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*);
3528cd7fcbbSMatthew G. Knepley     DMLabel          label;
3538cd7fcbbSMatthew G. Knepley     const char      *labelname;
354254c1ad2SMatthew G. Knepley     const PetscInt  *ids;
355254c1ad2SMatthew G. Knepley     PetscInt         numids, i;
356254c1ad2SMatthew G. Knepley     void            *ctx;
357254c1ad2SMatthew G. Knepley 
3588cd7fcbbSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, NULL, NULL, &labelname, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr);
3598cd7fcbbSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
360254c1ad2SMatthew G. Knepley     for (i = 0; i < numids; ++i) {
361254c1ad2SMatthew G. Knepley       IS              faceIS;
362254c1ad2SMatthew G. Knepley       const PetscInt *faces;
363254c1ad2SMatthew G. Knepley       PetscInt        numFaces, f;
364254c1ad2SMatthew G. Knepley 
365254c1ad2SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
366254c1ad2SMatthew G. Knepley       if (!faceIS) continue; /* No points with that id on this process */
367254c1ad2SMatthew G. Knepley       ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
368254c1ad2SMatthew G. Knepley       ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
369254c1ad2SMatthew G. Knepley       for (f = 0; f < numFaces; ++f) {
370254c1ad2SMatthew G. Knepley         const PetscInt     face = faces[f], *cells;
371254c1ad2SMatthew G. Knepley         const FaceGeom    *fg;
372254c1ad2SMatthew G. Knepley 
373254c1ad2SMatthew G. Knepley         if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
374254c1ad2SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
375254c1ad2SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
376254c1ad2SMatthew G. Knepley         if (Grad) {
377254c1ad2SMatthew G. Knepley           const CellGeom    *cg;
378254c1ad2SMatthew G. Knepley           const PetscScalar *cx, *cgrad;
3791475bf87SMatthew G. Knepley           PetscScalar       *xG;
3801475bf87SMatthew G. Knepley           PetscReal          dx[3];
381254c1ad2SMatthew G. Knepley           PetscInt           d;
382254c1ad2SMatthew G. Knepley 
383254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
384254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
385254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
386254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
387254c1ad2SMatthew G. Knepley           WaxpyD(dim, -1, cg->centroid, fg->centroid, dx);
388254c1ad2SMatthew G. Knepley           for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DotD(dim, &cgrad[d*dim], dx);
389254c1ad2SMatthew G. Knepley           ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
390254c1ad2SMatthew G. Knepley         } else {
391254c1ad2SMatthew G. Knepley           const PetscScalar *xI;
392254c1ad2SMatthew G. Knepley           PetscScalar       *xG;
393254c1ad2SMatthew G. Knepley 
394254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
395254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
396254c1ad2SMatthew G. Knepley           ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
397254c1ad2SMatthew G. Knepley         }
398254c1ad2SMatthew G. Knepley       }
399254c1ad2SMatthew G. Knepley       ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
400254c1ad2SMatthew G. Knepley       ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
401254c1ad2SMatthew G. Knepley     }
402254c1ad2SMatthew G. Knepley   }
403254c1ad2SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
404254c1ad2SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
405254c1ad2SMatthew G. Knepley   if (Grad) {
406254c1ad2SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
407254c1ad2SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
408254c1ad2SMatthew G. Knepley   }
409254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
410254c1ad2SMatthew G. Knepley }
411254c1ad2SMatthew G. Knepley 
412254c1ad2SMatthew G. Knepley #undef __FUNCT__
413254c1ad2SMatthew G. Knepley #define __FUNCT__ "TSComputeRHSFunction_DMPlex"
414254c1ad2SMatthew G. Knepley static PetscErrorCode TSComputeRHSFunction_DMPlex(TS ts, PetscReal time, Vec X, Vec F, void *ctx)
415254c1ad2SMatthew G. Knepley {
416254c1ad2SMatthew G. Knepley   DM                 dm;
417254c1ad2SMatthew G. Knepley   DMTS_Plex         *dmplexts = (DMTS_Plex *) ctx;
418254c1ad2SMatthew G. Knepley   void             (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *) = dmplexts->riemann;
419254c1ad2SMatthew G. Knepley   PetscFV            fvm;
420254c1ad2SMatthew G. Knepley   PetscLimiter       lim;
421254c1ad2SMatthew G. Knepley   Vec                faceGeometry = dmplexts->facegeom;
422254c1ad2SMatthew G. Knepley   Vec                cellGeometry = dmplexts->cellgeom;
423254c1ad2SMatthew G. Knepley   Vec                Grad = NULL, locGrad, locX;
424254c1ad2SMatthew G. Knepley   DM                 dmFace, dmCell;
4258cd7fcbbSMatthew G. Knepley   DMLabel            ghostLabel;
426254c1ad2SMatthew G. Knepley   PetscCellGeometry  fgeom, cgeom;
427254c1ad2SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
428254c1ad2SMatthew G. Knepley   PetscScalar       *grad, *f, *uL, *uR, *fluxL, *fluxR;
429254c1ad2SMatthew G. Knepley   PetscReal         *centroid, *normal, *vol, *cellPhi;
430254c1ad2SMatthew G. Knepley   PetscBool          computeGradients;
431254c1ad2SMatthew G. Knepley   PetscInt           Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface, cell, cStart, cEnd, cEndInterior;
432254c1ad2SMatthew G. Knepley   PetscErrorCode     ierr;
433254c1ad2SMatthew G. Knepley 
434254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
435254c1ad2SMatthew G. Knepley   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
436254c1ad2SMatthew G. Knepley   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
437254c1ad2SMatthew G. Knepley   PetscValidHeaderSpecific(F,VEC_CLASSID,5);
438254c1ad2SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
439254c1ad2SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
440254c1ad2SMatthew G. Knepley   ierr = VecZeroEntries(locX);CHKERRQ(ierr);
441254c1ad2SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
442254c1ad2SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
443254c1ad2SMatthew G. Knepley   ierr = VecZeroEntries(F);CHKERRQ(ierr);
4446dbbd306SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
4456dbbd306SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
4466dbbd306SMatthew G. Knepley   ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr);
447254c1ad2SMatthew G. Knepley   ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr);
4486dbbd306SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
449254c1ad2SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr);
450254c1ad2SMatthew G. Knepley   if (computeGradients) {
451254c1ad2SMatthew G. Knepley     ierr = DMGetGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr);
452254c1ad2SMatthew G. Knepley     ierr = VecZeroEntries(Grad);CHKERRQ(ierr);
453254c1ad2SMatthew G. Knepley     ierr = VecGetArray(Grad, &grad);CHKERRQ(ierr);
454254c1ad2SMatthew G. Knepley   }
4556dbbd306SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
4566dbbd306SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
4576dbbd306SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
4586dbbd306SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
4596dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
4606dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
4616dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr);
462254c1ad2SMatthew G. Knepley   /* Count faces and reconstruct gradients */
4636dbbd306SMatthew G. Knepley   for (face = fStart; face < fEnd; ++face) {
464254c1ad2SMatthew G. Knepley     const PetscInt    *cells;
465254c1ad2SMatthew G. Knepley     const FaceGeom    *fg;
466254c1ad2SMatthew G. Knepley     const PetscScalar *cx[2];
467254c1ad2SMatthew G. Knepley     PetscScalar       *cgrad[2];
4688cd7fcbbSMatthew G. Knepley     PetscBool          boundary;
4698cd7fcbbSMatthew G. Knepley     PetscInt           ghost, c, pd, d;
4706dbbd306SMatthew G. Knepley 
4716dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
4726dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
4736dbbd306SMatthew G. Knepley     ++numFaces;
474254c1ad2SMatthew G. Knepley     if (!computeGradients) continue;
4758cd7fcbbSMatthew G. Knepley     ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
4768cd7fcbbSMatthew G. Knepley     if (boundary) continue;
477254c1ad2SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
478254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
479254c1ad2SMatthew G. Knepley     for (c = 0; c < 2; ++c) {
480254c1ad2SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr);
481254c1ad2SMatthew G. Knepley       ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cells[c], grad, &cgrad[c]);CHKERRQ(ierr);
482254c1ad2SMatthew G. Knepley     }
483254c1ad2SMatthew G. Knepley     for (pd = 0; pd < pdim; ++pd) {
484254c1ad2SMatthew G. Knepley       PetscScalar delta = cx[1][pd] - cx[0][pd];
485254c1ad2SMatthew G. Knepley 
486254c1ad2SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
487254c1ad2SMatthew G. Knepley         if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta;
488254c1ad2SMatthew G. Knepley         if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta;
489254c1ad2SMatthew G. Knepley       }
490254c1ad2SMatthew G. Knepley     }
491254c1ad2SMatthew G. Knepley   }
492254c1ad2SMatthew G. Knepley   /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */
493254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
494254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
495254c1ad2SMatthew G. Knepley   ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
496254c1ad2SMatthew G. Knepley   for (cell = computeGradients && lim ? cStart : cEnd; cell < cEndInterior; ++cell) {
497254c1ad2SMatthew G. Knepley     const PetscInt    *faces;
498254c1ad2SMatthew G. Knepley     const PetscScalar *cx;
499254c1ad2SMatthew G. Knepley     const CellGeom    *cg;
500254c1ad2SMatthew G. Knepley     PetscScalar       *cgrad;
501254c1ad2SMatthew G. Knepley     PetscInt           coneSize, f, pd, d;
502254c1ad2SMatthew G. Knepley 
503254c1ad2SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
504254c1ad2SMatthew G. Knepley     ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
505254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr);
506254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr);
507254c1ad2SMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cell, grad, &cgrad);CHKERRQ(ierr);
508*71ea8452SMatthew G. Knepley     if (!cgrad) continue; /* Unowned overlap cell, we do not compute */
509254c1ad2SMatthew G. Knepley     /* Limiter will be minimum value over all neighbors */
510254c1ad2SMatthew G. Knepley     for (d = 0; d < pdim; ++d) cellPhi[d] = PETSC_MAX_REAL;
511254c1ad2SMatthew G. Knepley     for (f = 0; f < coneSize; ++f) {
512254c1ad2SMatthew G. Knepley       const PetscScalar *ncx;
513254c1ad2SMatthew G. Knepley       const CellGeom    *ncg;
514254c1ad2SMatthew G. Knepley       const PetscInt    *fcells;
5158cd7fcbbSMatthew G. Knepley       PetscInt           face = faces[f], ncell, ghost;
5161475bf87SMatthew G. Knepley       PetscReal          v[3];
5178cd7fcbbSMatthew G. Knepley       PetscBool          boundary;
518254c1ad2SMatthew G. Knepley 
519254c1ad2SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
5208cd7fcbbSMatthew G. Knepley       ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
5218cd7fcbbSMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
522254c1ad2SMatthew G. Knepley       ierr  = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr);
523254c1ad2SMatthew G. Knepley       ncell = cell == fcells[0] ? fcells[1] : fcells[0];
524254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr);
525254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr);
526254c1ad2SMatthew G. Knepley       WaxpyD(dim, -1, cg->centroid, ncg->centroid, v);
527254c1ad2SMatthew G. Knepley       for (d = 0; d < pdim; ++d) {
528254c1ad2SMatthew G. Knepley         /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */
5291475bf87SMatthew G. Knepley         PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DotD(dim, &cgrad[d*dim], v);
530254c1ad2SMatthew G. Knepley 
531254c1ad2SMatthew G. Knepley         ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr);
532254c1ad2SMatthew G. Knepley         cellPhi[d] = PetscMin(cellPhi[d], phi);
533254c1ad2SMatthew G. Knepley       }
534254c1ad2SMatthew G. Knepley     }
535254c1ad2SMatthew G. Knepley     /* Apply limiter to gradient */
536254c1ad2SMatthew G. Knepley     for (pd = 0; pd < pdim; ++pd)
537254c1ad2SMatthew G. Knepley       /* Scalar limiter applied to each component separately */
538254c1ad2SMatthew G. Knepley       for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd];
539254c1ad2SMatthew G. Knepley   }
540254c1ad2SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
541254c1ad2SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFVM_Static(dm, fvm, time, locX, Grad, dmplexts);CHKERRQ(ierr);
542254c1ad2SMatthew G. Knepley   if (computeGradients) {
543254c1ad2SMatthew G. Knepley     ierr = VecRestoreArray(Grad, &grad);CHKERRQ(ierr);
544254c1ad2SMatthew G. Knepley     ierr = DMGetLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr);
545254c1ad2SMatthew G. Knepley     ierr = DMGlobalToLocalBegin(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
546254c1ad2SMatthew G. Knepley     ierr = DMGlobalToLocalEnd(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
547254c1ad2SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr);
548254c1ad2SMatthew G. Knepley     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
5496dbbd306SMatthew G. Knepley   }
5506dbbd306SMatthew 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);
551254c1ad2SMatthew G. Knepley   /* Read out values */
5526dbbd306SMatthew G. Knepley   for (face = fStart, iface = 0; face < fEnd; ++face) {
5536dbbd306SMatthew G. Knepley     const PetscInt    *cells;
5546dbbd306SMatthew G. Knepley     const FaceGeom    *fg;
5556dbbd306SMatthew G. Knepley     const CellGeom    *cgL, *cgR;
556254c1ad2SMatthew G. Knepley     const PetscScalar *xL, *xR, *gL, *gR;
5576dbbd306SMatthew G. Knepley     PetscInt           ghost, d;
5586dbbd306SMatthew G. Knepley 
5596dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
5606dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
5616dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
5626dbbd306SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
5636dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr);
5646dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr);
5656dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr);
5666dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr);
567254c1ad2SMatthew G. Knepley     if (computeGradients) {
5681475bf87SMatthew G. Knepley       PetscReal dxL[3], dxR[3];
569254c1ad2SMatthew G. Knepley 
570254c1ad2SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr);
571254c1ad2SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr);
572254c1ad2SMatthew G. Knepley       WaxpyD(dim, -1, cgL->centroid, fg->centroid, dxL);
573254c1ad2SMatthew G. Knepley       WaxpyD(dim, -1, cgR->centroid, fg->centroid, dxR);
574254c1ad2SMatthew G. Knepley       for (d = 0; d < pdim; ++d) {
575254c1ad2SMatthew G. Knepley         uL[iface*pdim+d] = xL[d] + DotD(dim, &gL[d*dim], dxL);
576254c1ad2SMatthew G. Knepley         uR[iface*pdim+d] = xR[d] + DotD(dim, &gR[d*dim], dxR);
577254c1ad2SMatthew G. Knepley       }
578254c1ad2SMatthew G. Knepley     } else {
5796dbbd306SMatthew G. Knepley       for (d = 0; d < pdim; ++d) {
5806dbbd306SMatthew G. Knepley         uL[iface*pdim+d] = xL[d];
5816dbbd306SMatthew G. Knepley         uR[iface*pdim+d] = xR[d];
5826dbbd306SMatthew G. Knepley       }
583254c1ad2SMatthew G. Knepley     }
5846dbbd306SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
5856dbbd306SMatthew G. Knepley       centroid[iface*dim+d] = fg->centroid[d];
5866dbbd306SMatthew G. Knepley       normal[iface*dim+d]   = fg->normal[d];
5876dbbd306SMatthew G. Knepley     }
5886dbbd306SMatthew G. Knepley     vol[iface*2+0] = cgL->volume;
5896dbbd306SMatthew G. Knepley     vol[iface*2+1] = cgR->volume;
5906dbbd306SMatthew G. Knepley     ++iface;
5916dbbd306SMatthew G. Knepley   }
592254c1ad2SMatthew G. Knepley   if (computeGradients) {
593254c1ad2SMatthew G. Knepley     ierr = VecRestoreArrayRead(locGrad,&lgrad);CHKERRQ(ierr);
594254c1ad2SMatthew G. Knepley     ierr = DMRestoreLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr);
595254c1ad2SMatthew G. Knepley   }
5966dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr);
5976dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
5986dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
5996dbbd306SMatthew G. Knepley   fgeom.v0  = centroid;
6006dbbd306SMatthew G. Knepley   fgeom.n   = normal;
6016dbbd306SMatthew G. Knepley   cgeom.vol = vol;
602254c1ad2SMatthew G. Knepley   /* Riemann solve */
603254c1ad2SMatthew G. Knepley   ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, cgeom, uL, uR, riemann, fluxL, fluxR, dmplexts->rhsfunctionlocalctx);CHKERRQ(ierr);
604254c1ad2SMatthew G. Knepley   /* Insert fluxes */
6056dbbd306SMatthew G. Knepley   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
6066dbbd306SMatthew G. Knepley   for (face = fStart, iface = 0; face < fEnd; ++face) {
6076dbbd306SMatthew G. Knepley     const PetscInt *cells;
6086dbbd306SMatthew G. Knepley     PetscScalar    *fL, *fR;
6096dbbd306SMatthew G. Knepley     PetscInt        ghost, d;
6106dbbd306SMatthew G. Knepley 
6116dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
6126dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
6136dbbd306SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
6146dbbd306SMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr);
6156dbbd306SMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr);
6166dbbd306SMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
6176dbbd306SMatthew G. Knepley       if (fL) fL[d] -= fluxL[iface*pdim+d];
6186dbbd306SMatthew G. Knepley       if (fR) fR[d] += fluxR[iface*pdim+d];
6196dbbd306SMatthew G. Knepley     }
6206dbbd306SMatthew G. Knepley     ++iface;
6216dbbd306SMatthew G. Knepley   }
6226dbbd306SMatthew G. Knepley   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
6236dbbd306SMatthew G. Knepley   ierr = PetscFree7(centroid,normal,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr);
624254c1ad2SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
625254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
626254c1ad2SMatthew G. Knepley }
627254c1ad2SMatthew G. Knepley 
628254c1ad2SMatthew G. Knepley #undef __FUNCT__
629254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetRHSFunctionLocal"
630254c1ad2SMatthew G. Knepley /*@C
631254c1ad2SMatthew G. Knepley   DMPlexTSSetRHSFunctionLocal - set a local residual evaluation function
632254c1ad2SMatthew G. Knepley 
633254c1ad2SMatthew G. Knepley   Logically Collective
634254c1ad2SMatthew G. Knepley 
635254c1ad2SMatthew G. Knepley   Input Arguments:
636254c1ad2SMatthew G. Knepley + dm      - DM to associate callback with
637254c1ad2SMatthew G. Knepley . riemann - Riemann solver
638254c1ad2SMatthew G. Knepley - ctx     - optional context for Riemann solve
639254c1ad2SMatthew G. Knepley 
640254c1ad2SMatthew G. Knepley   Calling sequence for riemann:
641254c1ad2SMatthew G. Knepley 
642254c1ad2SMatthew G. Knepley $ riemann(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
643254c1ad2SMatthew G. Knepley 
644254c1ad2SMatthew G. Knepley + x    - The coordinates at a point on the interface
645254c1ad2SMatthew G. Knepley . n    - The normal vector to the interface
646254c1ad2SMatthew G. Knepley . uL   - The state vector to the left of the interface
647254c1ad2SMatthew G. Knepley . uR   - The state vector to the right of the interface
648254c1ad2SMatthew G. Knepley . flux - output array of flux through the interface
649254c1ad2SMatthew G. Knepley - ctx  - optional user context
650254c1ad2SMatthew G. Knepley 
651254c1ad2SMatthew G. Knepley   Level: beginner
652254c1ad2SMatthew G. Knepley 
653254c1ad2SMatthew G. Knepley .seealso: DMTSSetRHSFunctionLocal()
654254c1ad2SMatthew G. Knepley @*/
655254c1ad2SMatthew 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)
656254c1ad2SMatthew G. Knepley {
657254c1ad2SMatthew G. Knepley   DMTS           dmts;
658254c1ad2SMatthew G. Knepley   DMTS_Plex     *dmplexts;
659254c1ad2SMatthew G. Knepley   PetscFV        fvm;
660254c1ad2SMatthew G. Knepley   PetscInt       Nf;
661254c1ad2SMatthew G. Knepley   PetscBool      computeGradients;
662254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
663254c1ad2SMatthew G. Knepley 
664254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
665254c1ad2SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
666254c1ad2SMatthew G. Knepley   ierr = DMGetDMTSWrite(dm, &dmts);CHKERRQ(ierr);
667254c1ad2SMatthew G. Knepley   ierr = DMPlexTSGetContext(dm, dmts, &dmplexts);CHKERRQ(ierr);
668254c1ad2SMatthew G. Knepley   dmplexts->riemann             = riemann;
669254c1ad2SMatthew G. Knepley   dmplexts->rhsfunctionlocalctx = ctx;
670254c1ad2SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
671254c1ad2SMatthew G. Knepley   ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr);
672254c1ad2SMatthew G. Knepley   ierr = DMPlexTSSetupGeometry(dm, fvm, dmplexts);CHKERRQ(ierr);
673254c1ad2SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr);
674254c1ad2SMatthew G. Knepley   if (computeGradients) {ierr = DMPlexTSSetupGradient(dm, fvm, dmplexts);CHKERRQ(ierr);}
675254c1ad2SMatthew G. Knepley   ierr = DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMPlex, dmplexts);CHKERRQ(ierr);
6766dbbd306SMatthew G. Knepley   PetscFunctionReturn(0);
6776dbbd306SMatthew G. Knepley }
678