xref: /petsc/src/ts/utils/dmplexts.c (revision c49ccbb32729e871fff84de9e1f3d01f6ca029ee)
11faf85eaSMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2254c1ad2SMatthew G. Knepley #include <petsc-private/tsimpl.h>     /*I "petscts.h" I*/
3924a1b8fSMatthew G. Knepley #include <petscds.h>
46dbbd306SMatthew G. Knepley #include <petscfv.h>
56dbbd306SMatthew G. Knepley 
6254c1ad2SMatthew G. Knepley #undef __FUNCT__
7b6aca0f9SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometryFVM"
8c510411aSMatthew G. Knepley /*@
9b6aca0f9SMatthew G. Knepley   DMPlexTSGetGeometryFVM - Return precomputed geometric data
10c510411aSMatthew G. Knepley 
11c510411aSMatthew G. Knepley   Input Parameter:
12c510411aSMatthew G. Knepley . dm - The DM
13c510411aSMatthew G. Knepley 
14c510411aSMatthew G. Knepley   Output Parameters:
15c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry
16c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry
17c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
18c510411aSMatthew G. Knepley 
19c510411aSMatthew G. Knepley   Level: developer
20c510411aSMatthew G. Knepley 
21c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal()
22c510411aSMatthew G. Knepley @*/
23b6aca0f9SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
24a0ac79e7SMatthew G. Knepley {
25a0ac79e7SMatthew G. Knepley   DMTS           dmts;
261faf85eaSMatthew G. Knepley   PetscObject    obj;
27a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
28a0ac79e7SMatthew G. Knepley 
29a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
30924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
31a0ac79e7SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
32b6aca0f9SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
331faf85eaSMatthew G. Knepley   if (!obj) {
341faf85eaSMatthew G. Knepley     Vec cellgeom, facegeom;
351faf85eaSMatthew G. Knepley 
361faf85eaSMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr);
37b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
38b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
391faf85eaSMatthew G. Knepley     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
401faf85eaSMatthew G. Knepley     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
411faf85eaSMatthew G. Knepley   }
42b6aca0f9SMatthew G. Knepley   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
43b6aca0f9SMatthew G. Knepley   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
44113c68e6SMatthew G. Knepley   if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);}
45924a1b8fSMatthew G. Knepley   PetscFunctionReturn(0);
46924a1b8fSMatthew G. Knepley }
47924a1b8fSMatthew G. Knepley 
48924a1b8fSMatthew G. Knepley #undef __FUNCT__
49c5148223SMatthew G. Knepley #define __FUNCT__ "BuildGradientReconstruction"
50c5148223SMatthew G. Knepley static PetscErrorCode BuildGradientReconstruction(DM dm, PetscFV fvm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
51254c1ad2SMatthew G. Knepley {
528cd7fcbbSMatthew G. Knepley   DMLabel        ghostLabel;
53c5148223SMatthew G. Knepley   PetscScalar   *dx, *grad, **gref;
54c5148223SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, c, cEndInterior, maxNumFaces;
55254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
56254c1ad2SMatthew G. Knepley 
57254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
58c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
59254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
60c5148223SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
61254c1ad2SMatthew G. Knepley   ierr = DMPlexGetMaxSizes(dm, &maxNumFaces, NULL);CHKERRQ(ierr);
62c5148223SMatthew G. Knepley   ierr = PetscFVLeastSquaresSetMaxFaces(fvm, maxNumFaces);CHKERRQ(ierr);
63254c1ad2SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
64c5148223SMatthew G. Knepley   ierr = PetscMalloc3(maxNumFaces*dim, &dx, maxNumFaces*dim, &grad, maxNumFaces, &gref);CHKERRQ(ierr);
65254c1ad2SMatthew G. Knepley   for (c = cStart; c < cEndInterior; c++) {
66254c1ad2SMatthew G. Knepley     const PetscInt        *faces;
67c5148223SMatthew G. Knepley     PetscInt               numFaces, usedFaces, f, d;
681faf85eaSMatthew G. Knepley     const PetscFVCellGeom *cg;
698cd7fcbbSMatthew G. Knepley     PetscBool              boundary;
708cd7fcbbSMatthew G. Knepley     PetscInt               ghost;
718cd7fcbbSMatthew G. Knepley 
72254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, c, cgeom, &cg);CHKERRQ(ierr);
73c5148223SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, c, &numFaces);CHKERRQ(ierr);
74c5148223SMatthew G. Knepley     ierr = DMPlexGetCone(dm, c, &faces);CHKERRQ(ierr);
75c5148223SMatthew 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);
76c5148223SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
771faf85eaSMatthew G. Knepley       const PetscFVCellGeom *cg1;
781faf85eaSMatthew G. Knepley       PetscFVFaceGeom       *fg;
79254c1ad2SMatthew G. Knepley       const PetscInt        *fcells;
80254c1ad2SMatthew G. Knepley       PetscInt               ncell, side;
81254c1ad2SMatthew G. Knepley 
82254c1ad2SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
838cd7fcbbSMatthew G. Knepley       ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
848cd7fcbbSMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
85254c1ad2SMatthew G. Knepley       ierr  = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
86254c1ad2SMatthew G. Knepley       side  = (c != fcells[0]); /* c is on left=0 or right=1 of face */
87254c1ad2SMatthew G. Knepley       ncell = fcells[!side];    /* the neighbor */
88254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRef(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
89254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dmCell, ncell, cgeom, &cg1);CHKERRQ(ierr);
90c5148223SMatthew G. Knepley       for (d = 0; d < dim; ++d) dx[usedFaces*dim+d] = cg1->centroid[d] - cg->centroid[d];
91254c1ad2SMatthew G. Knepley       gref[usedFaces++] = fg->grad[side];  /* Gradient reconstruction term will go here */
92254c1ad2SMatthew G. Knepley     }
93254c1ad2SMatthew G. Knepley     if (!usedFaces) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Mesh contains isolated cell (no neighbors). Is it intentional?");
94c5148223SMatthew G. Knepley     ierr = PetscFVComputeGradient(fvm, usedFaces, dx, grad);CHKERRQ(ierr);
95c5148223SMatthew G. Knepley     for (f = 0, usedFaces = 0; f < numFaces; ++f) {
96254c1ad2SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, faces[f], &ghost);CHKERRQ(ierr);
978cd7fcbbSMatthew G. Knepley       ierr = DMPlexIsBoundaryPoint(dm, faces[f], &boundary);CHKERRQ(ierr);
988cd7fcbbSMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
99c5148223SMatthew G. Knepley       for (d = 0; d < dim; ++d) gref[usedFaces][d] = grad[usedFaces*dim+d];
100c5148223SMatthew G. Knepley       ++usedFaces;
101254c1ad2SMatthew G. Knepley     }
102254c1ad2SMatthew G. Knepley   }
103c5148223SMatthew G. Knepley   ierr = PetscFree3(dx, grad, gref);CHKERRQ(ierr);
104254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
105254c1ad2SMatthew G. Knepley }
106254c1ad2SMatthew G. Knepley 
107254c1ad2SMatthew G. Knepley #undef __FUNCT__
108*c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSSetupGradientFVM_Internal"
109*c49ccbb3SMatthew G. Knepley static PetscErrorCode DMPlexTSSetupGradientFVM_Internal(DM dm, PetscFV fvm, DM *dmGrad)
110254c1ad2SMatthew G. Knepley {
111*c49ccbb3SMatthew G. Knepley   DM             dmFace, dmCell;
112924a1b8fSMatthew G. Knepley   Vec            facegeom, cellgeom;
113254c1ad2SMatthew G. Knepley   PetscScalar   *fgeom, *cgeom;
114254c1ad2SMatthew G. Knepley   PetscSection   sectionGrad;
115254c1ad2SMatthew G. Knepley   PetscInt       dim, pdim, cStart, cEnd, cEndInterior, c;
116254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
117254c1ad2SMatthew G. Knepley 
118254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
119c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
120254c1ad2SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);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);
123*c49ccbb3SMatthew G. Knepley   /* Construct the interpolant corresponding to each face from the least-square solution over the cell neighborhood */
124b6aca0f9SMatthew G. Knepley   ierr = DMPlexTSGetGeometryFVM(dm, &facegeom, &cellgeom, NULL);CHKERRQ(ierr);
125924a1b8fSMatthew G. Knepley   ierr = VecGetDM(facegeom, &dmFace);CHKERRQ(ierr);
126924a1b8fSMatthew G. Knepley   ierr = VecGetDM(cellgeom, &dmCell);CHKERRQ(ierr);
127924a1b8fSMatthew G. Knepley   ierr = VecGetArray(facegeom, &fgeom);CHKERRQ(ierr);
128924a1b8fSMatthew G. Knepley   ierr = VecGetArray(cellgeom, &cgeom);CHKERRQ(ierr);
129c5148223SMatthew G. Knepley   ierr = BuildGradientReconstruction(dm, fvm, dmFace, fgeom, dmCell, cgeom);CHKERRQ(ierr);
130924a1b8fSMatthew G. Knepley   ierr = VecRestoreArray(facegeom, &fgeom);CHKERRQ(ierr);
131924a1b8fSMatthew G. Knepley   ierr = VecRestoreArray(cellgeom, &cgeom);CHKERRQ(ierr);
132254c1ad2SMatthew G. Knepley   /* Create storage for gradients */
133*c49ccbb3SMatthew G. Knepley   ierr = DMClone(dm, dmGrad);CHKERRQ(ierr);
134254c1ad2SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sectionGrad);CHKERRQ(ierr);
135254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetChart(sectionGrad, cStart, cEnd);CHKERRQ(ierr);
136254c1ad2SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {ierr = PetscSectionSetDof(sectionGrad, c, pdim*dim);CHKERRQ(ierr);}
137254c1ad2SMatthew G. Knepley   ierr = PetscSectionSetUp(sectionGrad);CHKERRQ(ierr);
138*c49ccbb3SMatthew G. Knepley   ierr = DMSetDefaultSection(*dmGrad, sectionGrad);CHKERRQ(ierr);
139254c1ad2SMatthew G. Knepley   ierr = PetscSectionDestroy(&sectionGrad);CHKERRQ(ierr);
140*c49ccbb3SMatthew G. Knepley   PetscFunctionReturn(0);
141*c49ccbb3SMatthew G. Knepley }
142*c49ccbb3SMatthew G. Knepley 
143*c49ccbb3SMatthew G. Knepley #undef __FUNCT__
144*c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM"
145*c49ccbb3SMatthew G. Knepley /*@C
146*c49ccbb3SMatthew G. Knepley   DMPlexTSGetGradientDM - Return gradient data layout
147*c49ccbb3SMatthew G. Knepley 
148*c49ccbb3SMatthew G. Knepley   Input Parameters:
149*c49ccbb3SMatthew G. Knepley + dm - The DM
150*c49ccbb3SMatthew G. Knepley - fv - The PetscFV
151*c49ccbb3SMatthew G. Knepley 
152*c49ccbb3SMatthew G. Knepley   Output Parameter:
153*c49ccbb3SMatthew G. Knepley . dmGrad - The layout for gradient values
154*c49ccbb3SMatthew G. Knepley 
155*c49ccbb3SMatthew G. Knepley   Level: developer
156*c49ccbb3SMatthew G. Knepley 
157*c49ccbb3SMatthew G. Knepley .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
158*c49ccbb3SMatthew G. Knepley @*/
159*c49ccbb3SMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
160*c49ccbb3SMatthew G. Knepley {
161*c49ccbb3SMatthew G. Knepley   DMTS           dmts;
162*c49ccbb3SMatthew G. Knepley   PetscObject    obj;
163*c49ccbb3SMatthew G. Knepley   PetscBool      computeGradients;
164*c49ccbb3SMatthew G. Knepley   PetscErrorCode ierr;
165*c49ccbb3SMatthew G. Knepley 
166*c49ccbb3SMatthew G. Knepley   PetscFunctionBegin;
167*c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
168*c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
169*c49ccbb3SMatthew G. Knepley   PetscValidPointer(dmGrad,3);
170*c49ccbb3SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
171*c49ccbb3SMatthew G. Knepley   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
172*c49ccbb3SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
173*c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
174*c49ccbb3SMatthew G. Knepley   if (!obj) {
175*c49ccbb3SMatthew G. Knepley     DM dmGrad;
176*c49ccbb3SMatthew G. Knepley 
177*c49ccbb3SMatthew G. Knepley     ierr = DMPlexTSSetupGradientFVM_Internal(dm, fv, &dmGrad);CHKERRQ(ierr);
178*c49ccbb3SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
179924a1b8fSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
180*c49ccbb3SMatthew G. Knepley   }
181*c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
182254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
183254c1ad2SMatthew G. Knepley }
184254c1ad2SMatthew G. Knepley 
185254c1ad2SMatthew G. Knepley #undef __FUNCT__
186254c1ad2SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFVM_Static"
187924a1b8fSMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValuesFVM_Static(DM dm, PetscFV fvm, PetscReal time, Vec locX, Vec Grad)
188254c1ad2SMatthew G. Knepley {
189924a1b8fSMatthew G. Knepley   Vec                faceGeometry, cellGeometry;
190924a1b8fSMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
191254c1ad2SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
192254c1ad2SMatthew G. Knepley   PetscScalar       *x, *fx;
193254c1ad2SMatthew G. Knepley   PetscInt           numBd, b, dim, pdim, fStart, fEnd;
194254c1ad2SMatthew G. Knepley   PetscErrorCode     ierr;
195254c1ad2SMatthew G. Knepley 
196254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
197c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
198b6aca0f9SMatthew G. Knepley   ierr = DMPlexTSGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
199*c49ccbb3SMatthew G. Knepley   ierr = DMPlexTSGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr);
200254c1ad2SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
201254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
202254c1ad2SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
203254c1ad2SMatthew G. Knepley   if (Grad) {
204254c1ad2SMatthew G. Knepley     ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
205254c1ad2SMatthew G. Knepley     ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
206254c1ad2SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
207254c1ad2SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
208254c1ad2SMatthew G. Knepley   }
209254c1ad2SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
210254c1ad2SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
211254c1ad2SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
212254c1ad2SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
2131475bf87SMatthew G. Knepley     PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*);
2148cd7fcbbSMatthew G. Knepley     DMLabel          label;
2158cd7fcbbSMatthew G. Knepley     const char      *labelname;
216254c1ad2SMatthew G. Knepley     const PetscInt  *ids;
217254c1ad2SMatthew G. Knepley     PetscInt         numids, i;
218254c1ad2SMatthew G. Knepley     void            *ctx;
219254c1ad2SMatthew G. Knepley 
2208cd7fcbbSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, NULL, NULL, &labelname, NULL, (void (**)()) &func, &numids, &ids, &ctx);CHKERRQ(ierr);
2218cd7fcbbSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
222254c1ad2SMatthew G. Knepley     for (i = 0; i < numids; ++i) {
223254c1ad2SMatthew G. Knepley       IS              faceIS;
224254c1ad2SMatthew G. Knepley       const PetscInt *faces;
225254c1ad2SMatthew G. Knepley       PetscInt        numFaces, f;
226254c1ad2SMatthew G. Knepley 
227254c1ad2SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
228254c1ad2SMatthew G. Knepley       if (!faceIS) continue; /* No points with that id on this process */
229254c1ad2SMatthew G. Knepley       ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
230254c1ad2SMatthew G. Knepley       ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
231254c1ad2SMatthew G. Knepley       for (f = 0; f < numFaces; ++f) {
232254c1ad2SMatthew G. Knepley         const PetscInt         face = faces[f], *cells;
2331faf85eaSMatthew G. Knepley         const PetscFVFaceGeom *fg;
234254c1ad2SMatthew G. Knepley 
235254c1ad2SMatthew G. Knepley         if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
236254c1ad2SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
237254c1ad2SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
238254c1ad2SMatthew G. Knepley         if (Grad) {
2391faf85eaSMatthew G. Knepley           const PetscFVCellGeom *cg;
240254c1ad2SMatthew G. Knepley           const PetscScalar     *cx, *cgrad;
2411475bf87SMatthew G. Knepley           PetscScalar           *xG;
2421475bf87SMatthew G. Knepley           PetscReal              dx[3];
243254c1ad2SMatthew G. Knepley           PetscInt               d;
244254c1ad2SMatthew G. Knepley 
245254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
246254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
247924a1b8fSMatthew G. Knepley           ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
248254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
2491faf85eaSMatthew G. Knepley           DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
2501faf85eaSMatthew G. Knepley           for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
251254c1ad2SMatthew G. Knepley           ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
252254c1ad2SMatthew G. Knepley         } else {
253254c1ad2SMatthew G. Knepley           const PetscScalar *xI;
254254c1ad2SMatthew G. Knepley           PetscScalar       *xG;
255254c1ad2SMatthew G. Knepley 
256254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
257254c1ad2SMatthew G. Knepley           ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
258254c1ad2SMatthew G. Knepley           ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
259254c1ad2SMatthew G. Knepley         }
260254c1ad2SMatthew G. Knepley       }
261254c1ad2SMatthew G. Knepley       ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
262254c1ad2SMatthew G. Knepley       ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
263254c1ad2SMatthew G. Knepley     }
264254c1ad2SMatthew G. Knepley   }
265254c1ad2SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
266254c1ad2SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
267254c1ad2SMatthew G. Knepley   if (Grad) {
268254c1ad2SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
269254c1ad2SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
270254c1ad2SMatthew G. Knepley   }
271254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
272254c1ad2SMatthew G. Knepley }
273254c1ad2SMatthew G. Knepley 
274254c1ad2SMatthew G. Knepley #undef __FUNCT__
275924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
276924a1b8fSMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *ctx)
277254c1ad2SMatthew G. Knepley {
278924a1b8fSMatthew G. Knepley   PetscDS            prob;
279924a1b8fSMatthew G. Knepley   void             (*riemann)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *);
280924a1b8fSMatthew G. Knepley   void              *rctx;
281254c1ad2SMatthew G. Knepley   PetscFV            fvm;
282254c1ad2SMatthew G. Knepley   PetscLimiter       lim;
283924a1b8fSMatthew G. Knepley   Vec                faceGeometry, cellGeometry;
284924a1b8fSMatthew G. Knepley   Vec                Grad = NULL, locGrad;
285924a1b8fSMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
2868cd7fcbbSMatthew G. Knepley   DMLabel            ghostLabel;
287b6aca0f9SMatthew G. Knepley   PetscFVFaceGeom   *fgeom;
288254c1ad2SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
289254c1ad2SMatthew G. Knepley   PetscScalar       *grad, *f, *uL, *uR, *fluxL, *fluxR;
290b6aca0f9SMatthew G. Knepley   PetscReal         *vol, *cellPhi;
291254c1ad2SMatthew G. Knepley   PetscBool          computeGradients;
292254c1ad2SMatthew G. Knepley   PetscInt           Nf, dim, pdim, fStart, fEnd, numFaces = 0, face, iface, cell, cStart, cEnd, cEndInterior;
293254c1ad2SMatthew G. Knepley   PetscErrorCode     ierr;
294254c1ad2SMatthew G. Knepley 
295254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
296924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
297924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(locX,VEC_CLASSID,3);
298254c1ad2SMatthew G. Knepley   PetscValidHeaderSpecific(F,VEC_CLASSID,5);
299b6aca0f9SMatthew G. Knepley   ierr = DMPlexTSGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
300924a1b8fSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
301924a1b8fSMatthew G. Knepley   ierr = PetscDSGetRiemannSolver(prob, 0, &riemann);CHKERRQ(ierr);
302924a1b8fSMatthew G. Knepley   ierr = PetscDSGetContext(prob, 0, &rctx);CHKERRQ(ierr);
303c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
3046dbbd306SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
3056dbbd306SMatthew G. Knepley   ierr = DMGetField(dm, 0, (PetscObject *) &fvm);CHKERRQ(ierr);
306*c49ccbb3SMatthew G. Knepley   ierr = DMPlexTSGetGradientDM(dm, fvm, &dmGrad);CHKERRQ(ierr);
307254c1ad2SMatthew G. Knepley   ierr = PetscFVGetLimiter(fvm, &lim);CHKERRQ(ierr);
3086dbbd306SMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
309254c1ad2SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fvm, &computeGradients);CHKERRQ(ierr);
310254c1ad2SMatthew G. Knepley   if (computeGradients) {
311924a1b8fSMatthew G. Knepley     ierr = DMGetGlobalVector(dmGrad, &Grad);CHKERRQ(ierr);
312254c1ad2SMatthew G. Knepley     ierr = VecZeroEntries(Grad);CHKERRQ(ierr);
313254c1ad2SMatthew G. Knepley     ierr = VecGetArray(Grad, &grad);CHKERRQ(ierr);
314254c1ad2SMatthew G. Knepley   }
3156dbbd306SMatthew G. Knepley   ierr = DMPlexGetLabel(dm, "ghost", &ghostLabel);CHKERRQ(ierr);
3166dbbd306SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
3176dbbd306SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
3186dbbd306SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
3196dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
3206dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
3216dbbd306SMatthew G. Knepley   ierr = VecGetArrayRead(locX, &x);CHKERRQ(ierr);
322254c1ad2SMatthew G. Knepley   /* Count faces and reconstruct gradients */
3236dbbd306SMatthew G. Knepley   for (face = fStart; face < fEnd; ++face) {
324254c1ad2SMatthew G. Knepley     const PetscInt        *cells;
3251faf85eaSMatthew G. Knepley     const PetscFVFaceGeom *fg;
326254c1ad2SMatthew G. Knepley     const PetscScalar     *cx[2];
327254c1ad2SMatthew G. Knepley     PetscScalar           *cgrad[2];
3288cd7fcbbSMatthew G. Knepley     PetscBool              boundary;
3298cd7fcbbSMatthew G. Knepley     PetscInt               ghost, c, pd, d;
3306dbbd306SMatthew G. Knepley 
3316dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
3326dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
3336dbbd306SMatthew G. Knepley     ++numFaces;
334254c1ad2SMatthew G. Knepley     if (!computeGradients) continue;
3358cd7fcbbSMatthew G. Knepley     ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
3368cd7fcbbSMatthew G. Knepley     if (boundary) continue;
337254c1ad2SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
338254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
339254c1ad2SMatthew G. Knepley     for (c = 0; c < 2; ++c) {
340254c1ad2SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dm, cells[c], x, &cx[c]);CHKERRQ(ierr);
341924a1b8fSMatthew G. Knepley       ierr = DMPlexPointGlobalRef(dmGrad, cells[c], grad, &cgrad[c]);CHKERRQ(ierr);
342254c1ad2SMatthew G. Knepley     }
343254c1ad2SMatthew G. Knepley     for (pd = 0; pd < pdim; ++pd) {
344254c1ad2SMatthew G. Knepley       PetscScalar delta = cx[1][pd] - cx[0][pd];
345254c1ad2SMatthew G. Knepley 
346254c1ad2SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
347254c1ad2SMatthew G. Knepley         if (cgrad[0]) cgrad[0][pd*dim+d] += fg->grad[0][d] * delta;
348254c1ad2SMatthew G. Knepley         if (cgrad[1]) cgrad[1][pd*dim+d] -= fg->grad[1][d] * delta;
349254c1ad2SMatthew G. Knepley       }
350254c1ad2SMatthew G. Knepley     }
351254c1ad2SMatthew G. Knepley   }
352254c1ad2SMatthew G. Knepley   /* Limit interior gradients (using cell-based loop because it generalizes better to vector limiters) */
353254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
354254c1ad2SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
355254c1ad2SMatthew G. Knepley   ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
356254c1ad2SMatthew G. Knepley   for (cell = computeGradients && lim ? cStart : cEnd; cell < cEndInterior; ++cell) {
357254c1ad2SMatthew G. Knepley     const PetscInt        *faces;
358254c1ad2SMatthew G. Knepley     const PetscScalar     *cx;
3591faf85eaSMatthew G. Knepley     const PetscFVCellGeom *cg;
360254c1ad2SMatthew G. Knepley     PetscScalar           *cgrad;
361254c1ad2SMatthew G. Knepley     PetscInt               coneSize, f, pd, d;
362254c1ad2SMatthew G. Knepley 
363254c1ad2SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr);
364254c1ad2SMatthew G. Knepley     ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
365254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cell, x, &cx);CHKERRQ(ierr);
366254c1ad2SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cell, cellgeom, &cg);CHKERRQ(ierr);
367924a1b8fSMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dmGrad, cell, grad, &cgrad);CHKERRQ(ierr);
36871ea8452SMatthew G. Knepley     if (!cgrad) continue; /* Unowned overlap cell, we do not compute */
369254c1ad2SMatthew G. Knepley     /* Limiter will be minimum value over all neighbors */
370254c1ad2SMatthew G. Knepley     for (d = 0; d < pdim; ++d) cellPhi[d] = PETSC_MAX_REAL;
371254c1ad2SMatthew G. Knepley     for (f = 0; f < coneSize; ++f) {
372254c1ad2SMatthew G. Knepley       const PetscScalar     *ncx;
3731faf85eaSMatthew G. Knepley       const PetscFVCellGeom *ncg;
374254c1ad2SMatthew G. Knepley       const PetscInt        *fcells;
3758cd7fcbbSMatthew G. Knepley       PetscInt               face = faces[f], ncell, ghost;
3761475bf87SMatthew G. Knepley       PetscReal              v[3];
3778cd7fcbbSMatthew G. Knepley       PetscBool              boundary;
378254c1ad2SMatthew G. Knepley 
379254c1ad2SMatthew G. Knepley       ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
3808cd7fcbbSMatthew G. Knepley       ierr = DMPlexIsBoundaryPoint(dm, face, &boundary);CHKERRQ(ierr);
3818cd7fcbbSMatthew G. Knepley       if ((ghost >= 0) || boundary) continue;
382254c1ad2SMatthew G. Knepley       ierr  = DMPlexGetSupport(dm, face, &fcells);CHKERRQ(ierr);
383254c1ad2SMatthew G. Knepley       ncell = cell == fcells[0] ? fcells[1] : fcells[0];
384254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dm, ncell, x, &ncx);CHKERRQ(ierr);
385254c1ad2SMatthew G. Knepley       ierr  = DMPlexPointLocalRead(dmCell, ncell, cellgeom, &ncg);CHKERRQ(ierr);
3861faf85eaSMatthew G. Knepley       DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, ncg->centroid, v);
387254c1ad2SMatthew G. Knepley       for (d = 0; d < pdim; ++d) {
388254c1ad2SMatthew G. Knepley         /* We use the symmetric slope limited form of Berger, Aftosmis, and Murman 2005 */
3891faf85eaSMatthew G. Knepley         PetscReal phi, flim = 0.5 * PetscRealPart(ncx[d] - cx[d]) / DMPlex_DotD_Internal(dim, &cgrad[d*dim], v);
390254c1ad2SMatthew G. Knepley 
391254c1ad2SMatthew G. Knepley         ierr = PetscLimiterLimit(lim, flim, &phi);CHKERRQ(ierr);
392254c1ad2SMatthew G. Knepley         cellPhi[d] = PetscMin(cellPhi[d], phi);
393254c1ad2SMatthew G. Knepley       }
394254c1ad2SMatthew G. Knepley     }
395254c1ad2SMatthew G. Knepley     /* Apply limiter to gradient */
396254c1ad2SMatthew G. Knepley     for (pd = 0; pd < pdim; ++pd)
397254c1ad2SMatthew G. Knepley       /* Scalar limiter applied to each component separately */
398254c1ad2SMatthew G. Knepley       for (d = 0; d < dim; ++d) cgrad[pd*dim+d] *= cellPhi[pd];
399254c1ad2SMatthew G. Knepley   }
400254c1ad2SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &cellPhi);CHKERRQ(ierr);
401924a1b8fSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValuesFVM_Static(dm, fvm, time, locX, Grad);CHKERRQ(ierr);
402254c1ad2SMatthew G. Knepley   if (computeGradients) {
403254c1ad2SMatthew G. Knepley     ierr = VecRestoreArray(Grad, &grad);CHKERRQ(ierr);
404924a1b8fSMatthew G. Knepley     ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
405924a1b8fSMatthew G. Knepley     ierr = DMGlobalToLocalBegin(dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
406924a1b8fSMatthew G. Knepley     ierr = DMGlobalToLocalEnd(dmGrad, Grad, INSERT_VALUES, locGrad);CHKERRQ(ierr);
407924a1b8fSMatthew G. Knepley     ierr = DMRestoreGlobalVector(dmGrad, &Grad);CHKERRQ(ierr);
408254c1ad2SMatthew G. Knepley     ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr);
4096dbbd306SMatthew G. Knepley   }
410b6aca0f9SMatthew G. Knepley   ierr = PetscMalloc6(numFaces,&fgeom,numFaces*2,&vol,numFaces*pdim,&uL,numFaces*pdim,&uR,numFaces*pdim,&fluxL,numFaces*pdim,&fluxR);CHKERRQ(ierr);
411254c1ad2SMatthew G. Knepley   /* Read out values */
4126dbbd306SMatthew G. Knepley   for (face = fStart, iface = 0; face < fEnd; ++face) {
4136dbbd306SMatthew G. Knepley     const PetscInt        *cells;
4141faf85eaSMatthew G. Knepley     const PetscFVFaceGeom *fg;
4151faf85eaSMatthew G. Knepley     const PetscFVCellGeom *cgL, *cgR;
416254c1ad2SMatthew G. Knepley     const PetscScalar     *xL, *xR, *gL, *gR;
4176dbbd306SMatthew G. Knepley     PetscInt               ghost, d;
4186dbbd306SMatthew G. Knepley 
4196dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
4206dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
4216dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
4226dbbd306SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
4236dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);CHKERRQ(ierr);
4246dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);CHKERRQ(ierr);
4256dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cells[0], x, &xL);CHKERRQ(ierr);
4266dbbd306SMatthew G. Knepley     ierr = DMPlexPointLocalRead(dm, cells[1], x, &xR);CHKERRQ(ierr);
427254c1ad2SMatthew G. Knepley     if (computeGradients) {
4281475bf87SMatthew G. Knepley       PetscReal dxL[3], dxR[3];
429254c1ad2SMatthew G. Knepley 
430924a1b8fSMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);CHKERRQ(ierr);
431924a1b8fSMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);CHKERRQ(ierr);
4321faf85eaSMatthew G. Knepley       DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL);
4331faf85eaSMatthew G. Knepley       DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR);
434254c1ad2SMatthew G. Knepley       for (d = 0; d < pdim; ++d) {
4351faf85eaSMatthew G. Knepley         uL[iface*pdim+d] = xL[d] + DMPlex_DotD_Internal(dim, &gL[d*dim], dxL);
4361faf85eaSMatthew G. Knepley         uR[iface*pdim+d] = xR[d] + DMPlex_DotD_Internal(dim, &gR[d*dim], dxR);
437254c1ad2SMatthew G. Knepley       }
438254c1ad2SMatthew G. Knepley     } else {
4396dbbd306SMatthew G. Knepley       for (d = 0; d < pdim; ++d) {
4406dbbd306SMatthew G. Knepley         uL[iface*pdim+d] = xL[d];
4416dbbd306SMatthew G. Knepley         uR[iface*pdim+d] = xR[d];
4426dbbd306SMatthew G. Knepley       }
443254c1ad2SMatthew G. Knepley     }
4446dbbd306SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
445b6aca0f9SMatthew G. Knepley       fgeom[iface].centroid[d] = fg->centroid[d];
446b6aca0f9SMatthew G. Knepley       fgeom[iface].normal[d]   = fg->normal[d];
4476dbbd306SMatthew G. Knepley     }
4486dbbd306SMatthew G. Knepley     vol[iface*2+0] = cgL->volume;
4496dbbd306SMatthew G. Knepley     vol[iface*2+1] = cgR->volume;
4506dbbd306SMatthew G. Knepley     ++iface;
4516dbbd306SMatthew G. Knepley   }
452254c1ad2SMatthew G. Knepley   if (computeGradients) {
453254c1ad2SMatthew G. Knepley     ierr = VecRestoreArrayRead(locGrad,&lgrad);CHKERRQ(ierr);
454924a1b8fSMatthew G. Knepley     ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr);
455254c1ad2SMatthew G. Knepley   }
4566dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(locX, &x);CHKERRQ(ierr);
4576dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
4586dbbd306SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
459254c1ad2SMatthew G. Knepley   /* Riemann solve */
460b6aca0f9SMatthew G. Knepley   ierr = PetscFVIntegrateRHSFunction(fvm, numFaces, Nf, &fvm, 0, fgeom, vol, uL, uR, riemann, fluxL, fluxR, rctx);CHKERRQ(ierr);
461254c1ad2SMatthew G. Knepley   /* Insert fluxes */
4626dbbd306SMatthew G. Knepley   ierr = VecGetArray(F, &f);CHKERRQ(ierr);
4636dbbd306SMatthew G. Knepley   for (face = fStart, iface = 0; face < fEnd; ++face) {
4646dbbd306SMatthew G. Knepley     const PetscInt *cells;
4656dbbd306SMatthew G. Knepley     PetscScalar    *fL, *fR;
4666dbbd306SMatthew G. Knepley     PetscInt        ghost, d;
4676dbbd306SMatthew G. Knepley 
4686dbbd306SMatthew G. Knepley     ierr = DMLabelGetValue(ghostLabel, face, &ghost);CHKERRQ(ierr);
4696dbbd306SMatthew G. Knepley     if (ghost >= 0) continue;
4706dbbd306SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
4716dbbd306SMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dm, cells[0], f, &fL);CHKERRQ(ierr);
4726dbbd306SMatthew G. Knepley     ierr = DMPlexPointGlobalRef(dm, cells[1], f, &fR);CHKERRQ(ierr);
4736dbbd306SMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
4746dbbd306SMatthew G. Knepley       if (fL) fL[d] -= fluxL[iface*pdim+d];
4756dbbd306SMatthew G. Knepley       if (fR) fR[d] += fluxR[iface*pdim+d];
4766dbbd306SMatthew G. Knepley     }
4776dbbd306SMatthew G. Knepley     ++iface;
4786dbbd306SMatthew G. Knepley   }
4796dbbd306SMatthew G. Knepley   ierr = VecRestoreArray(F, &f);CHKERRQ(ierr);
480b6aca0f9SMatthew G. Knepley   ierr = PetscFree6(fgeom,vol,uL,uR,fluxL,fluxR);CHKERRQ(ierr);
481254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
482254c1ad2SMatthew G. Knepley }
483254c1ad2SMatthew G. Knepley 
484254c1ad2SMatthew G. Knepley #undef __FUNCT__
48524cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
48624cdb843SMatthew G. Knepley /*@
48724cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
48824cdb843SMatthew G. Knepley 
48924cdb843SMatthew G. Knepley   Input Parameters:
49024cdb843SMatthew G. Knepley + dm - The mesh
49124cdb843SMatthew G. Knepley . t - The time
49224cdb843SMatthew G. Knepley . X  - Local solution
49324cdb843SMatthew G. Knepley . X_t - Local solution time derivative, or NULL
49424cdb843SMatthew G. Knepley - user - The user context
49524cdb843SMatthew G. Knepley 
49624cdb843SMatthew G. Knepley   Output Parameter:
49724cdb843SMatthew G. Knepley . F  - Local output vector
49824cdb843SMatthew G. Knepley 
49924cdb843SMatthew G. Knepley   Level: developer
50024cdb843SMatthew G. Knepley 
50124cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
50224cdb843SMatthew G. Knepley @*/
50324cdb843SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec X, Vec X_t, Vec F, void *user)
50424cdb843SMatthew G. Knepley {
50524cdb843SMatthew G. Knepley   PetscErrorCode ierr;
50624cdb843SMatthew G. Knepley 
50724cdb843SMatthew G. Knepley   PetscFunctionBegin;
50824cdb843SMatthew G. Knepley   ierr = DMPlexComputeResidualFEM_Internal(dm, X, X_t, F, user);CHKERRQ(ierr);
50924cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
51024cdb843SMatthew G. Knepley }
511