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