Actual source code: ex11.c
1: static char help[] = "Second Order TVD Finite Volume Example.\n";
We use a second order TVD finite volume method to evolve a system of PDEs. Our simple upwinded residual evaluation loops
over all mesh faces and uses a Riemann solver to produce the flux given the face geometry and cell values,
\begin{equation}
f_i = \mathrm{riemann}(\mathrm{phys}, p_\mathrm{centroid}, \hat n, x^L, x^R)
\end{equation}
and then update the cell values given the cell volume.
\begin{eqnarray}
f^L_i &-=& \frac{f_i}{vol^L} \\
f^R_i &+=& \frac{f_i}{vol^R}
\end{eqnarray}
As an example, we can consider the shallow water wave equation,
\begin{eqnarray}
h_t + \nabla\cdot \left( uh \right) &=& 0 \\
(uh)_t + \nabla\cdot \left( u\otimes uh + \frac{g h^2}{2} I \right) &=& 0
\end{eqnarray}
where $h$ is wave height, $u$ is wave velocity, and $g$ is the acceleration due to gravity.
A representative Riemann solver for the shallow water equations is given in the PhysicsRiemann_SW() function,
\begin{eqnarray}
f^{L,R}_h &=& uh^{L,R} \cdot \hat n \\
f^{L,R}_{uh} &=& \frac{f^{L,R}_h}{h^{L,R}} uh^{L,R} + g (h^{L,R})^2 \hat n \\
c^{L,R} &=& \sqrt{g h^{L,R}} \\
s &=& \max\left( \left|\frac{uh^L \cdot \hat n}{h^L}\right| + c^L, \left|\frac{uh^R \cdot \hat n}{h^R}\right| + c^R \right) \\
f_i &=& \frac{A_\mathrm{face}}{2} \left( f^L_i + f^R_i + s \left( x^L_i - x^R_i \right) \right)
\end{eqnarray}
where $c$ is the local gravity wave speed and $f_i$ is a Rusanov flux.
The more sophisticated residual evaluation in RHSFunctionLocal_LS() uses a least-squares fit to a quadratic polynomial
over a neighborhood of the given element.
The mesh is read in from an ExodusII file, usually generated by Cubit.
37: #include <petscdmplex.h>
38: #include <petscdmforest.h>
39: #include <petscds.h>
40: #include <petscts.h>
41: #include <petscsf.h>
43: #define DIM 2 /* Geometric dimension */
44: #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
46: static PetscFunctionList PhysicsList, PhysicsRiemannList_SW;
48: /* Represents continuum physical equations. */
49: typedef struct _n_Physics *Physics;
51: /* Physical model includes boundary conditions, initial conditions, and functionals of interest. It is
52: * discretization-independent, but its members depend on the scenario being solved. */
53: typedef struct _n_Model *Model;
55: /* 'User' implements a discretization of a continuous model. */
56: typedef struct _n_User *User;
57: typedef PetscErrorCode (*SolutionFunction)(Model,PetscReal,const PetscReal*,PetscScalar*,void*);
58: typedef PetscErrorCode (*SetUpBCFunction)(PetscDS,Physics);
59: typedef PetscErrorCode (*FunctionalFunction)(Model,PetscReal,const PetscReal*,const PetscScalar*,PetscReal*,void*);
60: typedef PetscErrorCode (*SetupFields)(Physics,PetscSection);
61: static PetscErrorCode ModelSolutionSetDefault(Model,SolutionFunction,void*);
62: static PetscErrorCode ModelFunctionalRegister(Model,const char*,PetscInt*,FunctionalFunction,void*);
63: static PetscErrorCode OutputVTK(DM,const char*,PetscViewer*);
65: struct FieldDescription {
66: const char *name;
67: PetscInt dof;
68: };
70: typedef struct _n_FunctionalLink *FunctionalLink;
71: struct _n_FunctionalLink {
72: char *name;
73: FunctionalFunction func;
74: void *ctx;
75: PetscInt offset;
76: FunctionalLink next;
77: };
79: struct _n_Physics {
80: PetscRiemannFunc riemann;
81: PetscInt dof; /* number of degrees of freedom per cell */
82: PetscReal maxspeed; /* kludge to pick initial time step, need to add monitoring and step control */
83: void *data;
84: PetscInt nfields;
85: const struct FieldDescription *field_desc;
86: };
88: struct _n_Model {
89: MPI_Comm comm; /* Does not do collective communicaton, but some error conditions can be collective */
90: Physics physics;
91: FunctionalLink functionalRegistry;
92: PetscInt maxComputed;
93: PetscInt numMonitored;
94: FunctionalLink *functionalMonitored;
95: PetscInt numCall;
96: FunctionalLink *functionalCall;
97: SolutionFunction solution;
98: SetUpBCFunction setupbc;
99: void *solutionctx;
100: PetscReal maxspeed; /* estimate of global maximum speed (for CFL calculation) */
101: PetscReal bounds[2*DIM];
102: DMBoundaryType bcs[3];
103: PetscErrorCode (*errorIndicator)(PetscInt, PetscReal, PetscInt, const PetscScalar[], const PetscScalar[], PetscReal *, void *);
104: void *errorCtx;
105: };
107: struct _n_User {
108: PetscInt numSplitFaces;
109: PetscInt vtkInterval; /* For monitor */
110: char outputBasename[PETSC_MAX_PATH_LEN]; /* Basename for output files */
111: PetscInt monitorStepOffset;
112: Model model;
113: PetscBool vtkmon;
114: };
116: PETSC_STATIC_INLINE PetscReal DotDIMReal(const PetscReal *x,const PetscReal *y)
117: {
118: PetscInt i;
119: PetscReal prod=0.0;
121: for (i=0; i<DIM; i++) prod += x[i]*y[i];
122: return prod;
123: }
124: PETSC_STATIC_INLINE PetscReal NormDIM(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(DotDIMReal(x,x))); }
126: PETSC_STATIC_INLINE PetscReal Dot2Real(const PetscReal *x,const PetscReal *y) { return x[0]*y[0] + x[1]*y[1];}
127: PETSC_STATIC_INLINE PetscReal Norm2Real(const PetscReal *x) { return PetscSqrtReal(PetscAbsReal(Dot2Real(x,x)));}
128: PETSC_STATIC_INLINE void Normalize2Real(PetscReal *x) { PetscReal a = 1./Norm2Real(x); x[0] *= a; x[1] *= a; }
129: PETSC_STATIC_INLINE void Waxpy2Real(PetscReal a,const PetscReal *x,const PetscReal *y,PetscReal *w) { w[0] = a*x[0] + y[0]; w[1] = a*x[1] + y[1]; }
130: PETSC_STATIC_INLINE void Scale2Real(PetscReal a,const PetscReal *x,PetscReal *y) { y[0] = a*x[0]; y[1] = a*x[1]; }
132: /******************* Advect ********************/
133: typedef enum {ADVECT_SOL_TILTED,ADVECT_SOL_BUMP,ADVECT_SOL_BUMP_CAVITY} AdvectSolType;
134: static const char *const AdvectSolTypes[] = {"TILTED","BUMP","BUMP_CAVITY","AdvectSolType","ADVECT_SOL_",0};
135: typedef enum {ADVECT_SOL_BUMP_CONE,ADVECT_SOL_BUMP_COS} AdvectSolBumpType;
136: static const char *const AdvectSolBumpTypes[] = {"CONE","COS","AdvectSolBumpType","ADVECT_SOL_BUMP_",0};
138: typedef struct {
139: PetscReal wind[DIM];
140: } Physics_Advect_Tilted;
141: typedef struct {
142: PetscReal center[DIM];
143: PetscReal radius;
144: AdvectSolBumpType type;
145: } Physics_Advect_Bump;
147: typedef struct {
148: PetscReal inflowState;
149: AdvectSolType soltype;
150: union {
151: Physics_Advect_Tilted tilted;
152: Physics_Advect_Bump bump;
153: } sol;
154: struct {
155: PetscInt Solution;
156: PetscInt Error;
157: } functional;
158: } Physics_Advect;
160: static const struct FieldDescription PhysicsFields_Advect[] = {{"U",1},{NULL,0}};
162: static PetscErrorCode PhysicsBoundary_Advect_Inflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
163: {
164: Physics phys = (Physics)ctx;
165: Physics_Advect *advect = (Physics_Advect*)phys->data;
168: xG[0] = advect->inflowState;
169: return(0);
170: }
172: static PetscErrorCode PhysicsBoundary_Advect_Outflow(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
173: {
175: xG[0] = xI[0];
176: return(0);
177: }
179: static void PhysicsRiemann_Advect(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
180: {
181: Physics_Advect *advect = (Physics_Advect*)phys->data;
182: PetscReal wind[DIM],wn;
184: switch (advect->soltype) {
185: case ADVECT_SOL_TILTED: {
186: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
187: wind[0] = tilted->wind[0];
188: wind[1] = tilted->wind[1];
189: } break;
190: case ADVECT_SOL_BUMP:
191: wind[0] = -qp[1];
192: wind[1] = qp[0];
193: break;
194: case ADVECT_SOL_BUMP_CAVITY:
195: {
196: PetscInt i;
197: PetscReal comp2[3] = {0.,0.,0.}, rad2;
199: rad2 = 0.;
200: for (i = 0; i < dim; i++) {
201: comp2[i] = qp[i] * qp[i];
202: rad2 += comp2[i];
203: }
205: wind[0] = -qp[1];
206: wind[1] = qp[0];
207: if (rad2 > 1.) {
208: PetscInt maxI = 0;
209: PetscReal maxComp2 = comp2[0];
211: for (i = 1; i < dim; i++) {
212: if (comp2[i] > maxComp2) {
213: maxI = i;
214: maxComp2 = comp2[i];
215: }
216: }
217: wind[maxI] = 0.;
218: }
219: }
220: break;
221: default:
222: {
223: PetscInt i;
224: for (i = 0; i < DIM; ++i) wind[i] = 0.0;
225: }
226: /* default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for solution type %s",AdvectSolBumpTypes[advect->soltype]); */
227: }
228: wn = Dot2Real(wind, n);
229: flux[0] = (wn > 0 ? xL[0] : xR[0]) * wn;
230: }
232: static PetscErrorCode PhysicsSolution_Advect(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
233: {
234: Physics phys = (Physics)ctx;
235: Physics_Advect *advect = (Physics_Advect*)phys->data;
238: switch (advect->soltype) {
239: case ADVECT_SOL_TILTED: {
240: PetscReal x0[DIM];
241: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
242: Waxpy2Real(-time,tilted->wind,x,x0);
243: if (x0[1] > 0) u[0] = 1.*x[0] + 3.*x[1];
244: else u[0] = advect->inflowState;
245: } break;
246: case ADVECT_SOL_BUMP_CAVITY:
247: case ADVECT_SOL_BUMP: {
248: Physics_Advect_Bump *bump = &advect->sol.bump;
249: PetscReal x0[DIM],v[DIM],r,cost,sint;
250: cost = PetscCosReal(time);
251: sint = PetscSinReal(time);
252: x0[0] = cost*x[0] + sint*x[1];
253: x0[1] = -sint*x[0] + cost*x[1];
254: Waxpy2Real(-1,bump->center,x0,v);
255: r = Norm2Real(v);
256: switch (bump->type) {
257: case ADVECT_SOL_BUMP_CONE:
258: u[0] = PetscMax(1 - r/bump->radius,0);
259: break;
260: case ADVECT_SOL_BUMP_COS:
261: u[0] = 0.5 + 0.5*PetscCosReal(PetscMin(r/bump->radius,1)*PETSC_PI);
262: break;
263: }
264: } break;
265: default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown solution type");
266: }
267: return(0);
268: }
270: static PetscErrorCode PhysicsFunctional_Advect(Model mod,PetscReal time,const PetscReal *x,const PetscScalar *y,PetscReal *f,void *ctx)
271: {
272: Physics phys = (Physics)ctx;
273: Physics_Advect *advect = (Physics_Advect*)phys->data;
274: PetscScalar yexact[1] = {0.0};
278: PhysicsSolution_Advect(mod,time,x,yexact,phys);
279: f[advect->functional.Solution] = PetscRealPart(y[0]);
280: f[advect->functional.Error] = PetscAbsScalar(y[0]-yexact[0]);
281: return(0);
282: }
284: static PetscErrorCode SetUpBC_Advect(PetscDS prob, Physics phys)
285: {
287: const PetscInt inflowids[] = {100,200,300},outflowids[] = {101};
290: /* Register "canned" boundary conditions and defaults for where to apply. */
291: PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "inflow", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Inflow, NULL, ALEN(inflowids), inflowids, phys);
292: PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "outflow", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Advect_Outflow, NULL, ALEN(outflowids), outflowids, phys);
293: return(0);
294: }
296: static PetscErrorCode PhysicsCreate_Advect(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
297: {
298: Physics_Advect *advect;
302: phys->field_desc = PhysicsFields_Advect;
303: phys->riemann = (PetscRiemannFunc)PhysicsRiemann_Advect;
304: PetscNew(&advect);
305: phys->data = advect;
306: mod->setupbc = SetUpBC_Advect;
308: PetscOptionsHead(PetscOptionsObject,"Advect options");
309: {
310: PetscInt two = 2,dof = 1;
311: advect->soltype = ADVECT_SOL_TILTED;
312: PetscOptionsEnum("-advect_sol_type","solution type","",AdvectSolTypes,(PetscEnum)advect->soltype,(PetscEnum*)&advect->soltype,NULL);
313: switch (advect->soltype) {
314: case ADVECT_SOL_TILTED: {
315: Physics_Advect_Tilted *tilted = &advect->sol.tilted;
316: two = 2;
317: tilted->wind[0] = 0.0;
318: tilted->wind[1] = 1.0;
319: PetscOptionsRealArray("-advect_tilted_wind","background wind vx,vy","",tilted->wind,&two,NULL);
320: advect->inflowState = -2.0;
321: PetscOptionsRealArray("-advect_tilted_inflow","Inflow state","",&advect->inflowState,&dof,NULL);
322: phys->maxspeed = Norm2Real(tilted->wind);
323: } break;
324: case ADVECT_SOL_BUMP_CAVITY:
325: case ADVECT_SOL_BUMP: {
326: Physics_Advect_Bump *bump = &advect->sol.bump;
327: two = 2;
328: bump->center[0] = 2.;
329: bump->center[1] = 0.;
330: PetscOptionsRealArray("-advect_bump_center","location of center of bump x,y","",bump->center,&two,NULL);
331: bump->radius = 0.9;
332: PetscOptionsReal("-advect_bump_radius","radius of bump","",bump->radius,&bump->radius,NULL);
333: bump->type = ADVECT_SOL_BUMP_CONE;
334: PetscOptionsEnum("-advect_bump_type","type of bump","",AdvectSolBumpTypes,(PetscEnum)bump->type,(PetscEnum*)&bump->type,NULL);
335: phys->maxspeed = 3.; /* radius of mesh, kludge */
336: } break;
337: }
338: }
339: PetscOptionsTail();
340: /* Initial/transient solution with default boundary conditions */
341: ModelSolutionSetDefault(mod,PhysicsSolution_Advect,phys);
342: /* Register "canned" functionals */
343: ModelFunctionalRegister(mod,"Solution",&advect->functional.Solution,PhysicsFunctional_Advect,phys);
344: ModelFunctionalRegister(mod,"Error",&advect->functional.Error,PhysicsFunctional_Advect,phys);
345: mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
346: return(0);
347: }
349: /******************* Shallow Water ********************/
350: typedef struct {
351: PetscReal gravity;
352: PetscReal boundaryHeight;
353: struct {
354: PetscInt Height;
355: PetscInt Speed;
356: PetscInt Energy;
357: } functional;
358: } Physics_SW;
359: typedef struct {
360: PetscReal h;
361: PetscReal uh[DIM];
362: } SWNode;
363: typedef union {
364: SWNode swnode;
365: PetscReal vals[DIM+1];
366: } SWNodeUnion;
368: static const struct FieldDescription PhysicsFields_SW[] = {{"Height",1},{"Momentum",DIM},{NULL,0}};
370: /*
371: * h_t + div(uh) = 0
372: * (uh)_t + div (u\otimes uh + g h^2 / 2 I) = 0
373: *
374: * */
375: static PetscErrorCode SWFlux(Physics phys,const PetscReal *n,const SWNode *x,SWNode *f)
376: {
377: Physics_SW *sw = (Physics_SW*)phys->data;
378: PetscReal uhn,u[DIM];
379: PetscInt i;
382: Scale2Real(1./x->h,x->uh,u);
383: uhn = x->uh[0] * n[0] + x->uh[1] * n[1];
384: f->h = uhn;
385: for (i=0; i<DIM; i++) f->uh[i] = u[i] * uhn + sw->gravity * PetscSqr(x->h) * n[i];
386: return(0);
387: }
389: static PetscErrorCode PhysicsBoundary_SW_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
390: {
392: xG[0] = xI[0];
393: xG[1] = -xI[1];
394: xG[2] = -xI[2];
395: return(0);
396: }
398: static void PhysicsRiemann_SW_HLL(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
399: {
400: Physics_SW *sw = (Physics_SW *) phys->data;
401: PetscReal aL, aR;
402: PetscReal nn[DIM];
403: #if !defined(PETSC_USE_COMPLEX)
404: const SWNode *uL = (const SWNode *) xL, *uR = (const SWNode *) xR;
405: #else
406: SWNodeUnion uLreal, uRreal;
407: const SWNode *uL = &uLreal.swnode;
408: const SWNode *uR = &uRreal.swnode;
409: #endif
410: SWNodeUnion fL, fR;
411: PetscInt i;
412: PetscReal zero = 0.;
414: #if defined(PETSC_USE_COMPLEX)
415: uLreal.swnode.h = 0; uRreal.swnode.h = 0;
416: for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
417: for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
418: #endif
419: if (uL->h <= 0 || uR->h <= 0) {
420: for (i = 0; i < 1 + dim; i++) flux[i] = zero;
421: return;
422: } /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
423: nn[0] = n[0];
424: nn[1] = n[1];
425: Normalize2Real(nn);
426: SWFlux(phys, nn, uL, &(fL.swnode));
427: SWFlux(phys, nn, uR, &(fR.swnode));
428: /* gravity wave speed */
429: aL = PetscSqrtReal(sw->gravity * uL->h);
430: aR = PetscSqrtReal(sw->gravity * uR->h);
431: // Defining u_tilda and v_tilda as u and v
432: PetscReal u_L, u_R;
433: u_L = Dot2Real(uL->uh,nn)/uL->h;
434: u_R = Dot2Real(uR->uh,nn)/uR->h;
435: PetscReal sL, sR;
436: sL = PetscMin(u_L - aL, u_R - aR);
437: sR = PetscMax(u_L + aL, u_R + aR);
438: if (sL > zero) {
439: for (i = 0; i < dim + 1; i++) {
440: flux[i] = fL.vals[i] * Norm2Real(n);
441: }
442: } else if (sR < zero) {
443: for (i = 0; i < dim + 1; i++) {
444: flux[i] = fR.vals[i] * Norm2Real(n);
445: }
446: } else {
447: for (i = 0; i < dim + 1; i++) {
448: flux[i] = ((sR * fL.vals[i] - sL * fR.vals[i] + sR * sL * (xR[i] - xL[i])) / (sR - sL)) * Norm2Real(n);
449: }
450: }
451: }
453: static void PhysicsRiemann_SW_Rusanov(PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n, const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
454: {
455: Physics_SW *sw = (Physics_SW*)phys->data;
456: PetscReal cL,cR,speed;
457: PetscReal nn[DIM];
458: #if !defined(PETSC_USE_COMPLEX)
459: const SWNode *uL = (const SWNode*)xL,*uR = (const SWNode*)xR;
460: #else
461: SWNodeUnion uLreal, uRreal;
462: const SWNode *uL = &uLreal.swnode;
463: const SWNode *uR = &uRreal.swnode;
464: #endif
465: SWNodeUnion fL,fR;
466: PetscInt i;
467: PetscReal zero=0.;
469: #if defined(PETSC_USE_COMPLEX)
470: uLreal.swnode.h = 0; uRreal.swnode.h = 0;
471: for (i = 0; i < 1+dim; i++) uLreal.vals[i] = PetscRealPart(xL[i]);
472: for (i = 0; i < 1+dim; i++) uRreal.vals[i] = PetscRealPart(xR[i]);
473: #endif
474: if (uL->h < 0 || uR->h < 0) {for (i=0; i<1+dim; i++) flux[i] = zero/zero; return;} /* SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative"); */
475: nn[0] = n[0];
476: nn[1] = n[1];
477: Normalize2Real(nn);
478: SWFlux(phys,nn,uL,&(fL.swnode));
479: SWFlux(phys,nn,uR,&(fR.swnode));
480: cL = PetscSqrtReal(sw->gravity*uL->h);
481: cR = PetscSqrtReal(sw->gravity*uR->h); /* gravity wave speed */
482: speed = PetscMax(PetscAbsReal(Dot2Real(uL->uh,nn)/uL->h) + cL,PetscAbsReal(Dot2Real(uR->uh,nn)/uR->h) + cR);
483: for (i=0; i<1+dim; i++) flux[i] = (0.5*(fL.vals[i] + fR.vals[i]) + 0.5*speed*(xL[i] - xR[i])) * Norm2Real(n);
484: }
486: static PetscErrorCode PhysicsSolution_SW(Model mod,PetscReal time,const PetscReal *x,PetscScalar *u,void *ctx)
487: {
488: PetscReal dx[2],r,sigma;
491: if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);
492: dx[0] = x[0] - 1.5;
493: dx[1] = x[1] - 1.0;
494: r = Norm2Real(dx);
495: sigma = 0.5;
496: u[0] = 1 + 2*PetscExpReal(-PetscSqr(r)/(2*PetscSqr(sigma)));
497: u[1] = 0.0;
498: u[2] = 0.0;
499: return(0);
500: }
502: static PetscErrorCode PhysicsFunctional_SW(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
503: {
504: Physics phys = (Physics)ctx;
505: Physics_SW *sw = (Physics_SW*)phys->data;
506: const SWNode *x = (const SWNode*)xx;
507: PetscReal u[2];
508: PetscReal h;
511: h = x->h;
512: Scale2Real(1./x->h,x->uh,u);
513: f[sw->functional.Height] = h;
514: f[sw->functional.Speed] = Norm2Real(u) + PetscSqrtReal(sw->gravity*h);
515: f[sw->functional.Energy] = 0.5*(Dot2Real(x->uh,u) + sw->gravity*PetscSqr(h));
516: return(0);
517: }
519: static PetscErrorCode SetUpBC_SW(PetscDS prob,Physics phys)
520: {
522: const PetscInt wallids[] = {100,101,200,300};
524: PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_SW_Wall, NULL, ALEN(wallids), wallids, phys);
525: return(0);
526: }
528: static PetscErrorCode PhysicsCreate_SW(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
529: {
530: Physics_SW *sw;
531: char sw_riemann[64] = "rusanov";
535: phys->field_desc = PhysicsFields_SW;
536: PetscNew(&sw);
537: phys->data = sw;
538: mod->setupbc = SetUpBC_SW;
540: PetscFunctionListAdd(&PhysicsRiemannList_SW, "rusanov", PhysicsRiemann_SW_Rusanov);
541: PetscFunctionListAdd(&PhysicsRiemannList_SW, "hll", PhysicsRiemann_SW_HLL);
543: PetscOptionsHead(PetscOptionsObject,"SW options");
544: {
545: void (*PhysicsRiemann_SW)(PetscInt, PetscInt, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, PetscInt, const PetscScalar, PetscScalar *, Physics);
546: sw->gravity = 1.0;
547: PetscOptionsReal("-sw_gravity","Gravitational constant","",sw->gravity,&sw->gravity,NULL);
548: PetscOptionsFList("-sw_riemann","Riemann solver","",PhysicsRiemannList_SW,sw_riemann,sw_riemann,sizeof sw_riemann,NULL);
549: PetscFunctionListFind(PhysicsRiemannList_SW,sw_riemann,&PhysicsRiemann_SW);
550: phys->riemann = (PetscRiemannFunc) PhysicsRiemann_SW;
551: }
552: PetscOptionsTail();
553: phys->maxspeed = PetscSqrtReal(2.0*sw->gravity); /* Mach 1 for depth of 2 */
555: ModelSolutionSetDefault(mod,PhysicsSolution_SW,phys);
556: ModelFunctionalRegister(mod,"Height",&sw->functional.Height,PhysicsFunctional_SW,phys);
557: ModelFunctionalRegister(mod,"Speed",&sw->functional.Speed,PhysicsFunctional_SW,phys);
558: ModelFunctionalRegister(mod,"Energy",&sw->functional.Energy,PhysicsFunctional_SW,phys);
560: mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
562: return(0);
563: }
565: /******************* Euler Density Shock (EULER_IV_SHOCK,EULER_SS_SHOCK) ********************/
566: /* An initial-value and self-similar solutions of the compressible Euler equations */
567: /* Ravi Samtaney and D. I. Pullin */
568: /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
569: typedef enum {EULER_PAR_GAMMA,EULER_PAR_RHOR,EULER_PAR_AMACH,EULER_PAR_ITANA,EULER_PAR_SIZE} EulerParamIdx;
570: typedef enum {EULER_IV_SHOCK,EULER_SS_SHOCK,EULER_SHOCK_TUBE,EULER_LINEAR_WAVE} EulerType;
571: typedef struct {
572: PetscReal r;
573: PetscReal ru[DIM];
574: PetscReal E;
575: } EulerNode;
576: typedef union {
577: EulerNode eulernode;
578: PetscReal vals[DIM+2];
579: } EulerNodeUnion;
580: typedef PetscErrorCode (*EquationOfState)(const PetscReal*, const EulerNode*, PetscReal*);
581: typedef struct {
582: EulerType type;
583: PetscReal pars[EULER_PAR_SIZE];
584: EquationOfState sound;
585: struct {
586: PetscInt Density;
587: PetscInt Momentum;
588: PetscInt Energy;
589: PetscInt Pressure;
590: PetscInt Speed;
591: } monitor;
592: } Physics_Euler;
594: static const struct FieldDescription PhysicsFields_Euler[] = {{"Density",1},{"Momentum",DIM},{"Energy",1},{NULL,0}};
596: /* initial condition */
597: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx);
598: static PetscErrorCode PhysicsSolution_Euler(Model mod, PetscReal time, const PetscReal *x, PetscScalar *u, void *ctx)
599: {
600: PetscInt i;
601: Physics phys = (Physics)ctx;
602: Physics_Euler *eu = (Physics_Euler*)phys->data;
603: EulerNode *uu = (EulerNode*)u;
604: PetscReal p0,gamma,c;
606: if (time != 0.0) SETERRQ1(mod->comm,PETSC_ERR_SUP,"No solution known for time %g",(double)time);
608: for (i=0; i<DIM; i++) uu->ru[i] = 0.0; /* zero out initial velocity */
609: /* set E and rho */
610: gamma = eu->pars[EULER_PAR_GAMMA];
612: if (eu->type==EULER_IV_SHOCK || eu->type==EULER_SS_SHOCK) {
613: /******************* Euler Density Shock ********************/
614: /* On initial-value and self-similar solutions of the compressible Euler equations */
615: /* Ravi Samtaney and D. I. Pullin */
616: /* Phys. Fluids 8, 2650 (1996); http://dx.doi.org/10.1063/1.869050 */
617: /* initial conditions 1: left of shock, 0: left of discontinuity 2: right of discontinuity, */
618: p0 = 1.;
619: if (x[0] < 0.0 + x[1]*eu->pars[EULER_PAR_ITANA]) {
620: if (x[0] < mod->bounds[0]*0.5) { /* left of shock (1) */
621: PetscReal amach,rho,press,gas1,p1;
622: amach = eu->pars[EULER_PAR_AMACH];
623: rho = 1.;
624: press = p0;
625: p1 = press*(1.0+2.0*gamma/(gamma+1.0)*(amach*amach-1.0));
626: gas1 = (gamma-1.0)/(gamma+1.0);
627: uu->r = rho*(p1/press+gas1)/(gas1*p1/press+1.0);
628: uu->ru[0] = ((uu->r - rho)*PetscSqrtReal(gamma*press/rho)*amach);
629: uu->E = p1/(gamma-1.0) + .5/uu->r*uu->ru[0]*uu->ru[0];
630: }
631: else { /* left of discontinuity (0) */
632: uu->r = 1.; /* rho = 1 */
633: uu->E = p0/(gamma-1.0);
634: }
635: }
636: else { /* right of discontinuity (2) */
637: uu->r = eu->pars[EULER_PAR_RHOR];
638: uu->E = p0/(gamma-1.0);
639: }
640: }
641: else if (eu->type==EULER_SHOCK_TUBE) {
642: /* For (x<x0) set (rho,u,p)=(8,0,10) and for (x>x0) set (rho,u,p)=(1,0,1). Choose x0 to the midpoint of the domain in the x-direction. */
643: if (x[0] < 0.0) {
644: uu->r = 8.;
645: uu->E = 10./(gamma-1.);
646: }
647: else {
648: uu->r = 1.;
649: uu->E = 1./(gamma-1.);
650: }
651: }
652: else if (eu->type==EULER_LINEAR_WAVE) {
653: initLinearWave( uu, gamma, x, mod->bounds[1] - mod->bounds[0]);
654: }
655: else SETERRQ1(mod->comm,PETSC_ERR_SUP,"Unknown type %d",eu->type);
657: /* set phys->maxspeed: (mod->maxspeed = phys->maxspeed) in main; */
658: eu->sound(&gamma,uu,&c);
659: c = (uu->ru[0]/uu->r) + c;
660: if (c > phys->maxspeed) phys->maxspeed = c;
662: return(0);
663: }
665: static PetscErrorCode Pressure_PG(const PetscReal gamma,const EulerNode *x,PetscReal *p)
666: {
667: PetscReal ru2;
670: ru2 = DotDIMReal(x->ru,x->ru);
671: (*p)=(x->E - 0.5*ru2/x->r)*(gamma - 1.0); /* (E - rho V^2/2)(gamma-1) = e rho (gamma-1) */
672: return(0);
673: }
675: static PetscErrorCode SpeedOfSound_PG(const PetscReal *gamma, const EulerNode *x, PetscReal *c)
676: {
677: PetscReal p;
680: Pressure_PG(*gamma,x,&p);
681: if (p<0.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"negative pressure time %g -- NEED TO FIX!!!!!!",(double) p);
682: /* pars[EULER_PAR_GAMMA] = heat capacity ratio */
683: (*c)=PetscSqrtReal(*gamma * p / x->r);
684: return(0);
685: }
687: /*
688: * x = (rho,rho*(u_1),...,rho*e)^T
689: * x_t+div(f_1(x))+...+div(f_DIM(x)) = 0
690: *
691: * f_i(x) = u_i*x+(0,0,...,p,...,p*u_i)^T
692: *
693: */
694: static PetscErrorCode EulerFlux(Physics phys,const PetscReal *n,const EulerNode *x,EulerNode *f)
695: {
696: Physics_Euler *eu = (Physics_Euler*)phys->data;
697: PetscReal nu,p;
698: PetscInt i;
701: Pressure_PG(eu->pars[EULER_PAR_GAMMA],x,&p);
702: nu = DotDIMReal(x->ru,n);
703: f->r = nu; /* A rho u */
704: nu /= x->r; /* A u */
705: for (i=0; i<DIM; i++) f->ru[i] = nu * x->ru[i] + n[i]*p; /* r u^2 + p */
706: f->E = nu * (x->E + p); /* u(e+p) */
707: return(0);
708: }
710: /* PetscReal* => EulerNode* conversion */
711: static PetscErrorCode PhysicsBoundary_Euler_Wall(PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *a_xI, PetscScalar *a_xG, void *ctx)
712: {
713: PetscInt i;
714: const EulerNode *xI = (const EulerNode*)a_xI;
715: EulerNode *xG = (EulerNode*)a_xG;
716: Physics phys = (Physics)ctx;
717: Physics_Euler *eu = (Physics_Euler*)phys->data;
719: xG->r = xI->r; /* ghost cell density - same */
720: xG->E = xI->E; /* ghost cell energy - same */
721: if (n[1] != 0.) { /* top and bottom */
722: xG->ru[0] = xI->ru[0]; /* copy tang to wall */
723: xG->ru[1] = -xI->ru[1]; /* reflect perp to t/b wall */
724: }
725: else { /* sides */
726: for (i=0; i<DIM; i++) xG->ru[i] = xI->ru[i]; /* copy */
727: }
728: if (eu->type == EULER_LINEAR_WAVE) { /* debug */
729: #if 0
730: PetscPrintf(PETSC_COMM_WORLD,"%s coord=%g,%g\n",PETSC_FUNCTION_NAME,c[0],c[1]);
731: #endif
732: }
733: return(0);
734: }
735: int godunovflux( const PetscScalar *ul, const PetscScalar *ur, PetscScalar *flux, const PetscReal *nn, const int *ndim, const PetscReal *gamma);
736: /* PetscReal* => EulerNode* conversion */
737: static void PhysicsRiemann_Euler_Godunov( PetscInt dim, PetscInt Nf, const PetscReal *qp, const PetscReal *n,
738: const PetscScalar *xL, const PetscScalar *xR, PetscInt numConstants, const PetscScalar constants[], PetscScalar *flux, Physics phys)
739: {
740: Physics_Euler *eu = (Physics_Euler*)phys->data;
741: PetscReal cL,cR,speed,velL,velR,nn[DIM],s2;
742: PetscInt i;
743: PetscErrorCode ierr;
746: for (i=0,s2=0.; i<DIM; i++) {
747: nn[i] = n[i];
748: s2 += nn[i]*nn[i];
749: }
750: s2 = PetscSqrtReal(s2); /* |n|_2 = sum(n^2)^1/2 */
751: for (i=0.; i<DIM; i++) nn[i] /= s2;
752: if (0) { /* Rusanov */
753: const EulerNode *uL = (const EulerNode*)xL,*uR = (const EulerNode*)xR;
754: EulerNodeUnion fL,fR;
755: EulerFlux(phys,nn,uL,&(fL.eulernode));
756: EulerFlux(phys,nn,uR,&(fR.eulernode));
757: eu->sound(&eu->pars[EULER_PAR_GAMMA],uL,&cL);if (ierr) exit(13);
758: eu->sound(&eu->pars[EULER_PAR_GAMMA],uR,&cR);if (ierr) exit(14);
759: velL = DotDIMReal(uL->ru,nn)/uL->r;
760: velR = DotDIMReal(uR->ru,nn)/uR->r;
761: speed = PetscMax(velR + cR, velL + cL);
762: for (i=0; i<2+dim; i++) flux[i] = 0.5*((fL.vals[i]+fR.vals[i]) + speed*(xL[i] - xR[i]))*s2;
763: }
764: else {
765: int dim = DIM;
766: /* int iwave = */
767: godunovflux(xL, xR, flux, nn, &dim, &eu->pars[EULER_PAR_GAMMA]);
768: for (i=0; i<2+dim; i++) flux[i] *= s2;
769: }
770: PetscFunctionReturnVoid();
771: }
773: static PetscErrorCode PhysicsFunctional_Euler(Model mod,PetscReal time,const PetscReal *coord,const PetscScalar *xx,PetscReal *f,void *ctx)
774: {
775: Physics phys = (Physics)ctx;
776: Physics_Euler *eu = (Physics_Euler*)phys->data;
777: const EulerNode *x = (const EulerNode*)xx;
778: PetscReal p;
781: f[eu->monitor.Density] = x->r;
782: f[eu->monitor.Momentum] = NormDIM(x->ru);
783: f[eu->monitor.Energy] = x->E;
784: f[eu->monitor.Speed] = NormDIM(x->ru)/x->r;
785: Pressure_PG(eu->pars[EULER_PAR_GAMMA], x, &p);
786: f[eu->monitor.Pressure] = p;
787: return(0);
788: }
790: static PetscErrorCode SetUpBC_Euler(PetscDS prob,Physics phys)
791: {
792: PetscErrorCode ierr;
793: Physics_Euler *eu = (Physics_Euler *) phys->data;
794: if (eu->type == EULER_LINEAR_WAVE) {
795: const PetscInt wallids[] = {100,101};
796: PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, ALEN(wallids), wallids, phys);
797: }
798: else {
799: const PetscInt wallids[] = {100,101,200,300};
800: PetscDSAddBoundary(prob, DM_BC_NATURAL_RIEMANN, "wall", "Face Sets", 0, 0, NULL, (void (*)(void)) PhysicsBoundary_Euler_Wall, NULL, ALEN(wallids), wallids, phys);
801: }
802: return(0);
803: }
805: static PetscErrorCode PhysicsCreate_Euler(Model mod,Physics phys,PetscOptionItems *PetscOptionsObject)
806: {
807: Physics_Euler *eu;
808: PetscErrorCode ierr;
811: phys->field_desc = PhysicsFields_Euler;
812: phys->riemann = (PetscRiemannFunc) PhysicsRiemann_Euler_Godunov;
813: PetscNew(&eu);
814: phys->data = eu;
815: mod->setupbc = SetUpBC_Euler;
816: PetscOptionsHead(PetscOptionsObject,"Euler options");
817: {
818: PetscReal alpha;
819: char type[64] = "linear_wave";
820: PetscBool is;
821: mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_GHOSTED;
822: eu->pars[EULER_PAR_GAMMA] = 1.4;
823: eu->pars[EULER_PAR_AMACH] = 2.02;
824: eu->pars[EULER_PAR_RHOR] = 3.0;
825: eu->pars[EULER_PAR_ITANA] = 0.57735026918963; /* angle of Euler self similar (SS) shock */
826: PetscOptionsReal("-eu_gamma","Heat capacity ratio","",eu->pars[EULER_PAR_GAMMA],&eu->pars[EULER_PAR_GAMMA],NULL);
827: PetscOptionsReal("-eu_amach","Shock speed (Mach)","",eu->pars[EULER_PAR_AMACH],&eu->pars[EULER_PAR_AMACH],NULL);
828: PetscOptionsReal("-eu_rho2","Density right of discontinuity","",eu->pars[EULER_PAR_RHOR],&eu->pars[EULER_PAR_RHOR],NULL);
829: alpha = 60.;
830: PetscOptionsReal("-eu_alpha","Angle of discontinuity","",alpha,&alpha,NULL);
831: if (alpha<=0. || alpha>90.) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Alpha bust be > 0 and <= 90 (%g)",alpha);
832: eu->pars[EULER_PAR_ITANA] = 1./PetscTanReal( alpha * PETSC_PI / 180.0);
833: PetscOptionsString("-eu_type","Type of Euler test","",type,type,sizeof(type),NULL);
834: PetscStrcmp(type,"linear_wave", &is);
835: if (is) {
836: eu->type = EULER_LINEAR_WAVE;
837: mod->bcs[0] = mod->bcs[1] = mod->bcs[2] = DM_BOUNDARY_PERIODIC;
838: mod->bcs[1] = DM_BOUNDARY_GHOSTED; /* debug */
839: PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"linear_wave");
840: }
841: else {
842: if (DIM != 2) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DIM must be 2 unless linear wave test %s",type);
843: PetscStrcmp(type,"iv_shock", &is);
844: if (is) {
845: eu->type = EULER_IV_SHOCK;
846: PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"iv_shock");
847: }
848: else {
849: PetscStrcmp(type,"ss_shock", &is);
850: if (is) {
851: eu->type = EULER_SS_SHOCK;
852: PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"ss_shock");
853: }
854: else {
855: PetscStrcmp(type,"shock_tube", &is);
856: if (is) eu->type = EULER_SHOCK_TUBE;
857: else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unknown Euler type %s",type);
858: PetscPrintf(PETSC_COMM_WORLD,"%s set Euler type: %s\n",PETSC_FUNCTION_NAME,"shock_tube");
859: }
860: }
861: }
862: }
863: PetscOptionsTail();
864: eu->sound = SpeedOfSound_PG;
865: phys->maxspeed = 0.; /* will get set in solution */
866: ModelSolutionSetDefault(mod,PhysicsSolution_Euler,phys);
867: ModelFunctionalRegister(mod,"Speed",&eu->monitor.Speed,PhysicsFunctional_Euler,phys);
868: ModelFunctionalRegister(mod,"Energy",&eu->monitor.Energy,PhysicsFunctional_Euler,phys);
869: ModelFunctionalRegister(mod,"Density",&eu->monitor.Density,PhysicsFunctional_Euler,phys);
870: ModelFunctionalRegister(mod,"Momentum",&eu->monitor.Momentum,PhysicsFunctional_Euler,phys);
871: ModelFunctionalRegister(mod,"Pressure",&eu->monitor.Pressure,PhysicsFunctional_Euler,phys);
873: return(0);
874: }
876: static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscInt numComps, const PetscScalar u[], const PetscScalar grad[], PetscReal *error, void *ctx)
877: {
878: PetscReal err = 0.;
879: PetscInt i, j;
882: for (i = 0; i < numComps; i++) {
883: for (j = 0; j < dim; j++) {
884: err += PetscSqr(PetscRealPart(grad[i * dim + j]));
885: }
886: }
887: *error = volume * err;
888: return(0);
889: }
891: PetscErrorCode ConstructCellBoundary(DM dm, User user)
892: {
893: const char *name = "Cell Sets";
894: const char *bdname = "split faces";
895: IS regionIS, innerIS;
896: const PetscInt *regions, *cells;
897: PetscInt numRegions, innerRegion, numCells, c;
898: PetscInt cStart, cEnd, cEndInterior, fStart, fEnd;
899: PetscBool hasLabel;
903: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
904: DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
905: DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);
907: DMHasLabel(dm, name, &hasLabel);
908: if (!hasLabel) return(0);
909: DMGetLabelSize(dm, name, &numRegions);
910: if (numRegions != 2) return(0);
911: /* Get the inner id */
912: DMGetLabelIdIS(dm, name, ®ionIS);
913: ISGetIndices(regionIS, ®ions);
914: innerRegion = regions[0];
915: ISRestoreIndices(regionIS, ®ions);
916: ISDestroy(®ionIS);
917: /* Find the faces between cells in different regions, could call DMPlexCreateNeighborCSR() */
918: DMGetStratumIS(dm, name, innerRegion, &innerIS);
919: ISGetLocalSize(innerIS, &numCells);
920: ISGetIndices(innerIS, &cells);
921: DMCreateLabel(dm, bdname);
922: for (c = 0; c < numCells; ++c) {
923: const PetscInt cell = cells[c];
924: const PetscInt *faces;
925: PetscInt numFaces, f;
927: if ((cell < cStart) || (cell >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", cell);
928: DMPlexGetConeSize(dm, cell, &numFaces);
929: DMPlexGetCone(dm, cell, &faces);
930: for (f = 0; f < numFaces; ++f) {
931: const PetscInt face = faces[f];
932: const PetscInt *neighbors;
933: PetscInt nC, regionA, regionB;
935: if ((face < fStart) || (face >= fEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a face", face);
936: DMPlexGetSupportSize(dm, face, &nC);
937: if (nC != 2) continue;
938: DMPlexGetSupport(dm, face, &neighbors);
939: if ((neighbors[0] >= cEndInterior) || (neighbors[1] >= cEndInterior)) continue;
940: if ((neighbors[0] < cStart) || (neighbors[0] >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", neighbors[0]);
941: if ((neighbors[1] < cStart) || (neighbors[1] >= cEnd)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Got invalid point %d which is not a cell", neighbors[1]);
942: DMGetLabelValue(dm, name, neighbors[0], ®ionA);
943: DMGetLabelValue(dm, name, neighbors[1], ®ionB);
944: if (regionA < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[0]);
945: if (regionB < 0) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid label %s: Cell %d has no value", name, neighbors[1]);
946: if (regionA != regionB) {
947: DMSetLabelValue(dm, bdname, faces[f], 1);
948: }
949: }
950: }
951: ISRestoreIndices(innerIS, &cells);
952: ISDestroy(&innerIS);
953: {
954: DMLabel label;
956: DMGetLabel(dm, bdname, &label);
957: DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD);
958: }
959: return(0);
960: }
962: /* Right now, I have just added duplicate faces, which see both cells. We can
963: - Add duplicate vertices and decouple the face cones
964: - Disconnect faces from cells across the rotation gap
965: */
966: PetscErrorCode SplitFaces(DM *dmSplit, const char labelName[], User user)
967: {
968: DM dm = *dmSplit, sdm;
969: PetscSF sfPoint, gsfPoint;
970: PetscSection coordSection, newCoordSection;
971: Vec coordinates;
972: IS idIS;
973: const PetscInt *ids;
974: PetscInt *newpoints;
975: PetscInt dim, depth, maxConeSize, maxSupportSize, numLabels, numGhostCells;
976: PetscInt numFS, fs, pStart, pEnd, p, cEnd, cEndInterior, vStart, vEnd, v, fStart, fEnd, newf, d, l;
977: PetscBool hasLabel;
981: DMHasLabel(dm, labelName, &hasLabel);
982: if (!hasLabel) return(0);
983: DMCreate(PetscObjectComm((PetscObject)dm), &sdm);
984: DMSetType(sdm, DMPLEX);
985: DMGetDimension(dm, &dim);
986: DMSetDimension(sdm, dim);
988: DMGetLabelIdIS(dm, labelName, &idIS);
989: ISGetLocalSize(idIS, &numFS);
990: ISGetIndices(idIS, &ids);
992: user->numSplitFaces = 0;
993: for (fs = 0; fs < numFS; ++fs) {
994: PetscInt numBdFaces;
996: DMGetStratumSize(dm, labelName, ids[fs], &numBdFaces);
997: user->numSplitFaces += numBdFaces;
998: }
999: DMPlexGetChart(dm, &pStart, &pEnd);
1000: pEnd += user->numSplitFaces;
1001: DMPlexSetChart(sdm, pStart, pEnd);
1002: DMPlexGetGhostCellStratum(dm, &cEndInterior, NULL);
1003: DMPlexGetHeightStratum(dm, 0, NULL, &cEnd);
1004: numGhostCells = cEnd - cEndInterior;
1005: /* Set cone and support sizes */
1006: DMPlexGetDepth(dm, &depth);
1007: for (d = 0; d <= depth; ++d) {
1008: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
1009: for (p = pStart; p < pEnd; ++p) {
1010: PetscInt newp = p;
1011: PetscInt size;
1013: DMPlexGetConeSize(dm, p, &size);
1014: DMPlexSetConeSize(sdm, newp, size);
1015: DMPlexGetSupportSize(dm, p, &size);
1016: DMPlexSetSupportSize(sdm, newp, size);
1017: }
1018: }
1019: DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
1020: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1021: IS faceIS;
1022: const PetscInt *faces;
1023: PetscInt numFaces, f;
1025: DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
1026: ISGetLocalSize(faceIS, &numFaces);
1027: ISGetIndices(faceIS, &faces);
1028: for (f = 0; f < numFaces; ++f, ++newf) {
1029: PetscInt size;
1031: /* Right now I think that both faces should see both cells */
1032: DMPlexGetConeSize(dm, faces[f], &size);
1033: DMPlexSetConeSize(sdm, newf, size);
1034: DMPlexGetSupportSize(dm, faces[f], &size);
1035: DMPlexSetSupportSize(sdm, newf, size);
1036: }
1037: ISRestoreIndices(faceIS, &faces);
1038: ISDestroy(&faceIS);
1039: }
1040: DMSetUp(sdm);
1041: /* Set cones and supports */
1042: DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
1043: PetscMalloc1(PetscMax(maxConeSize, maxSupportSize), &newpoints);
1044: DMPlexGetChart(dm, &pStart, &pEnd);
1045: for (p = pStart; p < pEnd; ++p) {
1046: const PetscInt *points, *orientations;
1047: PetscInt size, i, newp = p;
1049: DMPlexGetConeSize(dm, p, &size);
1050: DMPlexGetCone(dm, p, &points);
1051: DMPlexGetConeOrientation(dm, p, &orientations);
1052: for (i = 0; i < size; ++i) newpoints[i] = points[i];
1053: DMPlexSetCone(sdm, newp, newpoints);
1054: DMPlexSetConeOrientation(sdm, newp, orientations);
1055: DMPlexGetSupportSize(dm, p, &size);
1056: DMPlexGetSupport(dm, p, &points);
1057: for (i = 0; i < size; ++i) newpoints[i] = points[i];
1058: DMPlexSetSupport(sdm, newp, newpoints);
1059: }
1060: PetscFree(newpoints);
1061: for (fs = 0, newf = fEnd; fs < numFS; ++fs) {
1062: IS faceIS;
1063: const PetscInt *faces;
1064: PetscInt numFaces, f;
1066: DMGetStratumIS(dm, labelName, ids[fs], &faceIS);
1067: ISGetLocalSize(faceIS, &numFaces);
1068: ISGetIndices(faceIS, &faces);
1069: for (f = 0; f < numFaces; ++f, ++newf) {
1070: const PetscInt *points;
1072: DMPlexGetCone(dm, faces[f], &points);
1073: DMPlexSetCone(sdm, newf, points);
1074: DMPlexGetSupport(dm, faces[f], &points);
1075: DMPlexSetSupport(sdm, newf, points);
1076: }
1077: ISRestoreIndices(faceIS, &faces);
1078: ISDestroy(&faceIS);
1079: }
1080: ISRestoreIndices(idIS, &ids);
1081: ISDestroy(&idIS);
1082: DMPlexStratify(sdm);
1083: /* Convert coordinates */
1084: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1085: DMGetCoordinateSection(dm, &coordSection);
1086: PetscSectionCreate(PetscObjectComm((PetscObject)dm), &newCoordSection);
1087: PetscSectionSetNumFields(newCoordSection, 1);
1088: PetscSectionSetFieldComponents(newCoordSection, 0, dim);
1089: PetscSectionSetChart(newCoordSection, vStart, vEnd);
1090: for (v = vStart; v < vEnd; ++v) {
1091: PetscSectionSetDof(newCoordSection, v, dim);
1092: PetscSectionSetFieldDof(newCoordSection, v, 0, dim);
1093: }
1094: PetscSectionSetUp(newCoordSection);
1095: DMSetCoordinateSection(sdm, PETSC_DETERMINE, newCoordSection);
1096: PetscSectionDestroy(&newCoordSection); /* relinquish our reference */
1097: DMGetCoordinatesLocal(dm, &coordinates);
1098: DMSetCoordinatesLocal(sdm, coordinates);
1099: /* Convert labels */
1100: DMGetNumLabels(dm, &numLabels);
1101: for (l = 0; l < numLabels; ++l) {
1102: const char *lname;
1103: PetscBool isDepth, isDim;
1105: DMGetLabelName(dm, l, &lname);
1106: PetscStrcmp(lname, "depth", &isDepth);
1107: if (isDepth) continue;
1108: PetscStrcmp(lname, "dim", &isDim);
1109: if (isDim) continue;
1110: DMCreateLabel(sdm, lname);
1111: DMGetLabelIdIS(dm, lname, &idIS);
1112: ISGetLocalSize(idIS, &numFS);
1113: ISGetIndices(idIS, &ids);
1114: for (fs = 0; fs < numFS; ++fs) {
1115: IS pointIS;
1116: const PetscInt *points;
1117: PetscInt numPoints;
1119: DMGetStratumIS(dm, lname, ids[fs], &pointIS);
1120: ISGetLocalSize(pointIS, &numPoints);
1121: ISGetIndices(pointIS, &points);
1122: for (p = 0; p < numPoints; ++p) {
1123: PetscInt newpoint = points[p];
1125: DMSetLabelValue(sdm, lname, newpoint, ids[fs]);
1126: }
1127: ISRestoreIndices(pointIS, &points);
1128: ISDestroy(&pointIS);
1129: }
1130: ISRestoreIndices(idIS, &ids);
1131: ISDestroy(&idIS);
1132: }
1133: {
1134: /* Convert pointSF */
1135: const PetscSFNode *remotePoints;
1136: PetscSFNode *gremotePoints;
1137: const PetscInt *localPoints;
1138: PetscInt *glocalPoints,*newLocation,*newRemoteLocation;
1139: PetscInt numRoots, numLeaves;
1140: PetscMPIInt size;
1142: MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
1143: DMGetPointSF(dm, &sfPoint);
1144: DMGetPointSF(sdm, &gsfPoint);
1145: DMPlexGetChart(dm,&pStart,&pEnd);
1146: PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
1147: if (numRoots >= 0) {
1148: PetscMalloc2(numRoots,&newLocation,pEnd-pStart,&newRemoteLocation);
1149: for (l=0; l<numRoots; l++) newLocation[l] = l; /* + (l >= cEnd ? numGhostCells : 0); */
1150: PetscSFBcastBegin(sfPoint, MPIU_INT, newLocation, newRemoteLocation,MPI_REPLACE);
1151: PetscSFBcastEnd(sfPoint, MPIU_INT, newLocation, newRemoteLocation,MPI_REPLACE);
1152: PetscMalloc1(numLeaves, &glocalPoints);
1153: PetscMalloc1(numLeaves, &gremotePoints);
1154: for (l = 0; l < numLeaves; ++l) {
1155: glocalPoints[l] = localPoints[l]; /* localPoints[l] >= cEnd ? localPoints[l] + numGhostCells : localPoints[l]; */
1156: gremotePoints[l].rank = remotePoints[l].rank;
1157: gremotePoints[l].index = newRemoteLocation[localPoints[l]];
1158: }
1159: PetscFree2(newLocation,newRemoteLocation);
1160: PetscSFSetGraph(gsfPoint, numRoots+numGhostCells, numLeaves, glocalPoints, PETSC_OWN_POINTER, gremotePoints, PETSC_OWN_POINTER);
1161: }
1162: DMDestroy(dmSplit);
1163: *dmSplit = sdm;
1164: }
1165: return(0);
1166: }
1168: PetscErrorCode CreatePartitionVec(DM dm, DM *dmCell, Vec *partition)
1169: {
1170: PetscSF sfPoint;
1171: PetscSection coordSection;
1172: Vec coordinates;
1173: PetscSection sectionCell;
1174: PetscScalar *part;
1175: PetscInt cStart, cEnd, c;
1176: PetscMPIInt rank;
1180: DMGetCoordinateSection(dm, &coordSection);
1181: DMGetCoordinatesLocal(dm, &coordinates);
1182: DMClone(dm, dmCell);
1183: DMGetPointSF(dm, &sfPoint);
1184: DMSetPointSF(*dmCell, sfPoint);
1185: DMSetCoordinateSection(*dmCell, PETSC_DETERMINE, coordSection);
1186: DMSetCoordinatesLocal(*dmCell, coordinates);
1187: MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
1188: PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionCell);
1189: DMPlexGetHeightStratum(*dmCell, 0, &cStart, &cEnd);
1190: PetscSectionSetChart(sectionCell, cStart, cEnd);
1191: for (c = cStart; c < cEnd; ++c) {
1192: PetscSectionSetDof(sectionCell, c, 1);
1193: }
1194: PetscSectionSetUp(sectionCell);
1195: DMSetLocalSection(*dmCell, sectionCell);
1196: PetscSectionDestroy(§ionCell);
1197: DMCreateLocalVector(*dmCell, partition);
1198: PetscObjectSetName((PetscObject)*partition, "partition");
1199: VecGetArray(*partition, &part);
1200: for (c = cStart; c < cEnd; ++c) {
1201: PetscScalar *p;
1203: DMPlexPointLocalRef(*dmCell, c, part, &p);
1204: p[0] = rank;
1205: }
1206: VecRestoreArray(*partition, &part);
1207: return(0);
1208: }
1210: PetscErrorCode CreateMassMatrix(DM dm, Vec *massMatrix, User user)
1211: {
1212: DM plex, dmMass, dmFace, dmCell, dmCoord;
1213: PetscSection coordSection;
1214: Vec coordinates, facegeom, cellgeom;
1215: PetscSection sectionMass;
1216: PetscScalar *m;
1217: const PetscScalar *fgeom, *cgeom, *coords;
1218: PetscInt vStart, vEnd, v;
1219: PetscErrorCode ierr;
1222: DMConvert(dm, DMPLEX, &plex);
1223: DMGetCoordinateSection(dm, &coordSection);
1224: DMGetCoordinatesLocal(dm, &coordinates);
1225: DMClone(dm, &dmMass);
1226: DMSetCoordinateSection(dmMass, PETSC_DETERMINE, coordSection);
1227: DMSetCoordinatesLocal(dmMass, coordinates);
1228: PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ionMass);
1229: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1230: PetscSectionSetChart(sectionMass, vStart, vEnd);
1231: for (v = vStart; v < vEnd; ++v) {
1232: PetscInt numFaces;
1234: DMPlexGetSupportSize(dmMass, v, &numFaces);
1235: PetscSectionSetDof(sectionMass, v, numFaces*numFaces);
1236: }
1237: PetscSectionSetUp(sectionMass);
1238: DMSetLocalSection(dmMass, sectionMass);
1239: PetscSectionDestroy(§ionMass);
1240: DMGetLocalVector(dmMass, massMatrix);
1241: VecGetArray(*massMatrix, &m);
1242: DMPlexGetGeometryFVM(plex, &facegeom, &cellgeom, NULL);
1243: VecGetDM(facegeom, &dmFace);
1244: VecGetArrayRead(facegeom, &fgeom);
1245: VecGetDM(cellgeom, &dmCell);
1246: VecGetArrayRead(cellgeom, &cgeom);
1247: DMGetCoordinateDM(dm, &dmCoord);
1248: VecGetArrayRead(coordinates, &coords);
1249: for (v = vStart; v < vEnd; ++v) {
1250: const PetscInt *faces;
1251: PetscFVFaceGeom *fgA, *fgB, *cg;
1252: PetscScalar *vertex;
1253: PetscInt numFaces, sides[2], f, g;
1255: DMPlexPointLocalRead(dmCoord, v, coords, &vertex);
1256: DMPlexGetSupportSize(dmMass, v, &numFaces);
1257: DMPlexGetSupport(dmMass, v, &faces);
1258: for (f = 0; f < numFaces; ++f) {
1259: sides[0] = faces[f];
1260: DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fgA);
1261: for (g = 0; g < numFaces; ++g) {
1262: const PetscInt *cells = NULL;
1263: PetscReal area = 0.0;
1264: PetscInt numCells;
1266: sides[1] = faces[g];
1267: DMPlexPointLocalRead(dmFace, faces[g], fgeom, &fgB);
1268: DMPlexGetJoin(dmMass, 2, sides, &numCells, &cells);
1269: if (numCells != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Invalid join for faces");
1270: DMPlexPointLocalRead(dmCell, cells[0], cgeom, &cg);
1271: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgA->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgA->centroid[0] - cg->centroid[0]));
1272: area += PetscAbsScalar((vertex[0] - cg->centroid[0])*(fgB->centroid[1] - cg->centroid[1]) - (vertex[1] - cg->centroid[1])*(fgB->centroid[0] - cg->centroid[0]));
1273: m[f*numFaces+g] = Dot2Real(fgA->normal, fgB->normal)*area*0.5;
1274: DMPlexRestoreJoin(dmMass, 2, sides, &numCells, &cells);
1275: }
1276: }
1277: }
1278: VecRestoreArrayRead(facegeom, &fgeom);
1279: VecRestoreArrayRead(cellgeom, &cgeom);
1280: VecRestoreArrayRead(coordinates, &coords);
1281: VecRestoreArray(*massMatrix, &m);
1282: DMDestroy(&dmMass);
1283: DMDestroy(&plex);
1284: return(0);
1285: }
1287: /* Behavior will be different for multi-physics or when using non-default boundary conditions */
1288: static PetscErrorCode ModelSolutionSetDefault(Model mod,SolutionFunction func,void *ctx)
1289: {
1291: mod->solution = func;
1292: mod->solutionctx = ctx;
1293: return(0);
1294: }
1296: static PetscErrorCode ModelFunctionalRegister(Model mod,const char *name,PetscInt *offset,FunctionalFunction func,void *ctx)
1297: {
1299: FunctionalLink link,*ptr;
1300: PetscInt lastoffset = -1;
1303: for (ptr=&mod->functionalRegistry; *ptr; ptr = &(*ptr)->next) lastoffset = (*ptr)->offset;
1304: PetscNew(&link);
1305: PetscStrallocpy(name,&link->name);
1306: link->offset = lastoffset + 1;
1307: link->func = func;
1308: link->ctx = ctx;
1309: link->next = NULL;
1310: *ptr = link;
1311: *offset = link->offset;
1312: return(0);
1313: }
1315: static PetscErrorCode ModelFunctionalSetFromOptions(Model mod,PetscOptionItems *PetscOptionsObject)
1316: {
1318: PetscInt i,j;
1319: FunctionalLink link;
1320: char *names[256];
1323: mod->numMonitored = ALEN(names);
1324: PetscOptionsStringArray("-monitor","list of functionals to monitor","",names,&mod->numMonitored,NULL);
1325: /* Create list of functionals that will be computed somehow */
1326: PetscMalloc1(mod->numMonitored,&mod->functionalMonitored);
1327: /* Create index of calls that we will have to make to compute these functionals (over-allocation in general). */
1328: PetscMalloc1(mod->numMonitored,&mod->functionalCall);
1329: mod->numCall = 0;
1330: for (i=0; i<mod->numMonitored; i++) {
1331: for (link=mod->functionalRegistry; link; link=link->next) {
1332: PetscBool match;
1333: PetscStrcasecmp(names[i],link->name,&match);
1334: if (match) break;
1335: }
1336: if (!link) SETERRQ1(mod->comm,PETSC_ERR_USER,"No known functional '%s'",names[i]);
1337: mod->functionalMonitored[i] = link;
1338: for (j=0; j<i; j++) {
1339: if (mod->functionalCall[j]->func == link->func && mod->functionalCall[j]->ctx == link->ctx) goto next_name;
1340: }
1341: mod->functionalCall[mod->numCall++] = link; /* Just points to the first link using the result. There may be more results. */
1342: next_name:
1343: PetscFree(names[i]);
1344: }
1346: /* Find out the maximum index of any functional computed by a function we will be calling (even if we are not using it) */
1347: mod->maxComputed = -1;
1348: for (link=mod->functionalRegistry; link; link=link->next) {
1349: for (i=0; i<mod->numCall; i++) {
1350: FunctionalLink call = mod->functionalCall[i];
1351: if (link->func == call->func && link->ctx == call->ctx) {
1352: mod->maxComputed = PetscMax(mod->maxComputed,link->offset);
1353: }
1354: }
1355: }
1356: return(0);
1357: }
1359: static PetscErrorCode FunctionalLinkDestroy(FunctionalLink *link)
1360: {
1362: FunctionalLink l,next;
1365: if (!link) return(0);
1366: l = *link;
1367: *link = NULL;
1368: for (; l; l=next) {
1369: next = l->next;
1370: PetscFree(l->name);
1371: PetscFree(l);
1372: }
1373: return(0);
1374: }
1376: /* put the solution callback into a functional callback */
1377: static PetscErrorCode SolutionFunctional(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *modctx)
1378: {
1379: Model mod;
1382: mod = (Model) modctx;
1383: (*mod->solution)(mod, time, x, u, mod->solutionctx);
1384: return(0);
1385: }
1387: PetscErrorCode SetInitialCondition(DM dm, Vec X, User user)
1388: {
1389: PetscErrorCode (*func[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1390: void *ctx[1];
1391: Model mod = user->model;
1392: PetscErrorCode ierr;
1395: func[0] = SolutionFunctional;
1396: ctx[0] = (void *) mod;
1397: DMProjectFunction(dm,0.0,func,ctx,INSERT_ALL_VALUES,X);
1398: return(0);
1399: }
1401: static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1402: {
1406: PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer);
1407: PetscViewerSetType(*viewer, PETSCVIEWERVTK);
1408: PetscViewerFileSetName(*viewer, filename);
1409: return(0);
1410: }
1412: static PetscErrorCode MonitorVTK(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx)
1413: {
1414: User user = (User)ctx;
1415: DM dm, plex;
1416: PetscViewer viewer;
1417: char filename[PETSC_MAX_PATH_LEN],*ftable = NULL;
1418: PetscReal xnorm;
1422: PetscObjectSetName((PetscObject) X, "u");
1423: VecGetDM(X,&dm);
1424: VecNorm(X,NORM_INFINITY,&xnorm);
1426: if (stepnum >= 0) {
1427: stepnum += user->monitorStepOffset;
1428: }
1429: if (stepnum >= 0) { /* No summary for final time */
1430: Model mod = user->model;
1431: Vec cellgeom;
1432: PetscInt c,cStart,cEnd,fcount,i;
1433: size_t ftableused,ftablealloc;
1434: const PetscScalar *cgeom,*x;
1435: DM dmCell;
1436: DMLabel vtkLabel;
1437: PetscReal *fmin,*fmax,*fintegral,*ftmp;
1439: DMConvert(dm, DMPLEX, &plex);
1440: DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);
1441: fcount = mod->maxComputed+1;
1442: PetscMalloc4(fcount,&fmin,fcount,&fmax,fcount,&fintegral,fcount,&ftmp);
1443: for (i=0; i<fcount; i++) {
1444: fmin[i] = PETSC_MAX_REAL;
1445: fmax[i] = PETSC_MIN_REAL;
1446: fintegral[i] = 0;
1447: }
1448: VecGetDM(cellgeom,&dmCell);
1449: DMPlexGetSimplexOrBoxCells(dmCell,0,&cStart,&cEnd);
1450: VecGetArrayRead(cellgeom,&cgeom);
1451: VecGetArrayRead(X,&x);
1452: DMGetLabel(dm,"vtk",&vtkLabel);
1453: for (c = cStart; c < cEnd; ++c) {
1454: PetscFVCellGeom *cg;
1455: const PetscScalar *cx = NULL;
1456: PetscInt vtkVal = 0;
1458: /* not that these two routines as currently implemented work for any dm with a
1459: * localSection/globalSection */
1460: DMPlexPointLocalRead(dmCell,c,cgeom,&cg);
1461: DMPlexPointGlobalRead(dm,c,x,&cx);
1462: if (vtkLabel) {DMLabelGetValue(vtkLabel,c,&vtkVal);}
1463: if (!vtkVal || !cx) continue; /* ghost, or not a global cell */
1464: for (i=0; i<mod->numCall; i++) {
1465: FunctionalLink flink = mod->functionalCall[i];
1466: (*flink->func)(mod,time,cg->centroid,cx,ftmp,flink->ctx);
1467: }
1468: for (i=0; i<fcount; i++) {
1469: fmin[i] = PetscMin(fmin[i],ftmp[i]);
1470: fmax[i] = PetscMax(fmax[i],ftmp[i]);
1471: fintegral[i] += cg->volume * ftmp[i];
1472: }
1473: }
1474: VecRestoreArrayRead(cellgeom,&cgeom);
1475: VecRestoreArrayRead(X,&x);
1476: DMDestroy(&plex);
1477: MPI_Allreduce(MPI_IN_PLACE,fmin,fcount,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
1478: MPI_Allreduce(MPI_IN_PLACE,fmax,fcount,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1479: MPI_Allreduce(MPI_IN_PLACE,fintegral,fcount,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));
1481: ftablealloc = fcount * 100;
1482: ftableused = 0;
1483: PetscMalloc1(ftablealloc,&ftable);
1484: for (i=0; i<mod->numMonitored; i++) {
1485: size_t countused;
1486: char buffer[256],*p;
1487: FunctionalLink flink = mod->functionalMonitored[i];
1488: PetscInt id = flink->offset;
1489: if (i % 3) {
1490: PetscArraycpy(buffer," ",2);
1491: p = buffer + 2;
1492: } else if (i) {
1493: char newline[] = "\n";
1494: PetscMemcpy(buffer,newline,sizeof(newline)-1);
1495: p = buffer + sizeof(newline) - 1;
1496: } else {
1497: p = buffer;
1498: }
1499: PetscSNPrintfCount(p,sizeof buffer-(p-buffer),"%12s [%10.7g,%10.7g] int %10.7g",&countused,flink->name,(double)fmin[id],(double)fmax[id],(double)fintegral[id]);
1500: countused--;
1501: countused += p - buffer;
1502: if (countused > ftablealloc-ftableused-1) { /* reallocate */
1503: char *ftablenew;
1504: ftablealloc = 2*ftablealloc + countused;
1505: PetscMalloc(ftablealloc,&ftablenew);
1506: PetscArraycpy(ftablenew,ftable,ftableused);
1507: PetscFree(ftable);
1508: ftable = ftablenew;
1509: }
1510: PetscArraycpy(ftable+ftableused,buffer,countused);
1511: ftableused += countused;
1512: ftable[ftableused] = 0;
1513: }
1514: PetscFree4(fmin,fmax,fintegral,ftmp);
1516: PetscPrintf(PetscObjectComm((PetscObject)ts),"% 3D time %8.4g |x| %8.4g %s\n",stepnum,(double)time,(double)xnorm,ftable ? ftable : "");
1517: PetscFree(ftable);
1518: }
1519: if (user->vtkInterval < 1) return(0);
1520: if ((stepnum == -1) ^ (stepnum % user->vtkInterval == 0)) {
1521: if (stepnum == -1) { /* Final time is not multiple of normal time interval, write it anyway */
1522: TSGetStepNumber(ts,&stepnum);
1523: }
1524: PetscSNPrintf(filename,sizeof filename,"%s-%03D.vtu",user->outputBasename,stepnum);
1525: OutputVTK(dm,filename,&viewer);
1526: VecView(X,viewer);
1527: PetscViewerDestroy(&viewer);
1528: }
1529: return(0);
1530: }
1532: static PetscErrorCode initializeTS(DM dm, User user, TS *ts)
1533: {
1537: TSCreate(PetscObjectComm((PetscObject)dm), ts);
1538: TSSetType(*ts, TSSSP);
1539: TSSetDM(*ts, dm);
1540: if (user->vtkmon) {
1541: TSMonitorSet(*ts,MonitorVTK,user,NULL);
1542: }
1543: DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, user);
1544: DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFVM, user);
1545: TSSetMaxTime(*ts,2.0);
1546: TSSetExactFinalTime(*ts,TS_EXACTFINALTIME_STEPOVER);
1547: return(0);
1548: }
1550: static PetscErrorCode adaptToleranceFVM(PetscFV fvm, TS ts, Vec sol, VecTagger refineTag, VecTagger coarsenTag, User user, TS *tsNew, Vec *solNew)
1551: {
1552: DM dm, gradDM, plex, cellDM, adaptedDM = NULL;
1553: Vec cellGeom, faceGeom;
1554: PetscBool isForest, computeGradient;
1555: Vec grad, locGrad, locX, errVec;
1556: PetscInt cStart, cEnd, c, dim, nRefine, nCoarsen;
1557: PetscReal minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2], minInd, maxInd, time;
1558: PetscScalar *errArray;
1559: const PetscScalar *pointVals;
1560: const PetscScalar *pointGrads;
1561: const PetscScalar *pointGeom;
1562: DMLabel adaptLabel = NULL;
1563: IS refineIS, coarsenIS;
1564: PetscErrorCode ierr;
1567: TSGetTime(ts,&time);
1568: VecGetDM(sol, &dm);
1569: DMGetDimension(dm,&dim);
1570: PetscFVGetComputeGradients(fvm,&computeGradient);
1571: PetscFVSetComputeGradients(fvm,PETSC_TRUE);
1572: DMIsForest(dm, &isForest);
1573: DMConvert(dm, DMPLEX, &plex);
1574: DMPlexGetDataFVM(plex, fvm, &cellGeom, &faceGeom, &gradDM);
1575: DMCreateLocalVector(plex,&locX);
1576: DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, faceGeom, cellGeom, NULL);
1577: DMGlobalToLocalBegin(plex, sol, INSERT_VALUES, locX);
1578: DMGlobalToLocalEnd (plex, sol, INSERT_VALUES, locX);
1579: DMCreateGlobalVector(gradDM, &grad);
1580: DMPlexReconstructGradientsFVM(plex, locX, grad);
1581: DMCreateLocalVector(gradDM, &locGrad);
1582: DMGlobalToLocalBegin(gradDM, grad, INSERT_VALUES, locGrad);
1583: DMGlobalToLocalEnd(gradDM, grad, INSERT_VALUES, locGrad);
1584: VecDestroy(&grad);
1585: DMPlexGetSimplexOrBoxCells(plex,0,&cStart,&cEnd);
1586: VecGetArrayRead(locGrad,&pointGrads);
1587: VecGetArrayRead(cellGeom,&pointGeom);
1588: VecGetArrayRead(locX,&pointVals);
1589: VecGetDM(cellGeom,&cellDM);
1590: DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);
1591: VecCreateMPI(PetscObjectComm((PetscObject)plex),cEnd-cStart,PETSC_DETERMINE,&errVec);
1592: VecSetUp(errVec);
1593: VecGetArray(errVec,&errArray);
1594: for (c = cStart; c < cEnd; c++) {
1595: PetscReal errInd = 0.;
1596: PetscScalar *pointGrad;
1597: PetscScalar *pointVal;
1598: PetscFVCellGeom *cg;
1600: DMPlexPointLocalRead(gradDM,c,pointGrads,&pointGrad);
1601: DMPlexPointLocalRead(cellDM,c,pointGeom,&cg);
1602: DMPlexPointLocalRead(plex,c,pointVals,&pointVal);
1604: (user->model->errorIndicator)(dim,cg->volume,user->model->physics->dof,pointVal,pointGrad,&errInd,user->model->errorCtx);
1605: errArray[c-cStart] = errInd;
1606: minMaxInd[0] = PetscMin(minMaxInd[0],errInd);
1607: minMaxInd[1] = PetscMax(minMaxInd[1],errInd);
1608: }
1609: VecRestoreArray(errVec,&errArray);
1610: VecRestoreArrayRead(locX,&pointVals);
1611: VecRestoreArrayRead(cellGeom,&pointGeom);
1612: VecRestoreArrayRead(locGrad,&pointGrads);
1613: VecDestroy(&locGrad);
1614: VecDestroy(&locX);
1615: DMDestroy(&plex);
1617: VecTaggerComputeIS(refineTag,errVec,&refineIS);
1618: VecTaggerComputeIS(coarsenTag,errVec,&coarsenIS);
1619: ISGetSize(refineIS,&nRefine);
1620: ISGetSize(coarsenIS,&nCoarsen);
1621: if (nRefine) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_REFINE,refineIS);}
1622: if (nCoarsen) {DMLabelSetStratumIS(adaptLabel,DM_ADAPT_COARSEN,coarsenIS);}
1623: ISDestroy(&coarsenIS);
1624: ISDestroy(&refineIS);
1625: VecDestroy(&errVec);
1627: PetscFVSetComputeGradients(fvm,computeGradient);
1628: minMaxInd[1] = -minMaxInd[1];
1629: MPI_Allreduce(minMaxInd,minMaxIndGlobal,2,MPIU_REAL,MPI_MIN,PetscObjectComm((PetscObject)dm));
1630: minInd = minMaxIndGlobal[0];
1631: maxInd = -minMaxIndGlobal[1];
1632: PetscInfo2(ts, "error indicator range (%E, %E)\n", minInd, maxInd);
1633: if (nRefine || nCoarsen) { /* at least one cell is over the refinement threshold */
1634: DMAdaptLabel(dm,adaptLabel,&adaptedDM);
1635: }
1636: DMLabelDestroy(&adaptLabel);
1637: if (adaptedDM) {
1638: PetscInfo2(ts, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nRefine, nCoarsen);
1639: if (tsNew) {initializeTS(adaptedDM, user, tsNew);}
1640: if (solNew) {
1641: DMCreateGlobalVector(adaptedDM, solNew);
1642: PetscObjectSetName((PetscObject) *solNew, "solution");
1643: DMForestTransferVec(dm, sol, adaptedDM, *solNew, PETSC_TRUE, time);
1644: }
1645: if (isForest) {DMForestSetAdaptivityForest(adaptedDM,NULL);} /* clear internal references to the previous dm */
1646: DMDestroy(&adaptedDM);
1647: } else {
1648: if (tsNew) *tsNew = NULL;
1649: if (solNew) *solNew = NULL;
1650: }
1651: return(0);
1652: }
1654: int main(int argc, char **argv)
1655: {
1656: MPI_Comm comm;
1657: PetscDS prob;
1658: PetscFV fvm;
1659: PetscLimiter limiter = NULL, noneLimiter = NULL;
1660: User user;
1661: Model mod;
1662: Physics phys;
1663: DM dm, plex;
1664: PetscReal ftime, cfl, dt, minRadius;
1665: PetscInt dim, nsteps;
1666: TS ts;
1667: TSConvergedReason reason;
1668: Vec X;
1669: PetscViewer viewer;
1670: PetscBool simplex = PETSC_FALSE, vtkCellGeom, splitFaces, useAMR;
1671: PetscInt overlap, adaptInterval;
1672: char filename[PETSC_MAX_PATH_LEN] = "";
1673: char physname[256] = "advect";
1674: VecTagger refineTag = NULL, coarsenTag = NULL;
1675: PetscErrorCode ierr;
1677: PetscInitialize(&argc, &argv, (char*) 0, help);if (ierr) return ierr;
1678: comm = PETSC_COMM_WORLD;
1680: PetscNew(&user);
1681: PetscNew(&user->model);
1682: PetscNew(&user->model->physics);
1683: mod = user->model;
1684: phys = mod->physics;
1685: mod->comm = comm;
1686: useAMR = PETSC_FALSE;
1687: adaptInterval = 1;
1689: /* Register physical models to be available on the command line */
1690: PetscFunctionListAdd(&PhysicsList,"advect" ,PhysicsCreate_Advect);
1691: PetscFunctionListAdd(&PhysicsList,"sw" ,PhysicsCreate_SW);
1692: PetscFunctionListAdd(&PhysicsList,"euler" ,PhysicsCreate_Euler);
1695: PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Mesh Options","");
1696: {
1697: cfl = 0.9 * 4; /* default SSPRKS2 with s=5 stages is stable for CFL number s-1 */
1698: PetscOptionsReal("-ufv_cfl","CFL number per step","",cfl,&cfl,NULL);
1699: PetscOptionsString("-f","Exodus.II filename to read","",filename,filename,sizeof(filename),NULL);
1700: PetscOptionsBool("-simplex","Flag to use a simplex mesh","",simplex,&simplex,NULL);
1701: splitFaces = PETSC_FALSE;
1702: PetscOptionsBool("-ufv_split_faces","Split faces between cell sets","",splitFaces,&splitFaces,NULL);
1703: overlap = 1;
1704: PetscOptionsInt("-ufv_mesh_overlap","Number of cells to overlap partitions","",overlap,&overlap,NULL);
1705: user->vtkInterval = 1;
1706: PetscOptionsInt("-ufv_vtk_interval","VTK output interval (0 to disable)","",user->vtkInterval,&user->vtkInterval,NULL);
1707: user->vtkmon = PETSC_TRUE;
1708: PetscOptionsBool("-ufv_vtk_monitor","Use VTKMonitor routine","",user->vtkmon,&user->vtkmon,NULL);
1709: vtkCellGeom = PETSC_FALSE;
1710: PetscStrcpy(user->outputBasename, "ex11");
1711: PetscOptionsString("-ufv_vtk_basename","VTK output basename","",user->outputBasename,user->outputBasename,sizeof(user->outputBasename),NULL);
1712: PetscOptionsBool("-ufv_vtk_cellgeom","Write cell geometry (for debugging)","",vtkCellGeom,&vtkCellGeom,NULL);
1713: PetscOptionsBool("-ufv_use_amr","use local adaptive mesh refinement","",useAMR,&useAMR,NULL);
1714: PetscOptionsInt("-ufv_adapt_interval","time steps between AMR","",adaptInterval,&adaptInterval,NULL);
1715: }
1716: PetscOptionsEnd();
1718: if (useAMR) {
1719: VecTaggerBox refineBox, coarsenBox;
1721: refineBox.min = refineBox.max = PETSC_MAX_REAL;
1722: coarsenBox.min = coarsenBox.max = PETSC_MIN_REAL;
1724: VecTaggerCreate(comm,&refineTag);
1725: PetscObjectSetOptionsPrefix((PetscObject)refineTag,"refine_");
1726: VecTaggerSetType(refineTag,VECTAGGERABSOLUTE);
1727: VecTaggerAbsoluteSetBox(refineTag,&refineBox);
1728: VecTaggerSetFromOptions(refineTag);
1729: VecTaggerSetUp(refineTag);
1730: PetscObjectViewFromOptions((PetscObject)refineTag,NULL,"-tag_view");
1732: VecTaggerCreate(comm,&coarsenTag);
1733: PetscObjectSetOptionsPrefix((PetscObject)coarsenTag,"coarsen_");
1734: VecTaggerSetType(coarsenTag,VECTAGGERABSOLUTE);
1735: VecTaggerAbsoluteSetBox(coarsenTag,&coarsenBox);
1736: VecTaggerSetFromOptions(coarsenTag);
1737: VecTaggerSetUp(coarsenTag);
1738: PetscObjectViewFromOptions((PetscObject)coarsenTag,NULL,"-tag_view");
1739: }
1741: PetscOptionsBegin(comm,NULL,"Unstructured Finite Volume Physics Options","");
1742: {
1743: PetscErrorCode (*physcreate)(Model,Physics,PetscOptionItems*);
1744: PetscOptionsFList("-physics","Physics module to solve","",PhysicsList,physname,physname,sizeof physname,NULL);
1745: PetscFunctionListFind(PhysicsList,physname,&physcreate);
1746: PetscMemzero(phys,sizeof(struct _n_Physics));
1747: (*physcreate)(mod,phys,PetscOptionsObject);
1748: /* Count number of fields and dofs */
1749: for (phys->nfields=0,phys->dof=0; phys->field_desc[phys->nfields].name; phys->nfields++) phys->dof += phys->field_desc[phys->nfields].dof;
1750: if (phys->dof <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set dof",physname);
1751: ModelFunctionalSetFromOptions(mod,PetscOptionsObject);
1752: }
1753: PetscOptionsEnd();
1755: /* Create mesh */
1756: {
1757: size_t len,i;
1758: for (i = 0; i < DIM; i++) { mod->bounds[2*i] = 0.; mod->bounds[2*i+1] = 1.;};
1759: PetscStrlen(filename,&len);
1760: dim = DIM;
1761: if (!len) { /* a null name means just do a hex box */
1762: PetscInt cells[3] = {1, 1, 1}; /* coarse mesh is one cell; refine from there */
1763: PetscBool flg1, flg2, skew = PETSC_FALSE;
1764: PetscInt nret1 = DIM;
1765: PetscInt nret2 = 2*DIM;
1766: PetscOptionsBegin(comm,NULL,"Rectangular mesh options","");
1767: PetscOptionsIntArray("-grid_size","number of cells in each direction","",cells,&nret1,&flg1);
1768: PetscOptionsRealArray("-grid_bounds","bounds of the mesh in each direction (i.e., x_min,x_max,y_min,y_max","",mod->bounds,&nret2,&flg2);
1769: PetscOptionsBool("-grid_skew_60","Skew grid for 60 degree shock mesh","",skew,&skew,NULL);
1770: PetscOptionsEnd();
1771: if (flg1) {
1772: dim = nret1;
1773: if (dim != DIM) SETERRQ1(comm,PETSC_ERR_ARG_SIZ,"Dim wrong size %D in -grid_size",dim);
1774: }
1775: DMPlexCreateBoxMesh(comm, dim, simplex, cells, NULL, NULL, mod->bcs, PETSC_TRUE, &dm);
1776: if (flg2) {
1777: PetscInt dimEmbed, i;
1778: PetscInt nCoords;
1779: PetscScalar *coords;
1780: Vec coordinates;
1782: DMGetCoordinatesLocal(dm,&coordinates);
1783: DMGetCoordinateDim(dm,&dimEmbed);
1784: VecGetLocalSize(coordinates,&nCoords);
1785: if (nCoords % dimEmbed) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Coordinate vector the wrong size");
1786: VecGetArray(coordinates,&coords);
1787: for (i = 0; i < nCoords; i += dimEmbed) {
1788: PetscInt j;
1790: PetscScalar *coord = &coords[i];
1791: for (j = 0; j < dimEmbed; j++) {
1792: coord[j] = mod->bounds[2 * j] + coord[j] * (mod->bounds[2 * j + 1] - mod->bounds[2 * j]);
1793: if (dim==2 && cells[1]==1 && j==0 && skew) {
1794: if (cells[0]==2 && i==8) {
1795: coord[j] = .57735026918963; /* hack to get 60 deg skewed mesh */
1796: }
1797: else if (cells[0]==3) {
1798: if (i==2 || i==10) coord[j] = mod->bounds[1]/4.;
1799: else if (i==4) coord[j] = mod->bounds[1]/2.;
1800: else if (i==12) coord[j] = 1.57735026918963*mod->bounds[1]/2.;
1801: }
1802: }
1803: }
1804: }
1805: VecRestoreArray(coordinates,&coords);
1806: DMSetCoordinatesLocal(dm,coordinates);
1807: }
1808: } else {
1809: DMPlexCreateFromFile(comm, filename, PETSC_TRUE, &dm);
1810: }
1811: }
1812: DMViewFromOptions(dm, NULL, "-orig_dm_view");
1813: DMGetDimension(dm, &dim);
1815: /* set up BCs, functions, tags */
1816: DMCreateLabel(dm, "Face Sets");
1818: mod->errorIndicator = ErrorIndicator_Simple;
1820: {
1821: DM dmDist;
1823: DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE);
1824: DMPlexDistribute(dm, overlap, NULL, &dmDist);
1825: if (dmDist) {
1826: DMDestroy(&dm);
1827: dm = dmDist;
1828: }
1829: }
1831: DMSetFromOptions(dm);
1833: {
1834: DM gdm;
1836: DMPlexConstructGhostCells(dm, NULL, NULL, &gdm);
1837: DMDestroy(&dm);
1838: dm = gdm;
1839: DMViewFromOptions(dm, NULL, "-dm_view");
1840: }
1841: if (splitFaces) {ConstructCellBoundary(dm, user);}
1842: SplitFaces(&dm, "split faces", user);
1844: PetscFVCreate(comm, &fvm);
1845: PetscFVSetFromOptions(fvm);
1846: PetscFVSetNumComponents(fvm, phys->dof);
1847: PetscFVSetSpatialDimension(fvm, dim);
1848: PetscObjectSetName((PetscObject) fvm,"");
1849: {
1850: PetscInt f, dof;
1851: for (f=0,dof=0; f < phys->nfields; f++) {
1852: PetscInt newDof = phys->field_desc[f].dof;
1854: if (newDof == 1) {
1855: PetscFVSetComponentName(fvm,dof,phys->field_desc[f].name);
1856: }
1857: else {
1858: PetscInt j;
1860: for (j = 0; j < newDof; j++) {
1861: char compName[256] = "Unknown";
1863: PetscSNPrintf(compName,sizeof(compName),"%s_%d",phys->field_desc[f].name,j);
1864: PetscFVSetComponentName(fvm,dof+j,compName);
1865: }
1866: }
1867: dof += newDof;
1868: }
1869: }
1870: /* FV is now structured with one field having all physics as components */
1871: DMAddField(dm, NULL, (PetscObject) fvm);
1872: DMCreateDS(dm);
1873: DMGetDS(dm, &prob);
1874: PetscDSSetRiemannSolver(prob, 0, user->model->physics->riemann);
1875: PetscDSSetContext(prob, 0, user->model->physics);
1876: (*mod->setupbc)(prob,phys);
1877: PetscDSSetFromOptions(prob);
1878: {
1879: char convType[256];
1880: PetscBool flg;
1882: PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
1883: PetscOptionsFList("-dm_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);
1884: PetscOptionsEnd();
1885: if (flg) {
1886: DM dmConv;
1888: DMConvert(dm,convType,&dmConv);
1889: if (dmConv) {
1890: DMViewFromOptions(dmConv, NULL, "-dm_conv_view");
1891: DMDestroy(&dm);
1892: dm = dmConv;
1893: DMSetFromOptions(dm);
1894: }
1895: }
1896: }
1898: initializeTS(dm, user, &ts);
1900: DMCreateGlobalVector(dm, &X);
1901: PetscObjectSetName((PetscObject) X, "solution");
1902: SetInitialCondition(dm, X, user);
1903: if (useAMR) {
1904: PetscInt adaptIter;
1906: /* use no limiting when reconstructing gradients for adaptivity */
1907: PetscFVGetLimiter(fvm, &limiter);
1908: PetscObjectReference((PetscObject) limiter);
1909: PetscLimiterCreate(PetscObjectComm((PetscObject) fvm), &noneLimiter);
1910: PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE);
1912: PetscFVSetLimiter(fvm, noneLimiter);
1913: for (adaptIter = 0; ; ++adaptIter) {
1914: PetscLogDouble bytes;
1915: TS tsNew = NULL;
1917: PetscMemoryGetCurrentUsage(&bytes);
1918: PetscInfo2(ts, "refinement loop %D: memory used %g\n", adaptIter, bytes);
1919: DMViewFromOptions(dm, NULL, "-initial_dm_view");
1920: VecViewFromOptions(X, NULL, "-initial_vec_view");
1921: #if 0
1922: if (viewInitial) {
1923: PetscViewer viewer;
1924: char buf[256];
1925: PetscBool isHDF5, isVTK;
1927: PetscViewerCreate(comm,&viewer);
1928: PetscViewerSetType(viewer,PETSCVIEWERVTK);
1929: PetscViewerSetOptionsPrefix(viewer,"initial_");
1930: PetscViewerSetFromOptions(viewer);
1931: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&isHDF5);
1932: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isVTK);
1933: if (isHDF5) {
1934: PetscSNPrintf(buf, 256, "ex11-initial-%d.h5", adaptIter);
1935: } else if (isVTK) {
1936: PetscSNPrintf(buf, 256, "ex11-initial-%d.vtu", adaptIter);
1937: PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);
1938: }
1939: PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
1940: PetscViewerFileSetName(viewer,buf);
1941: if (isHDF5) {
1942: DMView(dm,viewer);
1943: PetscViewerFileSetMode(viewer,FILE_MODE_UPDATE);
1944: }
1945: VecView(X,viewer);
1946: PetscViewerDestroy(&viewer);
1947: }
1948: #endif
1950: adaptToleranceFVM(fvm, ts, X, refineTag, coarsenTag, user, &tsNew, NULL);
1951: if (!tsNew) {
1952: break;
1953: } else {
1954: DMDestroy(&dm);
1955: VecDestroy(&X);
1956: TSDestroy(&ts);
1957: ts = tsNew;
1958: TSGetDM(ts,&dm);
1959: PetscObjectReference((PetscObject)dm);
1960: DMCreateGlobalVector(dm,&X);
1961: PetscObjectSetName((PetscObject) X, "solution");
1962: SetInitialCondition(dm, X, user);
1963: }
1964: }
1965: /* restore original limiter */
1966: PetscFVSetLimiter(fvm, limiter);
1967: }
1969: DMConvert(dm, DMPLEX, &plex);
1970: if (vtkCellGeom) {
1971: DM dmCell;
1972: Vec cellgeom, partition;
1974: DMPlexGetGeometryFVM(plex, NULL, &cellgeom, NULL);
1975: OutputVTK(dm, "ex11-cellgeom.vtk", &viewer);
1976: VecView(cellgeom, viewer);
1977: PetscViewerDestroy(&viewer);
1978: CreatePartitionVec(dm, &dmCell, &partition);
1979: OutputVTK(dmCell, "ex11-partition.vtk", &viewer);
1980: VecView(partition, viewer);
1981: PetscViewerDestroy(&viewer);
1982: VecDestroy(&partition);
1983: DMDestroy(&dmCell);
1984: }
1985: /* collect max maxspeed from all processes -- todo */
1986: DMPlexGetGeometryFVM(plex, NULL, NULL, &minRadius);
1987: DMDestroy(&plex);
1988: MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
1989: if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
1990: dt = cfl * minRadius / mod->maxspeed;
1991: TSSetTimeStep(ts,dt);
1992: TSSetFromOptions(ts);
1993: if (!useAMR) {
1994: TSSolve(ts,X);
1995: TSGetSolveTime(ts,&ftime);
1996: TSGetStepNumber(ts,&nsteps);
1997: } else {
1998: PetscReal finalTime;
1999: PetscInt adaptIter;
2000: TS tsNew = NULL;
2001: Vec solNew = NULL;
2003: TSGetMaxTime(ts,&finalTime);
2004: TSSetMaxSteps(ts,adaptInterval);
2005: TSSolve(ts,X);
2006: TSGetSolveTime(ts,&ftime);
2007: TSGetStepNumber(ts,&nsteps);
2008: for (adaptIter = 0;ftime < finalTime;adaptIter++) {
2009: PetscLogDouble bytes;
2011: PetscMemoryGetCurrentUsage(&bytes);
2012: PetscInfo2(ts, "AMR time step loop %D: memory used %g\n", adaptIter, bytes);
2013: PetscFVSetLimiter(fvm,noneLimiter);
2014: adaptToleranceFVM(fvm,ts,X,refineTag,coarsenTag,user,&tsNew,&solNew);
2015: PetscFVSetLimiter(fvm,limiter);
2016: if (tsNew) {
2017: PetscInfo(ts, "AMR used\n");
2018: DMDestroy(&dm);
2019: VecDestroy(&X);
2020: TSDestroy(&ts);
2021: ts = tsNew;
2022: X = solNew;
2023: TSSetFromOptions(ts);
2024: VecGetDM(X,&dm);
2025: PetscObjectReference((PetscObject)dm);
2026: DMConvert(dm, DMPLEX, &plex);
2027: DMPlexGetGeometryFVM(dm, NULL, NULL, &minRadius);
2028: DMDestroy(&plex);
2029: MPI_Allreduce(&phys->maxspeed,&mod->maxspeed,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)ts));
2030: if (mod->maxspeed <= 0) SETERRQ1(comm,PETSC_ERR_ARG_WRONGSTATE,"Physics '%s' did not set maxspeed",physname);
2031: dt = cfl * minRadius / mod->maxspeed;
2032: TSSetStepNumber(ts,nsteps);
2033: TSSetTime(ts,ftime);
2034: TSSetTimeStep(ts,dt);
2035: } else {
2036: PetscInfo(ts, "AMR not used\n");
2037: }
2038: user->monitorStepOffset = nsteps;
2039: TSSetMaxSteps(ts,nsteps+adaptInterval);
2040: TSSolve(ts,X);
2041: TSGetSolveTime(ts,&ftime);
2042: TSGetStepNumber(ts,&nsteps);
2043: }
2044: }
2045: TSGetConvergedReason(ts,&reason);
2046: PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,nsteps);
2047: TSDestroy(&ts);
2049: VecTaggerDestroy(&refineTag);
2050: VecTaggerDestroy(&coarsenTag);
2051: PetscFunctionListDestroy(&PhysicsList);
2052: PetscFunctionListDestroy(&PhysicsRiemannList_SW);
2053: FunctionalLinkDestroy(&user->model->functionalRegistry);
2054: PetscFree(user->model->functionalMonitored);
2055: PetscFree(user->model->functionalCall);
2056: PetscFree(user->model->physics->data);
2057: PetscFree(user->model->physics);
2058: PetscFree(user->model);
2059: PetscFree(user);
2060: VecDestroy(&X);
2061: PetscLimiterDestroy(&limiter);
2062: PetscLimiterDestroy(&noneLimiter);
2063: PetscFVDestroy(&fvm);
2064: DMDestroy(&dm);
2065: PetscFinalize();
2066: return ierr;
2067: }
2069: /* Godunov fluxs */
2070: PetscScalar cvmgp_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2071: {
2072: /* System generated locals */
2073: PetscScalar ret_val;
2075: if (PetscRealPart(*test) > 0.) {
2076: goto L10;
2077: }
2078: ret_val = *b;
2079: return ret_val;
2080: L10:
2081: ret_val = *a;
2082: return ret_val;
2083: } /* cvmgp_ */
2085: PetscScalar cvmgm_(PetscScalar *a, PetscScalar *b, PetscScalar *test)
2086: {
2087: /* System generated locals */
2088: PetscScalar ret_val;
2090: if (PetscRealPart(*test) < 0.) {
2091: goto L10;
2092: }
2093: ret_val = *b;
2094: return ret_val;
2095: L10:
2096: ret_val = *a;
2097: return ret_val;
2098: } /* cvmgm_ */
2100: int riem1mdt( PetscScalar *gaml, PetscScalar *gamr, PetscScalar *rl, PetscScalar *pl,
2101: PetscScalar *uxl, PetscScalar *rr, PetscScalar *pr,
2102: PetscScalar *uxr, PetscScalar *rstarl, PetscScalar *rstarr, PetscScalar *
2103: pstar, PetscScalar *ustar)
2104: {
2105: /* Initialized data */
2107: static PetscScalar smallp = 1e-8;
2109: /* System generated locals */
2110: int i__1;
2111: PetscScalar d__1, d__2;
2113: /* Local variables */
2114: static int i0;
2115: static PetscScalar cl, cr, wl, zl, wr, zr, pst, durl, skpr1, skpr2;
2116: static int iwave;
2117: static PetscScalar gascl4, gascr4, cstarl, dpstar, cstarr;
2118: /* static PetscScalar csqrl, csqrr, gascl1, gascl2, gascl3, gascr1, gascr2, gascr3; */
2119: static int iterno;
2120: static PetscScalar ustarl, ustarr, rarepr1, rarepr2;
2124: /* gascl1 = *gaml - 1.; */
2125: /* gascl2 = (*gaml + 1.) * .5; */
2126: /* gascl3 = gascl2 / *gaml; */
2127: gascl4 = 1. / (*gaml - 1.);
2129: /* gascr1 = *gamr - 1.; */
2130: /* gascr2 = (*gamr + 1.) * .5; */
2131: /* gascr3 = gascr2 / *gamr; */
2132: gascr4 = 1. / (*gamr - 1.);
2133: iterno = 10;
2134: /* find pstar: */
2135: cl = PetscSqrtScalar(*gaml * *pl / *rl);
2136: cr = PetscSqrtScalar(*gamr * *pr / *rr);
2137: wl = *rl * cl;
2138: wr = *rr * cr;
2139: /* csqrl = wl * wl; */
2140: /* csqrr = wr * wr; */
2141: *pstar = (wl * *pr + wr * *pl) / (wl + wr);
2142: *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2143: pst = *pl / *pr;
2144: skpr1 = cr * (pst - 1.) * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2145: d__1 = (*gamr - 1.) / (*gamr * 2.);
2146: rarepr2 = gascr4 * 2. * cr * (1. - PetscPowScalar(pst, d__1));
2147: pst = *pr / *pl;
2148: skpr2 = cl * (pst - 1.) * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2149: d__1 = (*gaml - 1.) / (*gaml * 2.);
2150: rarepr1 = gascl4 * 2. * cl * (1. - PetscPowScalar(pst, d__1));
2151: durl = *uxr - *uxl;
2152: if (PetscRealPart(*pr) < PetscRealPart(*pl)) {
2153: if (PetscRealPart(durl) >= PetscRealPart(rarepr1)) {
2154: iwave = 100;
2155: } else if (PetscRealPart(durl) <= PetscRealPart(-skpr1)) {
2156: iwave = 300;
2157: } else {
2158: iwave = 400;
2159: }
2160: } else {
2161: if (PetscRealPart(durl) >= PetscRealPart(rarepr2)) {
2162: iwave = 100;
2163: } else if (PetscRealPart(durl) <= PetscRealPart(-skpr2)) {
2164: iwave = 300;
2165: } else {
2166: iwave = 200;
2167: }
2168: }
2169: if (iwave == 100) {
2170: /* 1-wave: rarefaction wave, 3-wave: rarefaction wave */
2171: /* case (100) */
2172: i__1 = iterno;
2173: for (i0 = 1; i0 <= i__1; ++i0) {
2174: d__1 = *pstar / *pl;
2175: d__2 = 1. / *gaml;
2176: *rstarl = *rl * PetscPowScalar(d__1, d__2);
2177: cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2178: ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2179: zl = *rstarl * cstarl;
2180: d__1 = *pstar / *pr;
2181: d__2 = 1. / *gamr;
2182: *rstarr = *rr * PetscPowScalar(d__1, d__2);
2183: cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2184: ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2185: zr = *rstarr * cstarr;
2186: dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2187: *pstar -= dpstar;
2188: *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2189: if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2190: #if 0
2191: break;
2192: #endif
2193: }
2194: }
2195: /* 1-wave: shock wave, 3-wave: rarefaction wave */
2196: } else if (iwave == 200) {
2197: /* case (200) */
2198: i__1 = iterno;
2199: for (i0 = 1; i0 <= i__1; ++i0) {
2200: pst = *pstar / *pl;
2201: ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2202: zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2203: d__1 = *pstar / *pr;
2204: d__2 = 1. / *gamr;
2205: *rstarr = *rr * PetscPowScalar(d__1, d__2);
2206: cstarr = PetscSqrtScalar(*gamr * *pstar / *rstarr);
2207: zr = *rstarr * cstarr;
2208: ustarr = *uxr + gascr4 * 2. * (cstarr - cr);
2209: dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2210: *pstar -= dpstar;
2211: *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2212: if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2213: #if 0
2214: break;
2215: #endif
2216: }
2217: }
2218: /* 1-wave: shock wave, 3-wave: shock */
2219: } else if (iwave == 300) {
2220: /* case (300) */
2221: i__1 = iterno;
2222: for (i0 = 1; i0 <= i__1; ++i0) {
2223: pst = *pstar / *pl;
2224: ustarl = *uxl - (pst - 1.) * cl * PetscSqrtScalar(2. / (*gaml * (*gaml - 1. + (*gaml + 1.) * pst)));
2225: zl = *pl / cl * PetscSqrtScalar(*gaml * 2. * (*gaml - 1. + (*gaml + 1.) * pst)) * (*gaml - 1. + (*gaml + 1.) * pst) / (*gaml * 3. - 1. + (*gaml + 1.) * pst);
2226: pst = *pstar / *pr;
2227: ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2228: zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2229: dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2230: *pstar -= dpstar;
2231: *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2232: if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2233: #if 0
2234: break;
2235: #endif
2236: }
2237: }
2238: /* 1-wave: rarefaction wave, 3-wave: shock */
2239: } else if (iwave == 400) {
2240: /* case (400) */
2241: i__1 = iterno;
2242: for (i0 = 1; i0 <= i__1; ++i0) {
2243: d__1 = *pstar / *pl;
2244: d__2 = 1. / *gaml;
2245: *rstarl = *rl * PetscPowScalar(d__1, d__2);
2246: cstarl = PetscSqrtScalar(*gaml * *pstar / *rstarl);
2247: ustarl = *uxl - gascl4 * 2. * (cstarl - cl);
2248: zl = *rstarl * cstarl;
2249: pst = *pstar / *pr;
2250: ustarr = *uxr + (pst - 1.) * cr * PetscSqrtScalar(2. / (*gamr * (*gamr - 1. + (*gamr + 1.) * pst)));
2251: zr = *pr / cr * PetscSqrtScalar(*gamr * 2. * (*gamr - 1. + (*gamr + 1.) * pst)) * (*gamr - 1. + (*gamr + 1.) * pst) / (*gamr * 3. - 1. + (*gamr + 1.) * pst);
2252: dpstar = zl * zr * (ustarr - ustarl) / (zl + zr);
2253: *pstar -= dpstar;
2254: *pstar = PetscMax(PetscRealPart(*pstar),PetscRealPart(smallp));
2255: if (PetscAbsScalar(dpstar) / PetscRealPart(*pstar) <= 1e-8) {
2256: #if 0
2257: break;
2258: #endif
2259: }
2260: }
2261: }
2263: *ustar = (zl * ustarr + zr * ustarl) / (zl + zr);
2264: if (PetscRealPart(*pstar) > PetscRealPart(*pl)) {
2265: pst = *pstar / *pl;
2266: *rstarl = ((*gaml + 1.) * pst + *gaml - 1.) / ((*gaml - 1.) * pst + *
2267: gaml + 1.) * *rl;
2268: }
2269: if (PetscRealPart(*pstar) > PetscRealPart(*pr)) {
2270: pst = *pstar / *pr;
2271: *rstarr = ((*gamr + 1.) * pst + *gamr - 1.) / ((*gamr - 1.) * pst + *
2272: gamr + 1.) * *rr;
2273: }
2274: return iwave;
2275: }
2277: PetscScalar sign(PetscScalar x)
2278: {
2279: if (PetscRealPart(x) > 0) return 1.0;
2280: if (PetscRealPart(x) < 0) return -1.0;
2281: return 0.0;
2282: }
2283: /* Riemann Solver */
2284: /* -------------------------------------------------------------------- */
2285: int riemannsolver(PetscScalar *xcen, PetscScalar *xp,
2286: PetscScalar *dtt, PetscScalar *rl, PetscScalar *uxl, PetscScalar *pl,
2287: PetscScalar *utl, PetscScalar *ubl, PetscScalar *gaml, PetscScalar *rho1l,
2288: PetscScalar *rr, PetscScalar *uxr, PetscScalar *pr, PetscScalar *utr,
2289: PetscScalar *ubr, PetscScalar *gamr, PetscScalar *rho1r, PetscScalar *rx,
2290: PetscScalar *uxm, PetscScalar *px, PetscScalar *utx, PetscScalar *ubx,
2291: PetscScalar *gam, PetscScalar *rho1)
2292: {
2293: /* System generated locals */
2294: PetscScalar d__1, d__2;
2296: /* Local variables */
2297: static PetscScalar s, c0, p0, r0, u0, w0, x0, x2, ri, cx, sgn0, wsp0, gasc1, gasc2, gasc3, gasc4;
2298: static PetscScalar cstar, pstar, rstar, ustar, xstar, wspst, ushock, streng, rstarl, rstarr, rstars;
2299: int iwave;
2301: if (*rl == *rr && *pr == *pl && *uxl == *uxr && *gaml == *gamr) {
2302: *rx = *rl;
2303: *px = *pl;
2304: *uxm = *uxl;
2305: *gam = *gaml;
2306: x2 = *xcen + *uxm * *dtt;
2308: if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2309: *utx = *utr;
2310: *ubx = *ubr;
2311: *rho1 = *rho1r;
2312: } else {
2313: *utx = *utl;
2314: *ubx = *ubl;
2315: *rho1 = *rho1l;
2316: }
2317: return 0;
2318: }
2319: iwave = riem1mdt(gaml, gamr, rl, pl, uxl, rr, pr, uxr, &rstarl, &rstarr, &pstar, &ustar);
2321: x2 = *xcen + ustar * *dtt;
2322: d__1 = *xp - x2;
2323: sgn0 = sign(d__1);
2324: /* x is in 3-wave if sgn0 = 1 */
2325: /* x is in 1-wave if sgn0 = -1 */
2326: r0 = cvmgm_(rl, rr, &sgn0);
2327: p0 = cvmgm_(pl, pr, &sgn0);
2328: u0 = cvmgm_(uxl, uxr, &sgn0);
2329: *gam = cvmgm_(gaml, gamr, &sgn0);
2330: gasc1 = *gam - 1.;
2331: gasc2 = (*gam + 1.) * .5;
2332: gasc3 = gasc2 / *gam;
2333: gasc4 = 1. / (*gam - 1.);
2334: c0 = PetscSqrtScalar(*gam * p0 / r0);
2335: streng = pstar - p0;
2336: w0 = *gam * r0 * p0 * (gasc3 * streng / p0 + 1.);
2337: rstars = r0 / (1. - r0 * streng / w0);
2338: d__1 = p0 / pstar;
2339: d__2 = -1. / *gam;
2340: rstarr = r0 * PetscPowScalar(d__1, d__2);
2341: rstar = cvmgm_(&rstarr, &rstars, &streng);
2342: w0 = PetscSqrtScalar(w0);
2343: cstar = PetscSqrtScalar(*gam * pstar / rstar);
2344: wsp0 = u0 + sgn0 * c0;
2345: wspst = ustar + sgn0 * cstar;
2346: ushock = ustar + sgn0 * w0 / rstar;
2347: wspst = cvmgp_(&ushock, &wspst, &streng);
2348: wsp0 = cvmgp_(&ushock, &wsp0, &streng);
2349: x0 = *xcen + wsp0 * *dtt;
2350: xstar = *xcen + wspst * *dtt;
2351: /* using gas formula to evaluate rarefaction wave */
2352: /* ri : reiman invariant */
2353: ri = u0 - sgn0 * 2. * gasc4 * c0;
2354: cx = sgn0 * .5 * gasc1 / gasc2 * ((*xp - *xcen) / *dtt - ri);
2355: *uxm = ri + sgn0 * 2. * gasc4 * cx;
2356: s = p0 / PetscPowScalar(r0, *gam);
2357: d__1 = cx * cx / (*gam * s);
2358: *rx = PetscPowScalar(d__1, gasc4);
2359: *px = cx * cx * *rx / *gam;
2360: d__1 = sgn0 * (x0 - *xp);
2361: *rx = cvmgp_(rx, &r0, &d__1);
2362: d__1 = sgn0 * (x0 - *xp);
2363: *px = cvmgp_(px, &p0, &d__1);
2364: d__1 = sgn0 * (x0 - *xp);
2365: *uxm = cvmgp_(uxm, &u0, &d__1);
2366: d__1 = sgn0 * (xstar - *xp);
2367: *rx = cvmgm_(rx, &rstar, &d__1);
2368: d__1 = sgn0 * (xstar - *xp);
2369: *px = cvmgm_(px, &pstar, &d__1);
2370: d__1 = sgn0 * (xstar - *xp);
2371: *uxm = cvmgm_(uxm, &ustar, &d__1);
2372: if (PetscRealPart(*xp) >= PetscRealPart(x2)) {
2373: *utx = *utr;
2374: *ubx = *ubr;
2375: *rho1 = *rho1r;
2376: } else {
2377: *utx = *utl;
2378: *ubx = *ubl;
2379: *rho1 = *rho1l;
2380: }
2381: return iwave;
2382: }
2383: int godunovflux( const PetscScalar *ul, const PetscScalar *ur,
2384: PetscScalar *flux, const PetscReal *nn, const int *ndim,
2385: const PetscReal *gamma)
2386: {
2387: /* System generated locals */
2388: int i__1,iwave;
2389: PetscScalar d__1, d__2, d__3;
2391: /* Local variables */
2392: static int k;
2393: static PetscScalar bn[3], fn, ft, tg[3], pl, rl, pm, pr, rr, xp, ubl, ubm,
2394: ubr, dtt, unm, tmp, utl, utm, uxl, utr, uxr, gaml, gamm, gamr,
2395: xcen, rhom, rho1l, rho1m, rho1r;
2397: /* Function Body */
2398: xcen = 0.;
2399: xp = 0.;
2400: i__1 = *ndim;
2401: for (k = 1; k <= i__1; ++k) {
2402: tg[k - 1] = 0.;
2403: bn[k - 1] = 0.;
2404: }
2405: dtt = 1.;
2406: if (*ndim == 3) {
2407: if (nn[0] == 0. && nn[1] == 0.) {
2408: tg[0] = 1.;
2409: } else {
2410: tg[0] = -nn[1];
2411: tg[1] = nn[0];
2412: }
2413: /* tmp=dsqrt(tg(1)**2+tg(2)**2) */
2414: /* tg=tg/tmp */
2415: bn[0] = -nn[2] * tg[1];
2416: bn[1] = nn[2] * tg[0];
2417: bn[2] = nn[0] * tg[1] - nn[1] * tg[0];
2418: /* Computing 2nd power */
2419: d__1 = bn[0];
2420: /* Computing 2nd power */
2421: d__2 = bn[1];
2422: /* Computing 2nd power */
2423: d__3 = bn[2];
2424: tmp = PetscSqrtScalar(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
2425: i__1 = *ndim;
2426: for (k = 1; k <= i__1; ++k) {
2427: bn[k - 1] /= tmp;
2428: }
2429: } else if (*ndim == 2) {
2430: tg[0] = -nn[1];
2431: tg[1] = nn[0];
2432: /* tmp=dsqrt(tg(1)**2+tg(2)**2) */
2433: /* tg=tg/tmp */
2434: bn[0] = 0.;
2435: bn[1] = 0.;
2436: bn[2] = 1.;
2437: }
2438: rl = ul[0];
2439: rr = ur[0];
2440: uxl = 0.;
2441: uxr = 0.;
2442: utl = 0.;
2443: utr = 0.;
2444: ubl = 0.;
2445: ubr = 0.;
2446: i__1 = *ndim;
2447: for (k = 1; k <= i__1; ++k) {
2448: uxl += ul[k] * nn[k-1];
2449: uxr += ur[k] * nn[k-1];
2450: utl += ul[k] * tg[k - 1];
2451: utr += ur[k] * tg[k - 1];
2452: ubl += ul[k] * bn[k - 1];
2453: ubr += ur[k] * bn[k - 1];
2454: }
2455: uxl /= rl;
2456: uxr /= rr;
2457: utl /= rl;
2458: utr /= rr;
2459: ubl /= rl;
2460: ubr /= rr;
2462: gaml = *gamma;
2463: gamr = *gamma;
2464: /* Computing 2nd power */
2465: d__1 = uxl;
2466: /* Computing 2nd power */
2467: d__2 = utl;
2468: /* Computing 2nd power */
2469: d__3 = ubl;
2470: pl = (*gamma - 1.) * (ul[*ndim + 1] - rl * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2471: /* Computing 2nd power */
2472: d__1 = uxr;
2473: /* Computing 2nd power */
2474: d__2 = utr;
2475: /* Computing 2nd power */
2476: d__3 = ubr;
2477: pr = (*gamma - 1.) * (ur[*ndim + 1] - rr * .5 * (d__1 * d__1 + d__2 * d__2 + d__3 * d__3));
2478: rho1l = rl;
2479: rho1r = rr;
2481: iwave = riemannsolver(&xcen, &xp, &dtt, &rl, &uxl, &pl, &utl, &ubl, &gaml, &
2482: rho1l, &rr, &uxr, &pr, &utr, &ubr, &gamr, &rho1r, &rhom, &unm, &
2483: pm, &utm, &ubm, &gamm, &rho1m);
2485: flux[0] = rhom * unm;
2486: fn = rhom * unm * unm + pm;
2487: ft = rhom * unm * utm;
2488: /* flux(2)=fn*nn(1)+ft*nn(2) */
2489: /* flux(3)=fn*tg(1)+ft*tg(2) */
2490: flux[1] = fn * nn[0] + ft * tg[0];
2491: flux[2] = fn * nn[1] + ft * tg[1];
2492: /* flux(2)=rhom*unm*(unm)+pm */
2493: /* flux(3)=rhom*(unm)*utm */
2494: if (*ndim == 3) {
2495: flux[3] = rhom * unm * ubm;
2496: }
2497: flux[*ndim + 1] = (rhom * .5 * (unm * unm + utm * utm + ubm * ubm) + gamm / (gamm - 1.) * pm) * unm;
2498: return iwave;
2499: } /* godunovflux_ */
2501: /* Subroutine to set up the initial conditions for the */
2502: /* Shock Interface interaction or linear wave (Ravi Samtaney,Mark Adams). */
2503: /* ----------------------------------------------------------------------- */
2504: int projecteqstate(PetscReal wc[], const PetscReal ueq[], PetscReal lv[][3])
2505: {
2506: int j,k;
2507: /* Wc=matmul(lv,Ueq) 3 vars */
2508: for (k = 0; k < 3; ++k) {
2509: wc[k] = 0.;
2510: for (j = 0; j < 3; ++j) {
2511: wc[k] += lv[k][j]*ueq[j];
2512: }
2513: }
2514: return 0;
2515: }
2516: /* ----------------------------------------------------------------------- */
2517: int projecttoprim(PetscReal v[], const PetscReal wc[], PetscReal rv[][3])
2518: {
2519: int k,j;
2520: /* V=matmul(rv,WC) 3 vars */
2521: for (k = 0; k < 3; ++k) {
2522: v[k] = 0.;
2523: for (j = 0; j < 3; ++j) {
2524: v[k] += rv[k][j]*wc[j];
2525: }
2526: }
2527: return 0;
2528: }
2529: /* ---------------------------------------------------------------------- */
2530: int eigenvectors(PetscReal rv[][3], PetscReal lv[][3], const PetscReal ueq[], PetscReal gamma)
2531: {
2532: int j,k;
2533: PetscReal rho,csnd,p0;
2534: /* PetscScalar u; */
2536: for (k = 0; k < 3; ++k) for (j = 0; j < 3; ++j) { lv[k][j] = 0.; rv[k][j] = 0.; }
2537: rho = ueq[0];
2538: /* u = ueq[1]; */
2539: p0 = ueq[2];
2540: csnd = PetscSqrtReal(gamma * p0 / rho);
2541: lv[0][1] = rho * .5;
2542: lv[0][2] = -.5 / csnd;
2543: lv[1][0] = csnd;
2544: lv[1][2] = -1. / csnd;
2545: lv[2][1] = rho * .5;
2546: lv[2][2] = .5 / csnd;
2547: rv[0][0] = -1. / csnd;
2548: rv[1][0] = 1. / rho;
2549: rv[2][0] = -csnd;
2550: rv[0][1] = 1. / csnd;
2551: rv[0][2] = 1. / csnd;
2552: rv[1][2] = 1. / rho;
2553: rv[2][2] = csnd;
2554: return 0;
2555: }
2557: int initLinearWave(EulerNode *ux, const PetscReal gamma, const PetscReal coord[], const PetscReal Lx)
2558: {
2559: PetscReal p0,u0,wcp[3],wc[3];
2560: PetscReal lv[3][3];
2561: PetscReal vp[3];
2562: PetscReal rv[3][3];
2563: PetscReal eps, ueq[3], rho0, twopi;
2565: /* Function Body */
2566: twopi = 2.*PETSC_PI;
2567: eps = 1e-4; /* perturbation */
2568: rho0 = 1e3; /* density of water */
2569: p0 = 101325.; /* init pressure of 1 atm (?) */
2570: u0 = 0.;
2571: ueq[0] = rho0;
2572: ueq[1] = u0;
2573: ueq[2] = p0;
2574: /* Project initial state to characteristic variables */
2575: eigenvectors(rv, lv, ueq, gamma);
2576: projecteqstate(wc, ueq, lv);
2577: wcp[0] = wc[0];
2578: wcp[1] = wc[1];
2579: wcp[2] = wc[2] + eps * PetscCosReal(coord[0] * 2. * twopi / Lx);
2580: projecttoprim(vp, wcp, rv);
2581: ux->r = vp[0]; /* density */
2582: ux->ru[0] = vp[0] * vp[1]; /* x momentum */
2583: ux->ru[1] = 0.;
2584: #if defined DIM > 2
2585: if (dim>2) ux->ru[2] = 0.;
2586: #endif
2587: /* E = rho * e + rho * v^2/2 = p/(gam-1) + rho*v^2/2 */
2588: ux->E = vp[2]/(gamma - 1.) + 0.5*vp[0]*vp[1]*vp[1];
2589: return 0;
2590: }
2592: /*TEST
2594: # 2D Advection 0-10
2595: test:
2596: suffix: 0
2597: requires: exodusii
2598: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2600: test:
2601: suffix: 1
2602: requires: exodusii
2603: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2605: test:
2606: suffix: 2
2607: requires: exodusii
2608: nsize: 2
2609: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2611: test:
2612: suffix: 3
2613: requires: exodusii
2614: nsize: 2
2615: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo
2617: test:
2618: suffix: 4
2619: requires: exodusii
2620: nsize: 8
2621: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2623: test:
2624: suffix: 5
2625: requires: exodusii
2626: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw -ts_adapt_reject_safety 1
2628: test:
2629: suffix: 6
2630: requires: exodusii
2631: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/squaremotor-30.exo -ufv_split_faces
2633: test:
2634: suffix: 7
2635: requires: exodusii
2636: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2638: test:
2639: suffix: 8
2640: requires: exodusii
2641: nsize: 2
2642: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2644: test:
2645: suffix: 9
2646: requires: exodusii
2647: nsize: 8
2648: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad-15.exo -dm_refine 1
2650: test:
2651: suffix: 10
2652: requires: exodusii
2653: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside-quad.exo
2655: # 2D Shallow water
2656: test:
2657: suffix: sw_0
2658: requires: exodusii
2659: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -bc_wall 100,101 -physics sw -ufv_cfl 5 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_time 1 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy
2661: test:
2662: suffix: sw_hll
2663: args: -ufv_vtk_interval 0 -bc_wall 1,2,3,4 -physics sw -ufv_cfl 3 -petscfv_type leastsquares -petsclimiter_type sin -ts_max_steps 5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -monitor height,energy -grid_bounds 0,5,0,5 -grid_size 25,25 -sw_riemann hll
2665: # 2D Advection: p4est
2666: test:
2667: suffix: p4est_advec_2d
2668: requires: p4est
2669: args: -ufv_vtk_interval 0 -f -dm_type p4est -dm_forest_minimum_refinement 1 -dm_forest_initial_refinement 2 -dm_p4est_refine_pattern hash -dm_forest_maximum_refinement 5
2671: # Advection in a box
2672: test:
2673: suffix: adv_2d_quad_0
2674: args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2676: test:
2677: suffix: adv_2d_quad_1
2678: args: -ufv_vtk_interval 0 -dm_refine 3 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2679: timeoutfactor: 3
2681: test:
2682: suffix: adv_2d_quad_p4est_0
2683: requires: p4est
2684: args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2686: test:
2687: suffix: adv_2d_quad_p4est_1
2688: requires: p4est
2689: args: -ufv_vtk_interval 0 -dm_refine 5 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2690: timeoutfactor: 3
2692: test:
2693: suffix: adv_2d_quad_p4est_adapt_0
2694: requires: p4est !__float128 #broken for quad precision
2695: args: -ufv_vtk_interval 0 -dm_refine 3 -dm_type p4est -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1 -ufv_use_amr -refine_vec_tagger_box 0.005,inf -coarsen_vec_tagger_box 0,1.e-5 -petscfv_type leastsquares -ts_max_time 0.01
2696: timeoutfactor: 3
2698: test:
2699: suffix: adv_2d_tri_0
2700: requires: triangle
2701: TODO: how did this ever get in main when there is no support for this
2702: args: -ufv_vtk_interval 0 -simplex -dm_refine 3 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2704: test:
2705: suffix: adv_2d_tri_1
2706: requires: triangle
2707: TODO: how did this ever get in main when there is no support for this
2708: args: -ufv_vtk_interval 0 -simplex -dm_refine 5 -dm_plex_separate_marker -grid_bounds -0.5,0.5,-0.5,0.5 -bc_inflow 1,2,4 -bc_outflow 3 -advect_sol_type bump -advect_bump_center 0.25,0 -advect_bump_radius 0.1
2710: test:
2711: suffix: adv_0
2712: requires: exodusii
2713: args: -ufv_vtk_interval 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -bc_inflow 100,101,200 -bc_outflow 201
2715: test:
2716: suffix: shock_0
2717: requires: p4est !single !complex
2718: args: -ufv_vtk_interval 0 -monitor density,energy -f -grid_size 2,1 -grid_bounds -1,1.,0.,1 -bc_wall 1,2,3,4 -dm_type p4est -dm_forest_partition_overlap 1 -dm_forest_maximum_refinement 6 -dm_forest_minimum_refinement 2 -dm_forest_initial_refinement 2 -ufv_use_amr -refine_vec_tagger_box 0.5,inf -coarsen_vec_tagger_box 0,1.e-2 -refine_tag_view -coarsen_tag_view -physics euler -eu_type iv_shock -ufv_cfl 10 -eu_alpha 60. -grid_skew_60 -eu_gamma 1.4 -eu_amach 2.02 -eu_rho2 3. -petscfv_type leastsquares -petsclimiter_type minmod -petscfv_compute_gradients 0 -ts_max_time 0.5 -ts_ssp_type rks2 -ts_ssp_nstages 10 -ufv_vtk_basename ${wPETSC_DIR}/ex11
2719: timeoutfactor: 3
2721: # Test GLVis visualization of PetscFV fields
2722: test:
2723: suffix: glvis_adv_2d_tet
2724: args: -ufv_vtk_interval 0 -ts_monitor_solution glvis: -ts_max_steps 0 -ufv_vtk_monitor 0 -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/square_periodic.msh -dm_plex_gmsh_periodic 0
2726: test:
2727: suffix: glvis_adv_2d_quad
2728: args: -ufv_vtk_interval 0 -ts_monitor_solution glvis: -ts_max_steps 0 -ufv_vtk_monitor 0 -dm_refine 5 -dm_plex_separate_marker -bc_inflow 1,2,4 -bc_outflow 3
2730: test:
2731: suffix: tut_1
2732: requires: exodusii
2733: nsize: 1
2734: args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo
2736: test:
2737: suffix: tut_2
2738: requires: exodusii
2739: nsize: 1
2740: args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/sevenside.exo -ts_type rosw
2742: test:
2743: suffix: tut_3
2744: requires: exodusii
2745: nsize: 4
2746: args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -monitor Error -advect_sol_type bump -petscfv_type leastsquares -petsclimiter_type sin
2748: test:
2749: suffix: tut_4
2750: requires: exodusii
2751: nsize: 4
2752: args: -f ${wPETSC_DIR}/share/petsc/datafiles/meshes/annulus-20.exo -physics sw -monitor Height,Energy -petscfv_type leastsquares -petsclimiter_type minmod
2754: TEST*/