xref: /petsc/src/ts/utils/dmplexts.c (revision 8cd7fcbbde4ae99b6c1fac7f008f82577400d6fc)
1 #include <petscdmplex.h>          /*I "petscdmplex.h" I*/
2 #include <petsc-private/tsimpl.h>   /*I "petscts.h" I*/
3 #include <petscfv.h>
4 
5 /* TODO Move LS stuff to dtfv.c */
6 #include <petscblaslapack.h>
7 
8 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];}
9 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;}
10 PETSC_STATIC_INLINE PetscReal NormD(PetscInt dim, const PetscScalar *x) {return PetscSqrtReal(PetscAbsScalar(DotD(dim,x,x)));}
11 
12 typedef struct {
13   PetscBool setupGeom; /* Flag for geometry setup */
14   PetscBool setupGrad; /* Flag for gradient calculation setup */
15   Vec       facegeom;  /* FaceGeom struct for each face */
16   Vec       cellgeom;  /* CellGeom struct for each cell */
17   DM        dmGrad;    /* Layout for the gradient data */
18   PetscReal minradius; /* Minimum distance from centroid to face */
19   void    (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *);
20   void     *rhsfunctionlocalctx;
21 } DMTS_Plex;
22 
23 #undef __FUNCT__
24 #define __FUNCT__ "DMTSDestroy_Plex"
25 static PetscErrorCode DMTSDestroy_Plex(DMTS dmts)
26 {
27   DMTS_Plex     *dmplexts = (DMTS_Plex *) dmts->data;
28   PetscErrorCode ierr;
29 
30   PetscFunctionBegin;
31   ierr = DMDestroy(&dmplexts->dmGrad);CHKERRQ(ierr);
32   ierr = VecDestroy(&dmplexts->facegeom);CHKERRQ(ierr);
33   ierr = VecDestroy(&dmplexts->cellgeom);CHKERRQ(ierr);
34   ierr = PetscFree(dmts->data);CHKERRQ(ierr);
35   PetscFunctionReturn(0);
36 }
37 
38 #undef __FUNCT__
39 #define __FUNCT__ "DMTSDuplicate_Plex"
40 static PetscErrorCode DMTSDuplicate_Plex(DMTS olddmts, DMTS dmts)
41 {
42   PetscErrorCode ierr;
43 
44   PetscFunctionBegin;
45   ierr = PetscNewLog(dmts, (DMTS_Plex **) &dmts->data);CHKERRQ(ierr);
46   if (olddmts->data) {ierr = PetscMemcpy(dmts->data, olddmts->data, sizeof(DMTS_Plex));CHKERRQ(ierr);}
47   PetscFunctionReturn(0);
48 }
49 
50 #undef __FUNCT__
51 #define __FUNCT__ "DMPlexTSGetContext"
52 static PetscErrorCode DMPlexTSGetContext(DM dm, DMTS dmts, DMTS_Plex **dmplexts)
53 {
54   PetscErrorCode ierr;
55 
56   PetscFunctionBegin;
57   *dmplexts = NULL;
58   if (!dmts->data) {
59     ierr = PetscNewLog(dm, (DMTS_Plex **) &dmts->data);CHKERRQ(ierr);
60     dmts->ops->destroy   = DMTSDestroy_Plex;
61     dmts->ops->duplicate = DMTSDuplicate_Plex;
62   }
63   *dmplexts = (DMTS_Plex *) dmts->data;
64   PetscFunctionReturn(0);
65 }
66 
67 #undef __FUNCT__
68 #define __FUNCT__ "DMPlexTSGetGeometry"
69 PetscErrorCode DMPlexTSGetGeometry(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
70 {
71   DMTS           dmts;
72   DMTS_Plex     *dmplexts;
73   PetscErrorCode ierr;
74 
75   PetscFunctionBegin;
76   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
77   ierr = DMPlexTSGetContext(dm, dmts, &dmplexts);CHKERRQ(ierr);
78   if (facegeom)  *facegeom  = dmplexts->facegeom;
79   if (cellgeom)  *cellgeom  = dmplexts->cellgeom;
80   if (minRadius) *minRadius = dmplexts->minradius;
81   PetscFunctionReturn(0);
82 }
83 
84 #undef __FUNCT__
85 #define __FUNCT__ "DMPlexTSSetupGeometry"
86 static PetscErrorCode DMPlexTSSetupGeometry(DM dm, PetscFV fvm, DMTS_Plex *dmplexts)
87 {
88   DM             dmFace, dmCell;
89   DMLabel        ghostLabel;
90   PetscSection   sectionFace, sectionCell;
91   PetscSection   coordSection;
92   Vec            coordinates;
93   PetscReal      minradius;
94   PetscScalar   *fgeom, *cgeom;
95   PetscInt       dim, cStart, cEnd, cEndInterior, c, fStart, fEnd, f;
96   PetscErrorCode ierr;
97 
98   PetscFunctionBegin;
99   if (dmplexts->setupGeom) PetscFunctionReturn(0);
100   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
101   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
102   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
103   /* Make cell centroids and volumes */
104   ierr = DMClone(dm, &dmCell);CHKERRQ(ierr);
105   ierr = DMSetCoordinateSection(dmCell, coordSection);CHKERRQ(ierr);
106   ierr = DMSetCoordinatesLocal(dmCell, coordinates);CHKERRQ(ierr);
107   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionCell);CHKERRQ(ierr);
108   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
109   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
110   ierr = PetscSectionSetChart(sectionCell, cStart, cEnd);CHKERRQ(ierr);
111   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionCell, c, sizeof(CellGeom)/sizeof(PetscScalar));CHKERRQ(ierr);}
112   ierr = PetscSectionSetUp(sectionCell);CHKERRQ(ierr);
113   ierr = DMSetDefaultSection(dmCell, sectionCell);CHKERRQ(ierr);
114   ierr = PetscSectionDestroy(&sectionCell);CHKERRQ(ierr);
115   ierr = DMCreateLocalVector(dmCell, &dmplexts->cellgeom);CHKERRQ(ierr);
116   ierr = VecGetArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
117   for (c = cStart; c < cEndInterior; ++c) {
118     CellGeom *cg;
119 
120     ierr = DMPlexPointLocalRef(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
121     ierr = PetscMemzero(cg, sizeof(*cg));CHKERRQ(ierr);
122     ierr = DMPlexComputeCellGeometryFVM(dmCell, c, &cg->volume, cg->centroid, NULL);CHKERRQ(ierr);
123   }
124   /* Compute face normals and minimum cell radius */
125   ierr = DMClone(dm, &dmFace);CHKERRQ(ierr);
126   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionFace);CHKERRQ(ierr);
127   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
128   ierr = PetscSectionSetChart(sectionFace, fStart, fEnd);CHKERRQ(ierr);
129   for (f = fStart; f < fEnd; ++f) {ierr = PetscSectionSetDof(sectionFace, f, sizeof(FaceGeom)/sizeof(PetscScalar));CHKERRQ(ierr);}
130   ierr = PetscSectionSetUp(sectionFace);CHKERRQ(ierr);
131   ierr = DMSetDefaultSection(dmFace, sectionFace);CHKERRQ(ierr);
132   ierr = PetscSectionDestroy(&sectionFace);CHKERRQ(ierr);
133   ierr = DMCreateLocalVector(dmFace, &dmplexts->facegeom);CHKERRQ(ierr);
134   ierr = VecGetArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
135   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
136   minradius = PETSC_MAX_REAL;
137   for (f = fStart; f < fEnd; ++f) {
138     FaceGeom *fg;
139     PetscReal area;
140     PetscInt  ghost, d;
141 
142     ierr = DMLabelGetValue(ghostLabel, f, &ghost);CHKERRQ(ierr);
143     if (ghost >= 0) continue;
144     ierr = DMPlexPointLocalRef(dmFace, f, fgeom, &fg);CHKERRQ(ierr);
145     ierr = DMPlexComputeCellGeometryFVM(dm, f, &area, fg->centroid, fg->normal);CHKERRQ(ierr);
146     for (d = 0; d < dim; ++d) fg->normal[d] *= area;
147     /* Flip face orientation if necessary to match ordering in support, and Update minimum radius */
148     {
149       CellGeom       *cL, *cR;
150       const PetscInt *cells;
151       PetscReal      *lcentroid, *rcentroid;
152       PetscScalar     v[3];
153 
154       ierr = DMPlexGetSupport(dm, f, &cells);CHKERRQ(ierr);
155       ierr = DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cL);CHKERRQ(ierr);
156       ierr = DMPlexPointLocalRead(dmCell, cells[1], cgeom, &cR);CHKERRQ(ierr);
157       lcentroid = cells[0] >= cEndInterior ? fg->centroid : cL->centroid;
158       rcentroid = cells[1] >= cEndInterior ? fg->centroid : cR->centroid;
159       WaxpyD(dim, -1, lcentroid, rcentroid, v);
160       if (DotD(dim, fg->normal, v) < 0) {
161         for (d = 0; d < dim; ++d) fg->normal[d] = -fg->normal[d];
162       }
163       if (DotD(dim, fg->normal, v) <= 0) {
164         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]);
165         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]);
166         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Direction for face %d could not be fixed", f);
167       }
168       if (cells[0] < cEndInterior) {
169         WaxpyD(dim, -1, fg->centroid, cL->centroid, v);
170         minradius = PetscMin(minradius, NormD(dim, v));
171       }
172       if (cells[1] < cEndInterior) {
173         WaxpyD(dim, -1, fg->centroid, cR->centroid, v);
174         minradius = PetscMin(minradius, NormD(dim, v));
175       }
176     }
177   }
178   ierr = MPI_Allreduce(&minradius, &dmplexts->minradius, 1, MPIU_SCALAR, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
179   /* Compute centroids of ghost cells */
180   for (c = cEndInterior; c < cEnd; ++c) {
181     FaceGeom       *fg;
182     const PetscInt *cone,    *support;
183     PetscInt        coneSize, supportSize, s;
184 
185     ierr = DMPlexGetConeSize(dmCell, c, &coneSize);CHKERRQ(ierr);
186     if (coneSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Ghost cell %d has cone size %d != 1", c, coneSize);
187     ierr = DMPlexGetCone(dmCell, c, &cone);CHKERRQ(ierr);
188     ierr = DMPlexGetSupportSize(dmCell, cone[0], &supportSize);CHKERRQ(ierr);
189     if (supportSize != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d has support size %d != 1", cone[0], supportSize);
190     ierr = DMPlexGetSupport(dmCell, cone[0], &support);CHKERRQ(ierr);
191     ierr = DMPlexPointLocalRef(dmFace, cone[0], fgeom, &fg);CHKERRQ(ierr);
192     for (s = 0; s < 2; ++s) {
193       /* Reflect ghost centroid across plane of face */
194       if (support[s] == c) {
195         const CellGeom *ci;
196         CellGeom       *cg;
197         PetscScalar     c2f[3], a;
198 
199         ierr = DMPlexPointLocalRead(dmCell, support[(s+1)%2], cgeom, &ci);CHKERRQ(ierr);
200         WaxpyD(dim, -1, ci->centroid, fg->centroid, c2f); /* cell to face centroid */
201         a    = DotD(dim, c2f, fg->normal)/DotD(dim, fg->normal, fg->normal);
202         ierr = DMPlexPointLocalRef(dmCell, support[s], cgeom, &cg);CHKERRQ(ierr);
203         WaxpyD(dim, 2*a, fg->normal, ci->centroid, cg->centroid);
204         cg->volume = ci->volume;
205       }
206     }
207   }
208   ierr = VecRestoreArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
209   ierr = VecRestoreArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
210   ierr = DMDestroy(&dmCell);CHKERRQ(ierr);
211   ierr = DMDestroy(&dmFace);CHKERRQ(ierr);
212   dmplexts->setupGeom = PETSC_TRUE;
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "PseudoInverse"
218 /* Overwrites A. Can only handle full-rank problems with m>=n */
219 static PetscErrorCode PseudoInverse(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
220 {
221   PetscBool      debug = PETSC_FALSE;
222   PetscErrorCode ierr;
223   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
224   PetscScalar    *R,*Q,*Aback,Alpha;
225 
226   PetscFunctionBegin;
227   if (debug) {
228     ierr = PetscMalloc1(m*n,&Aback);CHKERRQ(ierr);
229     ierr = PetscMemcpy(Aback,A,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
230   }
231 
232   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
233   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
234   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
235   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
236   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
237   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
238   ierr = PetscFPTrapPop();CHKERRQ(ierr);
239   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
240   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
241 
242   /* Extract an explicit representation of Q */
243   Q    = Ainv;
244   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
245   K    = N;                     /* full rank */
246   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
247   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
248 
249   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
250   Alpha = 1.0;
251   ldb   = lda;
252   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
253   /* Ainv is Q, overwritten with inverse */
254 
255   if (debug) {                      /* Check that pseudo-inverse worked */
256     PetscScalar Beta = 0.0;
257     PetscInt    ldc;
258     K   = N;
259     ldc = N;
260     BLASgemm_("ConjugateTranspose","Normal",&N,&K,&M,&Alpha,Ainv,&lda,Aback,&ldb,&Beta,work,&ldc);
261     ierr = PetscScalarView(n*n,work,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
262     ierr = PetscFree(Aback);CHKERRQ(ierr);
263   }
264   PetscFunctionReturn(0);
265 }
266 
267 #undef __FUNCT__
268 #define __FUNCT__ "PseudoInverseGetWorkRequired"
269 static PetscErrorCode PseudoInverseGetWorkRequired(PetscInt maxFaces, PetscInt *work)
270 {
271   PetscInt m,n,nrhs,minwork;
272 
273   PetscFunctionBegin;
274   m       = maxFaces;
275   n       = 3; /* spatial dimension */
276   nrhs    = maxFaces;
277   minwork = 3*PetscMin(m,n) + PetscMax(2*PetscMin(m,n), PetscMax(PetscMax(m,n), nrhs)); /* required by LAPACK */
278   *work   = 5*minwork;          /* We can afford to be extra generous */
279   PetscFunctionReturn(0);
280 }
281 
282 #undef __FUNCT__
283 #define __FUNCT__ "PseudoInverseSVD"
284 /* Overwrites A. Can handle degenerate problems and m<n. */
285 static PetscErrorCode PseudoInverseSVD(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
286 {
287   PetscBool      debug = PETSC_FALSE;
288   PetscErrorCode ierr;
289   PetscInt       i,j,maxmn;
290   PetscBLASInt   M,N,nrhs,lda,ldb,irank,ldwork,info;
291   PetscScalar    rcond,*tmpwork,*Brhs,*Aback;
292 
293   PetscFunctionBegin;
294   if (debug) {
295     ierr = PetscMalloc1(m*n,&Aback);CHKERRQ(ierr);
296     ierr = PetscMemcpy(Aback,A,m*n*sizeof(PetscScalar));CHKERRQ(ierr);
297   }
298 
299   /* initialize to identity */
300   tmpwork = Ainv;
301   Brhs = work;
302   maxmn = PetscMax(m,n);
303   for (j=0; j<maxmn; j++) {
304     for (i=0; i<maxmn; i++) Brhs[i + j*maxmn] = 1.0*(i == j);
305   }
306 
307   ierr  = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
308   ierr  = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
309   nrhs  = M;
310   ierr  = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
311   ierr  = PetscBLASIntCast(maxmn,&ldb);CHKERRQ(ierr);
312   ierr  = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
313   rcond = -1;
314   ierr  = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
315   LAPACKgelss_(&M,&N,&nrhs,A,&lda,Brhs,&ldb,tau,&rcond,&irank,tmpwork,&ldwork,&info);
316   ierr = PetscFPTrapPop();CHKERRQ(ierr);
317   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGELSS error");
318   /* The following check should be turned into a diagnostic as soon as someone wants to do this intentionally */
319   if (irank < PetscMin(M,N)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Rank deficient least squares fit, indicates an isolated cell with two colinear points");
320 
321   /* Brhs shaped (M,nrhs) column-major coldim=mstride was overwritten by Ainv shaped (N,nrhs) column-major coldim=maxmn.
322    * Here we transpose to (N,nrhs) row-major rowdim=mstride. */
323   for (i=0; i<n; i++) {
324     for (j=0; j<nrhs; j++) Ainv[i*mstride+j] = Brhs[i + j*maxmn];
325   }
326   PetscFunctionReturn(0);
327 }
328 
329 #undef __FUNCT__
330 #define __FUNCT__ "BuildLeastSquares"
331 /* Build least squares gradient reconstruction operators */
332 static PetscErrorCode BuildLeastSquares(DM dm,PetscInt cEndInterior,DM dmFace,PetscScalar *fgeom,DM dmCell,PetscScalar *cgeom)
333 {
334   DMLabel        ghostLabel;
335   PetscScalar   *B,*Binv,*work,*tau,**gref;
336   PetscInt       dim, c,cStart,cEnd,maxNumFaces,worksize;
337   PetscBool      useSVD = PETSC_TRUE;
338   PetscErrorCode ierr;
339 
340   PetscFunctionBegin;
341   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
342   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
343   ierr = DMPlexGetMaxSizes(dm,&maxNumFaces,NULL);CHKERRQ(ierr);
344   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
345   ierr = PseudoInverseGetWorkRequired(maxNumFaces,&worksize);CHKERRQ(ierr);
346   ierr = PetscMalloc5(maxNumFaces*dim,&B,worksize,&Binv,worksize,&work,maxNumFaces,&tau,maxNumFaces,&gref);CHKERRQ(ierr);
347   for (c = cStart; c < cEndInterior; c++) {
348     const PetscInt *faces;
349     PetscInt       numFaces,usedFaces,f,i,j;
350     const CellGeom *cg;
351     PetscBool       boundary;
352     PetscInt        ghost;
353 
354     ierr = DMPlexGetConeSize(dm,c,&numFaces);CHKERRQ(ierr);
355     if (numFaces < dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cell %D has only %D faces, not enough for gradient reconstruction",c,numFaces);
356     ierr = DMPlexGetCone(dm,c,&faces);CHKERRQ(ierr);
357     ierr = DMPlexPointLocalRead(dmCell,c,cgeom,&cg);CHKERRQ(ierr);
358     for (f=0,usedFaces=0; f<numFaces; f++) {
359       const PetscInt *fcells;
360       PetscInt       ncell,side;
361       FaceGeom       *fg;
362       const CellGeom *cg1;
363 
364       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
365       ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
366       if ((ghost >= 0) || boundary) continue;
367       ierr  = DMPlexGetSupport(dm,faces[f],&fcells);CHKERRQ(ierr);
368       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
369       ncell = fcells[!side];   /* the neighbor */
370       ierr  = DMPlexPointLocalRef(dmFace,faces[f],fgeom,&fg);CHKERRQ(ierr);
371       ierr  = DMPlexPointLocalRead(dmCell,ncell,cgeom,&cg1);CHKERRQ(ierr);
372       for (j=0; j<dim; j++) B[j*numFaces+usedFaces] = cg1->centroid[j] - cg->centroid[j];
373       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
374     }
375     if (!usedFaces) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Mesh contains isolated cell (no neighbors). Is it intentional?");
376     /* Overwrites B with garbage, returns Binv in row-major format */
377     if (useSVD) {
378       ierr = PseudoInverseSVD(usedFaces,numFaces,dim,B,Binv,tau,worksize,work);CHKERRQ(ierr);
379     } else {
380       ierr = PseudoInverse(usedFaces,numFaces,dim,B,Binv,tau,worksize,work);CHKERRQ(ierr);
381     }
382     for (f=0,i=0; f<numFaces; f++) {
383       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
384       ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
385       if ((ghost >= 0) || boundary) continue;
386       for (j=0; j<dim; j++) gref[i][j] = Binv[j*numFaces+i];
387       i++;
388     }
389 
390     if (0) {
391       PetscReal grad[2] = {0,0};
392       for (f=0; f<numFaces; f++) {
393         const PetscInt *fcells;
394         const CellGeom *cg1;
395         const FaceGeom *fg;
396         ierr = DMPlexGetSupport(dm,faces[f],&fcells);CHKERRQ(ierr);
397         ierr = DMPlexPointLocalRead(dmFace,faces[f],fgeom,&fg);CHKERRQ(ierr);
398         for (i=0; i<2; i++) {
399           if (fcells[i] == c) continue;
400           ierr = DMPlexPointLocalRead(dmCell,fcells[i],cgeom,&cg1);CHKERRQ(ierr);
401           PetscScalar du = cg1->centroid[0] + 3*cg1->centroid[1] - (cg->centroid[0] + 3*cg->centroid[1]);
402           grad[0] += fg->grad[!i][0] * du;
403           grad[1] += fg->grad[!i][1] * du;
404         }
405       }
406       printf("cell[%d] grad (%g,%g)\n",c,grad[0],grad[1]);
407     }
408   }
409   ierr = PetscFree5(B,Binv,work,tau,gref);CHKERRQ(ierr);
410   PetscFunctionReturn(0);
411 }
412 
413 #undef __FUNCT__
414 #define __FUNCT__ "DMPlexTSSetupGradient"
415 static PetscErrorCode DMPlexTSSetupGradient(DM dm, PetscFV fvm, DMTS_Plex *dmplexts)
416 {
417   DM             dmFace, dmCell;
418   PetscScalar   *fgeom, *cgeom;
419   PetscSection   sectionGrad;
420   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
421   PetscErrorCode ierr;
422 
423   PetscFunctionBegin;
424   if (dmplexts->setupGrad) PetscFunctionReturn(0);
425   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
426   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
427   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
428   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
429   /* Construct the interpolant corresponding to each face from the leat-square solution over the cell neighborhood */
430   ierr = VecGetDM(dmplexts->facegeom, &dmFace);CHKERRQ(ierr);
431   ierr = VecGetDM(dmplexts->cellgeom, &dmCell);CHKERRQ(ierr);
432   ierr = VecGetArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
433   ierr = VecGetArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
434   ierr = BuildLeastSquares(dm, cEndInterior, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
435   ierr = VecRestoreArray(dmplexts->facegeom, &fgeom);CHKERRQ(ierr);
436   ierr = VecRestoreArray(dmplexts->cellgeom, &cgeom);CHKERRQ(ierr);
437   /* Create storage for gradients */
438   ierr = DMClone(dm, &dmplexts->dmGrad);CHKERRQ(ierr);
439   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
440   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
441   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
442   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
443   ierr = DMSetDefaultSection(dmplexts->dmGrad, sectionGrad);CHKERRQ(ierr);
444   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
445   dmplexts->setupGrad = PETSC_TRUE;
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM_Static"
451 static PetscErrorCode DMPlexInsertBoundaryValuesFVM_Static(DM dm, PetscFV fvm, PetscReal time, Vec locX, Vec Grad, DMTS_Plex *dmplexts)
452 {
453   Vec                faceGeometry = dmplexts->facegeom;
454   Vec                cellGeometry = dmplexts->cellgeom;
455   DM                 dmFace, dmCell;
456   const PetscScalar *facegeom, *cellgeom, *grad;
457   PetscScalar       *x, *fx;
458   PetscInt           numBd, b, dim, pdim, fStart, fEnd;
459   PetscErrorCode     ierr;
460 
461   PetscFunctionBegin;
462   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
463   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
464   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
465   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
466   if (Grad) {
467     ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
468     ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
469     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
470     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
471   }
472   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
473   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
474   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
475   for (b = 0; b < numBd; ++b) {
476     PetscErrorCode (*func)(PetscReal,const PetscScalar*,const PetscScalar*,const PetscScalar*,PetscScalar*,void*);
477     DMLabel          label;
478     const char      *labelname;
479     const PetscInt  *ids;
480     PetscInt         numids, i;
481     void            *ctx;
482 
483     ierr = DMPlexGetBoundary(dm, b, NULL, NULL, &labelname, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr);
484     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
485     for (i = 0; i < numids; ++i) {
486       IS              faceIS;
487       const PetscInt *faces;
488       PetscInt        numFaces, f;
489 
490       ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
491       if (!faceIS) continue; /* No points with that id on this process */
492       ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
493       ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
494       for (f = 0; f < numFaces; ++f) {
495         const PetscInt     face = faces[f], *cells;
496         const FaceGeom    *fg;
497 
498         if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
499         ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
500         ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
501         if (Grad) {
502           const CellGeom    *cg;
503           const PetscScalar *cx, *cgrad;
504           PetscScalar       *xG, dx[3];
505           PetscInt           d;
506 
507           ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
508           ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
509           ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
510           ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
511           WaxpyD(dim, -1, cg->centroid, fg->centroid, dx);
512           for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DotD(dim, &cgrad[d*dim], dx);
513           ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
514         } else {
515           const PetscScalar *xI;
516           PetscScalar       *xG;
517 
518           ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
519           ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
520           ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
521         }
522       }
523       ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
524       ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
525     }
526   }
527   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
528   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
529   if (Grad) {
530     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
531     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
532   }
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "TSComputeRHSFunction_DMPlex"
538 static PetscErrorCode TSComputeRHSFunction_DMPlex(TS ts, PetscReal time, Vec X, Vec F, void *ctx)
539 {
540   DM                 dm;
541   DMTS_Plex         *dmplexts = (DMTS_Plex *) ctx;
542   void             (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *) = dmplexts->riemann;
543   PetscFV            fvm;
544   PetscLimiter       lim;
545   Vec                faceGeometry = dmplexts->facegeom;
546   Vec                cellGeometry = dmplexts->cellgeom;
547   Vec                Grad = NULL, locGrad, locX;
548   DM                 dmFace, dmCell;
549   DMLabel            ghostLabel;
550   PetscCellGeometry  fgeom, cgeom;
551   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
552   PetscScalar       *grad, *f, *uL, *uR, *fluxL, *fluxR;
553   PetscReal         *centroid, *normal, *vol, *cellPhi;
554   PetscBool          computeGradients;
555   PetscInt           Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface, cell, cStart, cEnd, cEndInterior;
556   PetscErrorCode     ierr;
557 
558   PetscFunctionBegin;
559   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
560   PetscValidHeaderSpecific(X,VEC_CLASSID,3);
561   PetscValidHeaderSpecific(F,VEC_CLASSID,5);
562   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
563   ierr = DMGetLocalVector(dm, &locX);CHKERRQ(ierr);
564   ierr = VecZeroEntries(locX);CHKERRQ(ierr);
565   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
566   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX);CHKERRQ(ierr);
567   ierr = VecZeroEntries(F);CHKERRQ(ierr);
568   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
569   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
570   ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr);
571   ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr);
572   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
573   ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr);
574   if (computeGradients) {
575     ierr = DMGetGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr);
576     ierr = VecZeroEntries(Grad);CHKERRQ(ierr);
577     ierr = VecGetArray(Grad, &grad);CHKERRQ(ierr);
578   }
579   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
580   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
581   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
582   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
583   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
584   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
585   ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr);
586   /* Count faces and reconstruct gradients */
587   for (face = fStart; face < fEnd; ++face) {
588     const PetscInt    *cells;
589     const FaceGeom    *fg;
590     const PetscScalar *cx[2];
591     PetscScalar       *cgrad[2];
592     PetscBool          boundary;
593     PetscInt           ghost, c, pd, d;
594 
595     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
596     if (ghost >= 0) continue;
597     ++numFaces;
598     if (!computeGradients) continue;
599     ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
600     if (boundary) continue;
601     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
602     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
603     for (c = 0; c < 2; ++c) {
604       ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr);
605       ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cells[c], grad, &cgrad[c]);CHKERRQ(ierr);
606     }
607     for (pd = 0; pd < pdim; ++pd) {
608       PetscScalar delta = cx[1][pd] - cx[0][pd];
609 
610       for (d = 0; d < dim; ++d) {
611         if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta;
612         if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta;
613       }
614     }
615   }
616   /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */
617   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
618   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
619   ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
620   for (cell = computeGradients && lim ? cStart : cEnd; cell < cEndInterior; ++cell) {
621     const PetscInt    *faces;
622     const PetscScalar *cx;
623     const CellGeom    *cg;
624     PetscScalar       *cgrad;
625     PetscInt           coneSize, f, pd, d;
626 
627     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
628     ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
629     ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr);
630     ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr);
631     ierr = DMPlexPointGlobalRef(dmplexts->dmGrad, cell, grad, &cgrad);CHKERRQ(ierr);
632     if (!cgrad) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Supposedly ghost cell %d, but this should be impossible", cell);
633     /* Limiter will be minimum value over all neighbors */
634     for (d = 0; d < pdim; ++d) cellPhi[d] = PETSC_MAX_REAL;
635     for (f = 0; f < coneSize; ++f) {
636       const PetscScalar *ncx;
637       const CellGeom    *ncg;
638       const PetscInt    *fcells;
639       PetscInt           face = faces[f], ncell, ghost;
640       PetscScalar        v[3];
641       PetscBool          boundary;
642 
643       ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
644       ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
645       if ((ghost >= 0) || boundary) continue;
646       ierr  = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr);
647       ncell = cell == fcells[0] ? fcells[1] : fcells[0];
648       ierr  = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr);
649       ierr  = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr);
650       WaxpyD(dim, -1, cg->centroid, ncg->centroid, v);
651       for (d = 0; d < pdim; ++d) {
652         /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */
653         PetscScalar phi, flim = 0.5 * (ncx[d] - cx[d]) / DotD(dim, &cgrad[d*dim], v);
654 
655         ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr);
656         cellPhi[d] = PetscMin(cellPhi[d], phi);
657       }
658     }
659     /* Apply limiter to gradient */
660     for (pd = 0; pd < pdim; ++pd)
661       /* Scalar limiter applied to each component separately */
662       for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd];
663   }
664   ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
665   ierr = DMPlexInsertBoundaryValuesFVM_Static(dm, fvm, time, locX, Grad, dmplexts);CHKERRQ(ierr);
666   if (computeGradients) {
667     ierr = VecRestoreArray(Grad, &grad);CHKERRQ(ierr);
668     ierr = DMGetLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr);
669     ierr = DMGlobalToLocalBegin(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
670     ierr = DMGlobalToLocalEnd(dmplexts->dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
671     ierr = DMRestoreGlobalVector(dmplexts->dmGrad, &Grad);CHKERRQ(ierr);
672     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
673   }
674   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);
675   /* Read out values */
676   for (face = fStart, iface = 0; face < fEnd; ++face) {
677     const PetscInt    *cells;
678     const FaceGeom    *fg;
679     const CellGeom    *cgL, *cgR;
680     const PetscScalar *xL, *xR, *gL, *gR;
681     PetscInt           ghost, d;
682 
683     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
684     if (ghost >= 0) continue;
685     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
686     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
687     ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr);
688     ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr);
689     ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr);
690     ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr);
691     if (computeGradients) {
692       PetscScalar dxL[3], dxR[3];
693 
694       ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr);
695       ierr = DMPlexPointLocalRead(dmplexts->dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr);
696       WaxpyD(dim, -1, cgL->centroid, fg->centroid, dxL);
697       WaxpyD(dim, -1, cgR->centroid, fg->centroid, dxR);
698       for (d = 0; d < pdim; ++d) {
699         uL[iface*pdim+d] = xL[d] + DotD(dim, &gL[d*dim], dxL);
700         uR[iface*pdim+d] = xR[d] + DotD(dim, &gR[d*dim], dxR);
701       }
702     } else {
703       for (d = 0; d < pdim; ++d) {
704         uL[iface*pdim+d] = xL[d];
705         uR[iface*pdim+d] = xR[d];
706       }
707     }
708     for (d = 0; d < dim; ++d) {
709       centroid[iface*dim+d] = fg->centroid[d];
710       normal[iface*dim+d]   = fg->normal[d];
711     }
712     vol[iface*2+0] = cgL->volume;
713     vol[iface*2+1] = cgR->volume;
714     ++iface;
715   }
716   if (computeGradients) {
717     ierr = VecRestoreArrayRead(locGrad,&lgrad);CHKERRQ(ierr);
718     ierr = DMRestoreLocalVector(dmplexts->dmGrad, &locGrad);CHKERRQ(ierr);
719   }
720   ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr);
721   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
722   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
723   fgeom.v0  = centroid;
724   fgeom.n   = normal;
725   cgeom.vol = vol;
726   /* Riemann solve */
727   ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, cgeom, uL, uR, riemann, fluxL, fluxR, dmplexts->rhsfunctionlocalctx);CHKERRQ(ierr);
728   /* Insert fluxes */
729   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
730   for (face = fStart, iface = 0; face < fEnd; ++face) {
731     const PetscInt *cells;
732     PetscScalar    *fL, *fR;
733     PetscInt        ghost, d;
734 
735     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
736     if (ghost >= 0) continue;
737     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
738     ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr);
739     ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr);
740     for (d = 0; d < pdim; ++d) {
741       if (fL) fL[d] -= fluxL[iface*pdim+d];
742       if (fR) fR[d] += fluxR[iface*pdim+d];
743     }
744     ++iface;
745   }
746   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
747   ierr = PetscFree7(centroid,normal,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr);
748   ierr = DMRestoreLocalVector(dm, &locX);CHKERRQ(ierr);
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "DMPlexTSSetRHSFunctionLocal"
754 /*@C
755   DMPlexTSSetRHSFunctionLocal - set a local residual evaluation function
756 
757   Logically Collective
758 
759   Input Arguments:
760 + dm      - DM to associate callback with
761 . riemann - Riemann solver
762 - ctx     - optional context for Riemann solve
763 
764   Calling sequence for riemann:
765 
766 $ riemann(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
767 
768 + x    - The coordinates at a point on the interface
769 . n    - The normal vector to the interface
770 . uL   - The state vector to the left of the interface
771 . uR   - The state vector to the right of the interface
772 . flux - output array of flux through the interface
773 - ctx  - optional user context
774 
775   Level: beginner
776 
777 .seealso: DMTSSetRHSFunctionLocal()
778 @*/
779 PetscErrorCode DMPlexTSSetRHSFunctionLocal(DM dm, void (*riemann)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx), void *ctx)
780 {
781   DMTS           dmts;
782   DMTS_Plex     *dmplexts;
783   PetscFV        fvm;
784   PetscInt       Nf;
785   PetscBool      computeGradients;
786   PetscErrorCode ierr;
787 
788   PetscFunctionBegin;
789   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
790   ierr = DMGetDMTSWrite(dm, &dmts);CHKERRQ(ierr);
791   ierr = DMPlexTSGetContext(dm, dmts, &dmplexts);CHKERRQ(ierr);
792   dmplexts->riemann             = riemann;
793   dmplexts->rhsfunctionlocalctx = ctx;
794   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
795   ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr);
796   ierr = DMPlexTSSetupGeometry(dm, fvm, dmplexts);CHKERRQ(ierr);
797   ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr);
798   if (computeGradients) {ierr = DMPlexTSSetupGradient(dm, fvm, dmplexts);CHKERRQ(ierr);}
799   ierr = DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMPlex, dmplexts);CHKERRQ(ierr);
800   PetscFunctionReturn(0);
801 }
802