Actual source code: ex9.c
petsc-dev 2014-02-02
1: static const char help[] = "1D periodic Finite Volume solver in slope-limiter form with semidiscrete time stepping.\n"
2: "Solves scalar and vector problems, choose the physical model with -physics\n"
3: " advection - Constant coefficient scalar advection\n"
4: " u_t + (a*u)_x = 0\n"
5: " burgers - Burgers equation\n"
6: " u_t + (u^2/2)_x = 0\n"
7: " traffic - Traffic equation\n"
8: " u_t + (u*(1-u))_x = 0\n"
9: " acoustics - Acoustic wave propagation\n"
10: " u_t + (c*z*v)_x = 0\n"
11: " v_t + (c/z*u)_x = 0\n"
12: " isogas - Isothermal gas dynamics\n"
13: " rho_t + (rho*u)_x = 0\n"
14: " (rho*u)_t + (rho*u^2 + c^2*rho)_x = 0\n"
15: " shallow - Shallow water equations\n"
16: " h_t + (h*u)_x = 0\n"
17: " (h*u)_t + (h*u^2 + g*h^2/2)_x = 0\n"
18: "Some of these physical models have multiple Riemann solvers, select these with -physics_xxx_riemann\n"
19: " exact - Exact Riemann solver which usually needs to perform a Newton iteration to connect\n"
20: " the states across shocks and rarefactions\n"
21: " roe - Linearized scheme, usually with an entropy fix inside sonic rarefactions\n"
22: "The systems provide a choice of reconstructions with -physics_xxx_reconstruct\n"
23: " characteristic - Limit the characteristic variables, this is usually preferred (default)\n"
24: " conservative - Limit the conservative variables directly, can cause undesired interaction of waves\n\n"
25: "A variety of limiters for high-resolution TVD limiters are available with -limit\n"
26: " upwind,minmod,superbee,mc,vanleer,vanalbada,koren,cada-torillhon (last two are nominally third order)\n"
27: " and non-TVD schemes lax-wendroff,beam-warming,fromm\n\n"
28: "To preserve the TVD property, one should time step with a strong stability preserving method.\n"
29: "The optimal high order explicit Runge-Kutta methods in TSSSP are recommended for non-stiff problems.\n\n"
30: "Several initial conditions can be chosen with -initial N\n\n"
31: "The problem size should be set with -da_grid_x M\n\n";
33: #include <petscts.h>
34: #include <petscdmda.h>
35: #include <petscdraw.h>
37: #include <petsc-private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */
39: PETSC_STATIC_INLINE PetscReal Sgn(PetscReal a) { return (a<0) ? -1 : 1; }
40: PETSC_STATIC_INLINE PetscReal Abs(PetscReal a) { return (a<0) ? 0 : a; }
41: PETSC_STATIC_INLINE PetscReal Sqr(PetscReal a) { return a*a; }
42: PETSC_STATIC_INLINE PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; }
43: PETSC_STATIC_INLINE PetscReal MinAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) < PetscAbs(b)) ? a : b; }
44: PETSC_STATIC_INLINE PetscReal MinMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMin(PetscAbs(a),PetscAbs(b)); }
45: PETSC_STATIC_INLINE PetscReal MaxMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMax(PetscAbs(a),PetscAbs(b)); }
46: PETSC_STATIC_INLINE PetscReal MinMod3(PetscReal a,PetscReal b,PetscReal c) {return (a*b<0 || a*c<0) ? 0 : Sgn(a)*PetscMin(PetscAbs(a),PetscMin(PetscAbs(b),PetscAbs(c))); }
48: PETSC_STATIC_INLINE PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin + fmod(range+fmod(a,range),range); }
51: /* ----------------------- Lots of limiters, these could go in a separate library ------------------------- */
52: typedef struct _LimitInfo {
53: PetscReal hx;
54: PetscInt m;
55: } *LimitInfo;
56: static void Limit_Upwind(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
57: {
58: PetscInt i;
59: for (i=0; i<info->m; i++) lmt[i] = 0;
60: }
61: static void Limit_LaxWendroff(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
62: {
63: PetscInt i;
64: for (i=0; i<info->m; i++) lmt[i] = jR[i];
65: }
66: static void Limit_BeamWarming(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
67: {
68: PetscInt i;
69: for (i=0; i<info->m; i++) lmt[i] = jL[i];
70: }
71: static void Limit_Fromm(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
72: {
73: PetscInt i;
74: for (i=0; i<info->m; i++) lmt[i] = 0.5*(jL[i] + jR[i]);
75: }
76: static void Limit_Minmod(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
77: {
78: PetscInt i;
79: for (i=0; i<info->m; i++) lmt[i] = MinMod2(jL[i],jR[i]);
80: }
81: static void Limit_Superbee(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
82: {
83: PetscInt i;
84: for (i=0; i<info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i],2*jR[i]),MinMod2(2*jL[i],jR[i]));
85: }
86: static void Limit_MC(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
87: {
88: PetscInt i;
89: for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],0.5*(jL[i]+jR[i]),2*jR[i]);
90: }
91: static void Limit_VanLeer(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
92: { /* phi = (t + abs(t)) / (1 + abs(t)) */
93: PetscInt i;
94: for (i=0; i<info->m; i++) lmt[i] = (jL[i]*Abs(jR[i]) + Abs(jL[i])*jR[i]) / (Abs(jL[i]) + Abs(jR[i]) + 1e-15);
95: }
96: static void Limit_VanAlbada(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
97: { /* phi = (t + t^2) / (1 + t^2) */
98: PetscInt i;
99: for (i=0; i<info->m; i++) lmt[i] = (jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15);
100: }
101: static void Limit_VanAlbadaTVD(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
102: { /* phi = (t + t^2) / (1 + t^2) */
103: PetscInt i;
104: for (i=0; i<info->m; i++) lmt[i] = (jL[i]*jR[i]<0) ? 0
105: : (jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15);
106: }
107: static void Limit_Koren(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
108: { /* phi = (t + 2*t^2) / (2 - t + 2*t^2) */
109: PetscInt i;
110: for (i=0; i<info->m; i++) lmt[i] = ((jL[i]*Sqr(jR[i]) + 2*Sqr(jL[i])*jR[i])
111: / (2*Sqr(jL[i]) - jL[i]*jR[i] + 2*Sqr(jR[i]) + 1e-15));
112: }
113: static void Limit_KorenSym(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
114: { /* Symmetric version of above */
115: PetscInt i;
116: for (i=0; i<info->m; i++) lmt[i] = (1.5*(jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i])
117: / (2*Sqr(jL[i]) - jL[i]*jR[i] + 2*Sqr(jR[i]) + 1e-15));
118: }
119: static void Limit_Koren3(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
120: { /* Eq 11 of Cada-Torrilhon 2009 */
121: PetscInt i;
122: for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i]);
123: }
125: static PetscReal CadaTorrilhonPhiHatR_Eq13(PetscReal L,PetscReal R)
126: {
127: return PetscMax(0,PetscMin((L+2*R)/3,
128: PetscMax(-0.5*L,
129: PetscMin(2*L,
130: PetscMin((L+2*R)/3,1.6*R)))));
131: }
132: static void Limit_CadaTorrilhon2(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
133: { /* Cada-Torrilhon 2009, Eq 13 */
134: PetscInt i;
135: for (i=0; i<info->m; i++) lmt[i] = CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]);
136: }
137: static void Limit_CadaTorrilhon3R(PetscReal r,LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
138: { /* Cada-Torrilhon 2009, Eq 22 */
139: /* They recommend 0.001 < r < 1, but larger values are more accurate in smooth regions */
140: const PetscReal eps = 1e-7,hx = info->hx;
141: PetscInt i;
142: for (i=0; i<info->m; i++) {
143: const PetscReal eta = (Sqr(jL[i]) + Sqr(jR[i])) / Sqr(r*hx);
144: lmt[i] = ((eta < 1-eps)
145: ? (jL[i] + 2*jR[i]) / 3
146: : ((eta > 1+eps)
147: ? CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i])
148: : 0.5*((1-(eta-1)/eps)*(jL[i]+2*jR[i])/3
149: + (1+(eta+1)/eps)*CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]))));
150: }
151: }
152: static void Limit_CadaTorrilhon3R0p1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
153: { Limit_CadaTorrilhon3R(0.1,info,jL,jR,lmt); }
154: static void Limit_CadaTorrilhon3R1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
155: { Limit_CadaTorrilhon3R(1,info,jL,jR,lmt); }
156: static void Limit_CadaTorrilhon3R10(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
157: { Limit_CadaTorrilhon3R(10,info,jL,jR,lmt); }
158: static void Limit_CadaTorrilhon3R100(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
159: { Limit_CadaTorrilhon3R(100,info,jL,jR,lmt); }
162: /* --------------------------------- Finite Volume data structures ----------------------------------- */
164: typedef enum {FVBC_PERIODIC, FVBC_OUTFLOW} FVBCType;
165: static const char *FVBCTypes[] = {"PERIODIC","OUTFLOW","FVBCType","FVBC_",0};
166: typedef PetscErrorCode (*RiemannFunction)(void*,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscReal*);
167: typedef PetscErrorCode (*ReconstructFunction)(void*,PetscInt,const PetscScalar*,PetscScalar*,PetscScalar*,PetscReal*);
169: typedef struct {
170: PetscErrorCode (*sample)(void*,PetscInt,FVBCType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*);
171: RiemannFunction riemann;
172: ReconstructFunction characteristic;
173: PetscErrorCode (*destroy)(void*);
174: void *user;
175: PetscInt dof;
176: char *fieldname[16];
177: } PhysicsCtx;
179: typedef struct {
180: void (*limit)(LimitInfo,const PetscScalar*,const PetscScalar*,PetscScalar*);
181: PhysicsCtx physics;
183: MPI_Comm comm;
184: char prefix[256];
186: /* Local work arrays */
187: PetscScalar *R,*Rinv; /* Characteristic basis, and it's inverse. COLUMN-MAJOR */
188: PetscScalar *cjmpLR; /* Jumps at left and right edge of cell, in characteristic basis, len=2*dof */
189: PetscScalar *cslope; /* Limited slope, written in characteristic basis */
190: PetscScalar *uLR; /* Solution at left and right of interface, conservative variables, len=2*dof */
191: PetscScalar *flux; /* Flux across interface */
192: PetscReal *speeds; /* Speeds of each wave */
194: PetscReal cfl_idt; /* Max allowable value of 1/Delta t */
195: PetscReal cfl;
196: PetscReal xmin,xmax;
197: PetscInt initial;
198: PetscBool exact;
199: FVBCType bctype;
200: } FVCtx;
203: /* Utility */
207: PetscErrorCode RiemannListAdd(PetscFunctionList *flist,const char *name,RiemannFunction rsolve)
208: {
212: PetscFunctionListAdd(flist,name,rsolve);
213: return(0);
214: }
218: PetscErrorCode RiemannListFind(PetscFunctionList flist,const char *name,RiemannFunction *rsolve)
219: {
223: PetscFunctionListFind(flist,name,rsolve);
224: if (!*rsolve) SETERRQ1(PETSC_COMM_SELF,1,"Riemann solver \"%s\" could not be found",name);
225: return(0);
226: }
230: PetscErrorCode ReconstructListAdd(PetscFunctionList *flist,const char *name,ReconstructFunction r)
231: {
235: PetscFunctionListAdd(flist,name,r);
236: return(0);
237: }
241: PetscErrorCode ReconstructListFind(PetscFunctionList flist,const char *name,ReconstructFunction *r)
242: {
246: PetscFunctionListFind(flist,name,r);
247: if (!*r) SETERRQ1(PETSC_COMM_SELF,1,"Reconstruction \"%s\" could not be found",name);
248: return(0);
249: }
252: /* --------------------------------- Physics ----------------------------------- */
253: /**
254: * Each physical model consists of Riemann solver and a function to determine the basis to use for reconstruction. These
255: * are set with the PhysicsCreate_XXX function which allocates private storage and sets these methods as well as the
256: * number of fields and their names, and a function to deallocate private storage.
257: **/
259: /* First a few functions useful to several different physics */
262: static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
263: {
264: PetscInt i,j;
267: for (i=0; i<m; i++) {
268: for (j=0; j<m; j++) Xi[i*m+j] = X[i*m+j] = (PetscScalar)(i==j);
269: speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */
270: }
271: return(0);
272: }
276: static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx)
277: {
281: PetscFree(vctx);
282: return(0);
283: }
287: /* --------------------------------- Advection ----------------------------------- */
289: typedef struct {
290: PetscReal a; /* advective velocity */
291: } AdvectCtx;
295: static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
296: {
297: AdvectCtx *ctx = (AdvectCtx*)vctx;
298: PetscReal speed;
301: speed = ctx->a;
302: flux[0] = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0];
303: *maxspeed = speed;
304: return(0);
305: }
309: static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
310: {
311: AdvectCtx *ctx = (AdvectCtx*)vctx;
314: X[0] = 1.;
315: Xi[0] = 1.;
316: speeds[0] = ctx->a;
317: return(0);
318: }
322: static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
323: {
324: AdvectCtx *ctx = (AdvectCtx*)vctx;
325: PetscReal a = ctx->a,x0;
328: switch (bctype) {
329: case FVBC_OUTFLOW: x0 = x-a*t; break;
330: case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break;
331: default: SETERRQ(PETSC_COMM_SELF,1,"unknown BCType");
332: }
333: switch (initial) {
334: case 0: u[0] = (x0 < 0) ? 1 : -1; break;
335: case 1: u[0] = (x0 < 0) ? -1 : 1; break;
336: case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
337: case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
338: case 4: u[0] = PetscAbs(x0); break;
339: case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
340: case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
341: default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
342: }
343: return(0);
344: }
348: static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
349: {
351: AdvectCtx *user;
354: PetscNew(&user);
355: ctx->physics.sample = PhysicsSample_Advect;
356: ctx->physics.riemann = PhysicsRiemann_Advect;
357: ctx->physics.characteristic = PhysicsCharacteristic_Advect;
358: ctx->physics.destroy = PhysicsDestroy_SimpleFree;
359: ctx->physics.user = user;
360: ctx->physics.dof = 1;
361: PetscStrallocpy("u",&ctx->physics.fieldname[0]);
362: user->a = 1;
363: PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");
364: {
365: PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL);
366: }
367: PetscOptionsEnd();
368: return(0);
369: }
373: /* --------------------------------- Burgers ----------------------------------- */
375: typedef struct {
376: PetscReal lxf_speed;
377: } BurgersCtx;
381: static PetscErrorCode PhysicsSample_Burgers(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
382: {
385: if (bctype == FVBC_PERIODIC && t > 0) SETERRQ(PETSC_COMM_SELF,1,"Exact solution not implemented for periodic");
386: switch (initial) {
387: case 0: u[0] = (x < 0) ? 1 : -1; break;
388: case 1:
389: if (x < -t) u[0] = -1;
390: else if (x < t) u[0] = x/t;
391: else u[0] = 1;
392: break;
393: case 2:
394: if (x < 0) u[0] = 0;
395: else if (x <= t) u[0] = x/t;
396: else if (x < 1+0.5*t) u[0] = 1;
397: else u[0] = 0;
398: break;
399: case 3:
400: if (x < 0.2*t) u[0] = 0.2;
401: else if (x < t) u[0] = x/t;
402: else u[0] = 1;
403: break;
404: case 4:
405: if (t > 0) SETERRQ(PETSC_COMM_SELF,1,"Only initial condition available");
406: u[0] = 0.7 + 0.3*PetscSinReal(2*PETSC_PI*((x-xmin)/(xmax-xmin)));
407: break;
408: case 5: /* Pure shock solution */
409: if (x < 0.5*t) u[0] = 1;
410: else u[0] = 0;
411: break;
412: default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
413: }
414: return(0);
415: }
419: static PetscErrorCode PhysicsRiemann_Burgers_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
420: {
423: if (uL[0] < uR[0]) { /* rarefaction */
424: flux[0] = (uL[0]*uR[0] < 0)
425: ? 0 /* sonic rarefaction */
426: : 0.5*PetscMin(PetscSqr(uL[0]),PetscSqr(uR[0]));
427: } else { /* shock */
428: flux[0] = 0.5*PetscMax(PetscSqr(uL[0]),PetscSqr(uR[0]));
429: }
430: *maxspeed = (PetscAbs(uL[0]) > PetscAbs(uR[0])) ? uL[0] : uR[0];
431: return(0);
432: }
436: static PetscErrorCode PhysicsRiemann_Burgers_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
437: {
438: PetscReal speed;
441: speed = 0.5*(uL[0] + uR[0]);
442: flux[0] = 0.25*(PetscSqr(uL[0]) + PetscSqr(uR[0])) - 0.5*PetscAbs(speed)*(uR[0]-uL[0]);
443: if (uL[0] <= 0 && 0 <= uR[0]) flux[0] = 0; /* Entropy fix for sonic rarefaction */
444: *maxspeed = speed;
445: return(0);
446: }
450: static PetscErrorCode PhysicsRiemann_Burgers_LxF(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
451: {
452: PetscReal c;
453: PetscScalar fL,fR;
456: c = ((BurgersCtx*)vctx)->lxf_speed;
457: fL = 0.5*PetscSqr(uL[0]);
458: fR = 0.5*PetscSqr(uR[0]);
459: flux[0] = 0.5*(fL + fR) - 0.5*c*(uR[0] - uL[0]);
460: *maxspeed = c;
461: return(0);
462: }
466: static PetscErrorCode PhysicsRiemann_Burgers_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
467: {
468: PetscReal c;
469: PetscScalar fL,fR;
472: c = PetscMax(PetscAbs(uL[0]),PetscAbs(uR[0]));
473: fL = 0.5*PetscSqr(uL[0]);
474: fR = 0.5*PetscSqr(uR[0]);
475: flux[0] = 0.5*(fL + fR) - 0.5*c*(uR[0] - uL[0]);
476: *maxspeed = c;
477: return(0);
478: }
482: static PetscErrorCode PhysicsCreate_Burgers(FVCtx *ctx)
483: {
484: BurgersCtx *user;
485: PetscErrorCode ierr;
486: RiemannFunction r;
487: PetscFunctionList rlist = 0;
488: char rname[256] = "exact";
491: PetscNew(&user);
493: ctx->physics.sample = PhysicsSample_Burgers;
494: ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
495: ctx->physics.destroy = PhysicsDestroy_SimpleFree;
496: ctx->physics.user = user;
497: ctx->physics.dof = 1;
499: PetscStrallocpy("u",&ctx->physics.fieldname[0]);
500: RiemannListAdd(&rlist,"exact", PhysicsRiemann_Burgers_Exact);
501: RiemannListAdd(&rlist,"roe", PhysicsRiemann_Burgers_Roe);
502: RiemannListAdd(&rlist,"lxf", PhysicsRiemann_Burgers_LxF);
503: RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Burgers_Rusanov);
504: PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");
505: {
506: PetscOptionsFList("-physics_burgers_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
507: }
508: PetscOptionsEnd();
509: RiemannListFind(rlist,rname,&r);
510: PetscFunctionListDestroy(&rlist);
511: ctx->physics.riemann = r;
513: /* *
514: * Hack to deal with LxF in semi-discrete form
515: * max speed is 1 for the basic initial conditions (where |u| <= 1)
516: * */
517: if (r == PhysicsRiemann_Burgers_LxF) user->lxf_speed = 1;
518: return(0);
519: }
523: /* --------------------------------- Traffic ----------------------------------- */
525: typedef struct {
526: PetscReal lxf_speed;
527: PetscReal a;
528: } TrafficCtx;
530: PETSC_STATIC_INLINE PetscScalar TrafficFlux(PetscScalar a,PetscScalar u) { return a*u*(1-u); }
534: static PetscErrorCode PhysicsSample_Traffic(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
535: {
536: PetscReal a = ((TrafficCtx*)vctx)->a;
539: if (bctype == FVBC_PERIODIC && t > 0) SETERRQ(PETSC_COMM_SELF,1,"Exact solution not implemented for periodic");
540: switch (initial) {
541: case 0:
542: u[0] = (-a*t < x) ? 2 : 0; break;
543: case 1:
544: if (x < PetscMin(2*a*t,0.5+a*t)) u[0] = -1;
545: else if (x < 1) u[0] = 0;
546: else u[0] = 1;
547: break;
548: case 2:
549: if (t > 0) SETERRQ(PETSC_COMM_SELF,1,"Only initial condition available");
550: u[0] = 0.7 + 0.3*PetscSinReal(2*PETSC_PI*((x-xmin)/(xmax-xmin)));
551: break;
552: default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
553: }
554: return(0);
555: }
559: static PetscErrorCode PhysicsRiemann_Traffic_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
560: {
561: PetscReal a = ((TrafficCtx*)vctx)->a;
564: if (uL[0] < uR[0]) {
565: flux[0] = PetscMin(TrafficFlux(a,uL[0]),
566: TrafficFlux(a,uR[0]));
567: } else {
568: flux[0] = (uR[0] < 0.5 && 0.5 < uL[0])
569: ? TrafficFlux(a,0.5)
570: : PetscMax(TrafficFlux(a,uL[0]),
571: TrafficFlux(a,uR[0]));
572: }
573: *maxspeed = a*MaxAbs(1-2*uL[0],1-2*uR[0]);
574: return(0);
575: }
579: static PetscErrorCode PhysicsRiemann_Traffic_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
580: {
581: PetscReal a = ((TrafficCtx*)vctx)->a;
582: PetscReal speed;
585: speed = a*(1 - (uL[0] + uR[0]));
586: flux[0] = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*PetscAbs(speed)*(uR[0]-uL[0]);
587: *maxspeed = speed;
588: return(0);
589: }
593: static PetscErrorCode PhysicsRiemann_Traffic_LxF(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
594: {
595: TrafficCtx *phys = (TrafficCtx*)vctx;
596: PetscReal a = phys->a;
597: PetscReal speed;
600: speed = a*(1 - (uL[0] + uR[0]));
601: flux[0] = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*phys->lxf_speed*(uR[0]-uL[0]);
602: *maxspeed = speed;
603: return(0);
604: }
608: static PetscErrorCode PhysicsRiemann_Traffic_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
609: {
610: PetscReal a = ((TrafficCtx*)vctx)->a;
611: PetscReal speed;
614: speed = a*PetscMax(PetscAbs(1-2*uL[0]),PetscAbs(1-2*uR[0]));
615: flux[0] = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*speed*(uR[0]-uL[0]);
616: *maxspeed = speed;
617: return(0);
618: }
622: static PetscErrorCode PhysicsCreate_Traffic(FVCtx *ctx)
623: {
624: PetscErrorCode ierr;
625: TrafficCtx *user;
626: RiemannFunction r;
627: PetscFunctionList rlist = 0;
628: char rname[256] = "exact";
631: PetscNew(&user);
632: ctx->physics.sample = PhysicsSample_Traffic;
633: ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
634: ctx->physics.destroy = PhysicsDestroy_SimpleFree;
635: ctx->physics.user = user;
636: ctx->physics.dof = 1;
638: PetscStrallocpy("density",&ctx->physics.fieldname[0]);
639: user->a = 0.5;
640: RiemannListAdd(&rlist,"exact", PhysicsRiemann_Traffic_Exact);
641: RiemannListAdd(&rlist,"roe", PhysicsRiemann_Traffic_Roe);
642: RiemannListAdd(&rlist,"lxf", PhysicsRiemann_Traffic_LxF);
643: RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Traffic_Rusanov);
644: PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Traffic","");
645: {
646: PetscOptionsReal("-physics_traffic_a","Flux = a*u*(1-u)","",user->a,&user->a,NULL);
647: PetscOptionsFList("-physics_traffic_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
648: }
649: PetscOptionsEnd();
651: RiemannListFind(rlist,rname,&r);
652: PetscFunctionListDestroy(&rlist);
654: ctx->physics.riemann = r;
656: /* *
657: * Hack to deal with LxF in semi-discrete form
658: * max speed is 3*a for the basic initial conditions (-1 <= u <= 2)
659: * */
660: if (r == PhysicsRiemann_Traffic_LxF) user->lxf_speed = 3*user->a;
661: return(0);
662: }
664: /* --------------------------------- Linear Acoustics ----------------------------------- */
666: /* Flux: u_t + (A u)_x
667: * z = sqrt(rho*bulk), c = sqrt(rho/bulk)
668: * Spectral decomposition: A = R * D * Rinv
669: * [ cz] = [-z z] [-c ] [-1/2z 1/2]
670: * [c/z ] = [ 1 1] [ c] [ 1/2z 1/2]
671: *
672: * We decompose this into the left-traveling waves Al = R * D^- Rinv
673: * and the right-traveling waves Ar = R * D^+ * Rinv
674: * Multiplying out these expressions produces the following two matrices
675: */
677: typedef struct {
678: PetscReal c; /* speed of sound: c = sqrt(bulk/rho) */
679: PetscReal z; /* impedence: z = sqrt(rho*bulk) */
680: } AcousticsCtx;
682: PETSC_STATIC_INLINE void AcousticsFlux(AcousticsCtx *ctx,const PetscScalar *u,PetscScalar *f)
683: {
684: f[0] = ctx->c*ctx->z*u[1];
685: f[1] = ctx->c/ctx->z*u[0];
686: }
690: static PetscErrorCode PhysicsCharacteristic_Acoustics(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
691: {
692: AcousticsCtx *phys = (AcousticsCtx*)vctx;
693: PetscReal z = phys->z,c = phys->c;
696: X[0*2+0] = -z;
697: X[0*2+1] = z;
698: X[1*2+0] = 1;
699: X[1*2+1] = 1;
700: Xi[0*2+0] = -1./(2*z);
701: Xi[0*2+1] = 1./2;
702: Xi[1*2+0] = 1./(2*z);
703: Xi[1*2+1] = 1./2;
704: speeds[0] = -c;
705: speeds[1] = c;
706: return(0);
707: }
711: static PetscErrorCode PhysicsSample_Acoustics_Initial(AcousticsCtx *phys,PetscInt initial,PetscReal xmin,PetscReal xmax,PetscReal x,PetscReal *u)
712: {
714: switch (initial) {
715: case 0:
716: u[0] = (PetscAbs((x - xmin)/(xmax - xmin) - 0.2) < 0.1) ? 1 : 0.5;
717: u[1] = (PetscAbs((x - xmin)/(xmax - xmin) - 0.7) < 0.1) ? 1 : -0.5;
718: break;
719: case 1:
720: u[0] = PetscCosReal(3 * 2*PETSC_PI*x/(xmax-xmin));
721: u[1] = PetscExpReal(-PetscSqr(x - (xmax + xmin)/2) / (2*PetscSqr(0.2*(xmax - xmin)))) - 0.5;
722: break;
723: default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
724: }
725: return(0);
726: }
730: static PetscErrorCode PhysicsSample_Acoustics(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
731: {
732: AcousticsCtx *phys = (AcousticsCtx*)vctx;
733: PetscReal c = phys->c;
734: PetscReal x0a,x0b,u0a[2],u0b[2],tmp[2];
735: PetscReal X[2][2],Xi[2][2],dummy[2];
739: switch (bctype) {
740: case FVBC_OUTFLOW:
741: x0a = x+c*t;
742: x0b = x-c*t;
743: break;
744: case FVBC_PERIODIC:
745: x0a = RangeMod(x+c*t,xmin,xmax);
746: x0b = RangeMod(x-c*t,xmin,xmax);
747: break;
748: default: SETERRQ(PETSC_COMM_SELF,1,"unknown BCType");
749: }
750: PhysicsSample_Acoustics_Initial(phys,initial,xmin,xmax,x0a,u0a);
751: PhysicsSample_Acoustics_Initial(phys,initial,xmin,xmax,x0b,u0b);
752: PhysicsCharacteristic_Acoustics(vctx,2,u,&X[0][0],&Xi[0][0],dummy);
753: tmp[0] = Xi[0][0]*u0a[0] + Xi[0][1]*u0a[1];
754: tmp[1] = Xi[1][0]*u0b[0] + Xi[1][1]*u0b[1];
755: u[0] = X[0][0]*tmp[0] + X[0][1]*tmp[1];
756: u[1] = X[1][0]*tmp[0] + X[1][1]*tmp[1];
757: return(0);
758: }
762: static PetscErrorCode PhysicsRiemann_Acoustics_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
763: {
764: AcousticsCtx *phys = (AcousticsCtx*)vctx;
765: PetscReal c = phys->c,z = phys->z;
766: PetscReal
767: Al[2][2] = {{-c/2 , c*z/2 },
768: {c/(2*z) , -c/2 }}, /* Left traveling waves */
769: Ar[2][2] = {{c/2 , c*z/2 },
770: {c/(2*z) , c/2 }}; /* Right traveling waves */
773: flux[0] = Al[0][0]*uR[0] + Al[0][1]*uR[1] + Ar[0][0]*uL[0] + Ar[0][1]*uL[1];
774: flux[1] = Al[1][0]*uR[0] + Al[1][1]*uR[1] + Ar[1][0]*uL[0] + Ar[1][1]*uL[1];
775: *maxspeed = c;
776: return(0);
777: }
781: static PetscErrorCode PhysicsCreate_Acoustics(FVCtx *ctx)
782: {
783: PetscErrorCode ierr;
784: AcousticsCtx *user;
785: PetscFunctionList rlist = 0,rclist = 0;
786: char rname[256] = "exact",rcname[256] = "characteristic";
789: PetscNew(&user);
790: ctx->physics.sample = PhysicsSample_Acoustics;
791: ctx->physics.destroy = PhysicsDestroy_SimpleFree;
792: ctx->physics.user = user;
793: ctx->physics.dof = 2;
795: PetscStrallocpy("u",&ctx->physics.fieldname[0]);
796: PetscStrallocpy("v",&ctx->physics.fieldname[1]);
798: user->c = 1;
799: user->z = 1;
801: RiemannListAdd(&rlist,"exact", PhysicsRiemann_Acoustics_Exact);
802: ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_Acoustics);
803: ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative);
804: PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for linear Acoustics","");
805: {
806: PetscOptionsReal("-physics_acoustics_c","c = sqrt(bulk/rho)","",user->c,&user->c,NULL);
807: PetscOptionsReal("-physics_acoustics_z","z = sqrt(bulk*rho)","",user->z,&user->z,NULL);
808: PetscOptionsFList("-physics_acoustics_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
809: PetscOptionsFList("-physics_acoustics_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);
810: }
811: PetscOptionsEnd();
812: RiemannListFind(rlist,rname,&ctx->physics.riemann);
813: ReconstructListFind(rclist,rcname,&ctx->physics.characteristic);
814: PetscFunctionListDestroy(&rlist);
815: PetscFunctionListDestroy(&rclist);
816: return(0);
817: }
819: /* --------------------------------- Isothermal Gas Dynamics ----------------------------------- */
821: typedef struct {
822: PetscReal acoustic_speed;
823: } IsoGasCtx;
825: PETSC_STATIC_INLINE void IsoGasFlux(PetscReal c,const PetscScalar *u,PetscScalar *f)
826: {
827: f[0] = u[1];
828: f[1] = PetscSqr(u[1])/u[0] + c*c*u[0];
829: }
833: static PetscErrorCode PhysicsSample_IsoGas(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
834: {
837: if (t > 0) SETERRQ(PETSC_COMM_SELF,1,"Exact solutions not implemented for t > 0");
838: switch (initial) {
839: case 0:
840: u[0] = (x < 0) ? 1 : 0.5;
841: u[1] = (x < 0) ? 1 : 0.7;
842: break;
843: case 1:
844: u[0] = 1+0.5*PetscSinReal(2*PETSC_PI*x);
845: u[1] = 1*u[0];
846: break;
847: default: SETERRQ(PETSC_COMM_SELF,1,"unknown initial condition");
848: }
849: return(0);
850: }
854: static PetscErrorCode PhysicsRiemann_IsoGas_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
855: {
856: IsoGasCtx *phys = (IsoGasCtx*)vctx;
857: PetscReal c = phys->acoustic_speed;
858: PetscScalar ubar,du[2],a[2],fL[2],fR[2],lam[2],ustar[2],R[2][2];
859: PetscInt i;
862: ubar = (uL[1]/PetscSqrtScalar(uL[0]) + uR[1]/PetscSqrtScalar(uR[0])) / (PetscSqrtScalar(uL[0]) + PetscSqrtScalar(uR[0]));
863: /* write fluxuations in characteristic basis */
864: du[0] = uR[0] - uL[0];
865: du[1] = uR[1] - uL[1];
866: a[0] = (1/(2*c)) * ((ubar + c)*du[0] - du[1]);
867: a[1] = (1/(2*c)) * ((-ubar + c)*du[0] + du[1]);
868: /* wave speeds */
869: lam[0] = ubar - c;
870: lam[1] = ubar + c;
871: /* Right eigenvectors */
872: R[0][0] = 1; R[0][1] = ubar-c;
873: R[1][0] = 1; R[1][1] = ubar+c;
874: /* Compute state in star region (between the 1-wave and 2-wave) */
875: for (i=0; i<2; i++) ustar[i] = uL[i] + a[0]*R[0][i];
876: if (uL[1]/uL[0] < c && c < ustar[1]/ustar[0]) { /* 1-wave is sonic rarefaction */
877: PetscScalar ufan[2];
878: ufan[0] = uL[0]*PetscExpScalar(uL[1]/(uL[0]*c) - 1);
879: ufan[1] = c*ufan[0];
880: IsoGasFlux(c,ufan,flux);
881: } else if (ustar[1]/ustar[0] < -c && -c < uR[1]/uR[0]) { /* 2-wave is sonic rarefaction */
882: PetscScalar ufan[2];
883: ufan[0] = uR[0]*PetscExpScalar(-uR[1]/(uR[0]*c) - 1);
884: ufan[1] = -c*ufan[0];
885: IsoGasFlux(c,ufan,flux);
886: } else { /* Centered form */
887: IsoGasFlux(c,uL,fL);
888: IsoGasFlux(c,uR,fR);
889: for (i=0; i<2; i++) {
890: PetscScalar absdu = PetscAbsScalar(lam[0])*a[0]*R[0][i] + PetscAbsScalar(lam[1])*a[1]*R[1][i];
891: flux[i] = 0.5*(fL[i]+fR[i]) - 0.5*absdu;
892: }
893: }
894: *maxspeed = MaxAbs(lam[0],lam[1]);
895: return(0);
896: }
900: static PetscErrorCode PhysicsRiemann_IsoGas_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
901: {
902: IsoGasCtx *phys = (IsoGasCtx*)vctx;
903: PetscReal c = phys->acoustic_speed;
904: PetscScalar ustar[2];
905: struct {PetscScalar rho,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
906: PetscInt i;
907: PetscErrorCode ierr;
910: if (!(L.rho > 0 && R.rho > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed density is negative");
911: {
912: /* Solve for star state */
913: PetscScalar res,tmp,rho = 0.5*(L.rho + R.rho); /* initial guess */
914: for (i=0; i<20; i++) {
915: PetscScalar fr,fl,dfr,dfl;
916: fl = (L.rho < rho)
917: ? (rho-L.rho)/PetscSqrtScalar(L.rho*rho) /* shock */
918: : PetscLogScalar(rho) - PetscLogScalar(L.rho); /* rarefaction */
919: fr = (R.rho < rho)
920: ? (rho-R.rho)/PetscSqrtScalar(R.rho*rho) /* shock */
921: : PetscLogScalar(rho) - PetscLogScalar(R.rho); /* rarefaction */
922: res = R.u-L.u + c*(fr+fl);
923: PetscIsInfOrNanScalar(res);
924: if (PetscAbsScalar(res) < 1e-10) {
925: star.rho = rho;
926: star.u = L.u - c*fl;
927: goto converged;
928: }
929: dfl = (L.rho < rho)
930: ? 1/PetscSqrtScalar(L.rho*rho)*(1 - 0.5*(rho-L.rho)/rho)
931: : 1/rho;
932: dfr = (R.rho < rho)
933: ? 1/PetscSqrtScalar(R.rho*rho)*(1 - 0.5*(rho-R.rho)/rho)
934: : 1/rho;
935: tmp = rho - res/(c*(dfr+dfl));
936: if (tmp <= 0) rho /= 2; /* Guard against Newton shooting off to a negative density */
937: else rho = tmp;
938: if (!((rho > 0) && PetscIsNormalScalar(rho))) SETERRQ1(PETSC_COMM_SELF,1,"non-normal iterate rho=%g",(double)PetscRealPart(rho));
939: }
940: SETERRQ1(PETSC_COMM_SELF,1,"Newton iteration for star.rho diverged after %D iterations",i);
941: }
942: converged:
943: if (L.u-c < 0 && 0 < star.u-c) { /* 1-wave is sonic rarefaction */
944: PetscScalar ufan[2];
945: ufan[0] = L.rho*PetscExpScalar(L.u/c - 1);
946: ufan[1] = c*ufan[0];
947: IsoGasFlux(c,ufan,flux);
948: } else if (star.u+c < 0 && 0 < R.u+c) { /* 2-wave is sonic rarefaction */
949: PetscScalar ufan[2];
950: ufan[0] = R.rho*PetscExpScalar(-R.u/c - 1);
951: ufan[1] = -c*ufan[0];
952: IsoGasFlux(c,ufan,flux);
953: } else if ((L.rho >= star.rho && L.u-c >= 0)
954: || (L.rho < star.rho && (star.rho*star.u-L.rho*L.u)/(star.rho-L.rho) > 0)) {
955: /* 1-wave is supersonic rarefaction, or supersonic shock */
956: IsoGasFlux(c,uL,flux);
957: } else if ((star.rho <= R.rho && R.u+c <= 0)
958: || (star.rho > R.rho && (R.rho*R.u-star.rho*star.u)/(R.rho-star.rho) < 0)) {
959: /* 2-wave is supersonic rarefaction or supersonic shock */
960: IsoGasFlux(c,uR,flux);
961: } else {
962: ustar[0] = star.rho;
963: ustar[1] = star.rho*star.u;
964: IsoGasFlux(c,ustar,flux);
965: }
966: *maxspeed = MaxAbs(MaxAbs(star.u-c,star.u+c),MaxAbs(L.u-c,R.u+c));
967: return(0);
968: }
972: static PetscErrorCode PhysicsRiemann_IsoGas_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
973: {
974: IsoGasCtx *phys = (IsoGasCtx*)vctx;
975: PetscScalar c = phys->acoustic_speed,fL[2],fR[2],s;
976: struct {PetscScalar rho,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};
979: if (!(L.rho > 0 && R.rho > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed density is negative");
980: IsoGasFlux(c,uL,fL);
981: IsoGasFlux(c,uR,fR);
982: s = PetscMax(PetscAbs(L.u),PetscAbs(R.u))+c;
983: flux[0] = 0.5*(fL[0] + fR[0]) + 0.5*s*(uL[0] - uR[0]);
984: flux[1] = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
985: *maxspeed = s;
986: return(0);
987: }
991: static PetscErrorCode PhysicsCharacteristic_IsoGas(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
992: {
993: IsoGasCtx *phys = (IsoGasCtx*)vctx;
994: PetscReal c = phys->acoustic_speed;
998: speeds[0] = u[1]/u[0] - c;
999: speeds[1] = u[1]/u[0] + c;
1000: X[0*2+0] = 1;
1001: X[0*2+1] = speeds[0];
1002: X[1*2+0] = 1;
1003: X[1*2+1] = speeds[1];
1005: PetscMemcpy(Xi,X,4*sizeof(X[0]));
1006: PetscKernel_A_gets_inverse_A_2(Xi,0);
1007: return(0);
1008: }
1012: static PetscErrorCode PhysicsCreate_IsoGas(FVCtx *ctx)
1013: {
1014: PetscErrorCode ierr;
1015: IsoGasCtx *user;
1016: PetscFunctionList rlist = 0,rclist = 0;
1017: char rname[256] = "exact",rcname[256] = "characteristic";
1020: PetscNew(&user);
1021: ctx->physics.sample = PhysicsSample_IsoGas;
1022: ctx->physics.destroy = PhysicsDestroy_SimpleFree;
1023: ctx->physics.user = user;
1024: ctx->physics.dof = 2;
1026: PetscStrallocpy("density",&ctx->physics.fieldname[0]);
1027: PetscStrallocpy("momentum",&ctx->physics.fieldname[1]);
1029: user->acoustic_speed = 1;
1031: RiemannListAdd(&rlist,"exact", PhysicsRiemann_IsoGas_Exact);
1032: RiemannListAdd(&rlist,"roe", PhysicsRiemann_IsoGas_Roe);
1033: RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_IsoGas_Rusanov);
1034: ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_IsoGas);
1035: ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative);
1036: PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for IsoGas","");
1037: {
1038: PetscOptionsReal("-physics_isogas_acoustic_speed","Acoustic speed","",user->acoustic_speed,&user->acoustic_speed,NULL);
1039: PetscOptionsFList("-physics_isogas_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
1040: PetscOptionsFList("-physics_isogas_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);
1041: }
1042: PetscOptionsEnd();
1043: RiemannListFind(rlist,rname,&ctx->physics.riemann);
1044: ReconstructListFind(rclist,rcname,&ctx->physics.characteristic);
1045: PetscFunctionListDestroy(&rlist);
1046: PetscFunctionListDestroy(&rclist);
1047: return(0);
1048: }
1052: /* --------------------------------- Shallow Water ----------------------------------- */
1054: typedef struct {
1055: PetscReal gravity;
1056: } ShallowCtx;
1058: PETSC_STATIC_INLINE void ShallowFlux(ShallowCtx *phys,const PetscScalar *u,PetscScalar *f)
1059: {
1060: f[0] = u[1];
1061: f[1] = PetscSqr(u[1])/u[0] + 0.5*phys->gravity*PetscSqr(u[0]);
1062: }
1066: static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
1067: {
1068: ShallowCtx *phys = (ShallowCtx*)vctx;
1069: PetscScalar g = phys->gravity,ustar[2],cL,cR,c,cstar;
1070: struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
1071: PetscInt i;
1072: PetscErrorCode ierr;
1075: if (!(L.h > 0 && R.h > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed thickness is negative");
1076: cL = PetscSqrtScalar(g*L.h);
1077: cR = PetscSqrtScalar(g*R.h);
1078: c = PetscMax(cL,cR);
1079: {
1080: /* Solve for star state */
1081: const PetscInt maxits = 50;
1082: PetscScalar tmp,res,res0=0,h0,h = 0.5*(L.h + R.h); /* initial guess */
1083: h0 = h;
1084: for (i=0; i<maxits; i++) {
1085: PetscScalar fr,fl,dfr,dfl;
1086: fl = (L.h < h)
1087: ? PetscSqrtScalar(0.5*g*(h*h - L.h*L.h)*(1/L.h - 1/h)) /* shock */
1088: : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*L.h); /* rarefaction */
1089: fr = (R.h < h)
1090: ? PetscSqrtScalar(0.5*g*(h*h - R.h*R.h)*(1/R.h - 1/h)) /* shock */
1091: : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*R.h); /* rarefaction */
1092: res = R.u - L.u + fr + fl;
1093: PetscIsInfOrNanScalar(res);
1094: if (PetscAbsScalar(res) < 1e-8 || (i > 0 && PetscAbsScalar(h-h0) < 1e-8)) {
1095: star.h = h;
1096: star.u = L.u - fl;
1097: goto converged;
1098: } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) { /* Line search */
1099: h = 0.8*h0 + 0.2*h;
1100: continue;
1101: }
1102: /* Accept the last step and take another */
1103: res0 = res;
1104: h0 = h;
1105: dfl = (L.h < h)
1106: ? 0.5/fl*0.5*g*(-L.h*L.h/(h*h) - 1 + 2*h/L.h)
1107: : PetscSqrtScalar(g/h);
1108: dfr = (R.h < h)
1109: ? 0.5/fr*0.5*g*(-R.h*R.h/(h*h) - 1 + 2*h/R.h)
1110: : PetscSqrtScalar(g/h);
1111: tmp = h - res/(dfr+dfl);
1112: if (tmp <= 0) h /= 2; /* Guard against Newton shooting off to a negative thickness */
1113: else h = tmp;
1114: if (!((h > 0) && PetscIsNormalScalar(h))) SETERRQ1(PETSC_COMM_SELF,1,"non-normal iterate h=%g",(double)h);
1115: }
1116: SETERRQ1(PETSC_COMM_SELF,1,"Newton iteration for star.h diverged after %D iterations",i);
1117: }
1118: converged:
1119: cstar = PetscSqrtScalar(g*star.h);
1120: if (L.u-cL < 0 && 0 < star.u-cstar) { /* 1-wave is sonic rarefaction */
1121: PetscScalar ufan[2];
1122: ufan[0] = 1/g*PetscSqr(L.u/3 + 2./3*cL);
1123: ufan[1] = PetscSqrtScalar(g*ufan[0])*ufan[0];
1124: ShallowFlux(phys,ufan,flux);
1125: } else if (star.u+cstar < 0 && 0 < R.u+cR) { /* 2-wave is sonic rarefaction */
1126: PetscScalar ufan[2];
1127: ufan[0] = 1/g*PetscSqr(R.u/3 - 2./3*cR);
1128: ufan[1] = -PetscSqrtScalar(g*ufan[0])*ufan[0];
1129: ShallowFlux(phys,ufan,flux);
1130: } else if ((L.h >= star.h && L.u-c >= 0)
1131: || (L.h<star.h && (star.h*star.u-L.h*L.u)/(star.h-L.h) > 0)) {
1132: /* 1-wave is right-travelling shock (supersonic) */
1133: ShallowFlux(phys,uL,flux);
1134: } else if ((star.h <= R.h && R.u+c <= 0)
1135: || (star.h>R.h && (R.h*R.u-star.h*star.h)/(R.h-star.h) < 0)) {
1136: /* 2-wave is left-travelling shock (supersonic) */
1137: ShallowFlux(phys,uR,flux);
1138: } else {
1139: ustar[0] = star.h;
1140: ustar[1] = star.h*star.u;
1141: ShallowFlux(phys,ustar,flux);
1142: }
1143: *maxspeed = MaxAbs(MaxAbs(star.u-cstar,star.u+cstar),MaxAbs(L.u-cL,R.u+cR));
1144: return(0);
1145: }
1149: static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
1150: {
1151: ShallowCtx *phys = (ShallowCtx*)vctx;
1152: PetscScalar g = phys->gravity,fL[2],fR[2],s;
1153: struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};
1156: if (!(L.h > 0 && R.h > 0)) SETERRQ(PETSC_COMM_SELF,1,"Reconstructed thickness is negative");
1157: ShallowFlux(phys,uL,fL);
1158: ShallowFlux(phys,uR,fR);
1159: s = PetscMax(PetscAbs(L.u)+PetscSqrtScalar(g*L.h),PetscAbs(R.u)+PetscSqrtScalar(g*R.h));
1160: flux[0] = 0.5*(fL[0] + fR[0]) + 0.5*s*(uL[0] - uR[0]);
1161: flux[1] = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
1162: *maxspeed = s;
1163: return(0);
1164: }
1168: static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
1169: {
1170: ShallowCtx *phys = (ShallowCtx*)vctx;
1171: PetscReal c;
1175: c = PetscSqrtScalar(u[0]*phys->gravity);
1176: speeds[0] = u[1]/u[0] - c;
1177: speeds[1] = u[1]/u[0] + c;
1178: X[0*2+0] = 1;
1179: X[0*2+1] = speeds[0];
1180: X[1*2+0] = 1;
1181: X[1*2+1] = speeds[1];
1183: PetscMemcpy(Xi,X,4*sizeof(X[0]));
1184: PetscKernel_A_gets_inverse_A_2(Xi,0);
1185: return(0);
1186: }
1190: static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx)
1191: {
1192: PetscErrorCode ierr;
1193: ShallowCtx *user;
1194: PetscFunctionList rlist = 0,rclist = 0;
1195: char rname[256] = "exact",rcname[256] = "characteristic";
1198: PetscNew(&user);
1199: /* Shallow water and Isothermal Gas dynamics are similar so we reuse initial conditions for now */
1200: ctx->physics.sample = PhysicsSample_IsoGas;
1201: ctx->physics.destroy = PhysicsDestroy_SimpleFree;
1202: ctx->physics.user = user;
1203: ctx->physics.dof = 2;
1205: PetscStrallocpy("density",&ctx->physics.fieldname[0]);
1206: PetscStrallocpy("momentum",&ctx->physics.fieldname[1]);
1208: user->gravity = 1;
1210: RiemannListAdd(&rlist,"exact", PhysicsRiemann_Shallow_Exact);
1211: RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Shallow_Rusanov);
1212: ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_Shallow);
1213: ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative);
1214: PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Shallow","");
1215: {
1216: PetscOptionsReal("-physics_shallow_gravity","Gravity","",user->gravity,&user->gravity,NULL);
1217: PetscOptionsFList("-physics_shallow_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL);
1218: PetscOptionsFList("-physics_shallow_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL);
1219: }
1220: PetscOptionsEnd();
1221: RiemannListFind(rlist,rname,&ctx->physics.riemann);
1222: ReconstructListFind(rclist,rcname,&ctx->physics.characteristic);
1223: PetscFunctionListDestroy(&rlist);
1224: PetscFunctionListDestroy(&rclist);
1225: return(0);
1226: }
1230: /* --------------------------------- Finite Volume Solver ----------------------------------- */
1234: static PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
1235: {
1236: FVCtx *ctx = (FVCtx*)vctx;
1238: PetscInt i,j,k,Mx,dof,xs,xm;
1239: PetscReal hx,cfl_idt = 0;
1240: PetscScalar *x,*f,*slope;
1241: Vec Xloc;
1242: DM da;
1245: TSGetDM(ts,&da);
1246: DMGetLocalVector(da,&Xloc);
1247: DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1248: hx = (ctx->xmax - ctx->xmin)/Mx;
1249: DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
1250: DMGlobalToLocalEnd (da,X,INSERT_VALUES,Xloc);
1252: VecZeroEntries(F);
1254: DMDAVecGetArray(da,Xloc,&x);
1255: DMDAVecGetArray(da,F,&f);
1256: DMDAGetArray(da,PETSC_TRUE,&slope);
1258: DMDAGetCorners(da,&xs,0,0,&xm,0,0);
1260: if (ctx->bctype == FVBC_OUTFLOW) {
1261: for (i=xs-2; i<0; i++) {
1262: for (j=0; j<dof; j++) x[i*dof+j] = x[j];
1263: }
1264: for (i=Mx; i<xs+xm+2; i++) {
1265: for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
1266: }
1267: }
1268: for (i=xs-1; i<xs+xm+1; i++) {
1269: struct _LimitInfo info;
1270: PetscScalar *cjmpL,*cjmpR;
1271: /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
1272: (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);
1273: /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
1274: PetscMemzero(ctx->cjmpLR,2*dof*sizeof(ctx->cjmpLR[0]));
1275: cjmpL = &ctx->cjmpLR[0];
1276: cjmpR = &ctx->cjmpLR[dof];
1277: for (j=0; j<dof; j++) {
1278: PetscScalar jmpL,jmpR;
1279: jmpL = x[(i+0)*dof+j] - x[(i-1)*dof+j];
1280: jmpR = x[(i+1)*dof+j] - x[(i+0)*dof+j];
1281: for (k=0; k<dof; k++) {
1282: cjmpL[k] += ctx->Rinv[k+j*dof] * jmpL;
1283: cjmpR[k] += ctx->Rinv[k+j*dof] * jmpR;
1284: }
1285: }
1286: /* Apply limiter to the left and right characteristic jumps */
1287: info.m = dof;
1288: info.hx = hx;
1289: (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
1290: for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
1291: for (j=0; j<dof; j++) {
1292: PetscScalar tmp = 0;
1293: for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof] * ctx->cslope[k];
1294: slope[i*dof+j] = tmp;
1295: }
1296: }
1298: for (i=xs; i<xs+xm+1; i++) {
1299: PetscReal maxspeed;
1300: PetscScalar *uL,*uR;
1301: uL = &ctx->uLR[0];
1302: uR = &ctx->uLR[dof];
1303: for (j=0; j<dof; j++) {
1304: uL[j] = x[(i-1)*dof+j] + slope[(i-1)*dof+j]*hx/2;
1305: uR[j] = x[(i-0)*dof+j] - slope[(i-0)*dof+j]*hx/2;
1306: }
1307: (*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed);
1308: cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
1310: if (i > xs) {
1311: for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx;
1312: }
1313: if (i < xs+xm) {
1314: for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx;
1315: }
1316: }
1318: DMDAVecRestoreArray(da,Xloc,&x);
1319: DMDAVecRestoreArray(da,F,&f);
1320: DMDARestoreArray(da,PETSC_TRUE,&slope);
1321: DMRestoreLocalVector(da,&Xloc);
1323: MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));
1324: if (0) {
1325: /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
1326: PetscReal dt,tnow;
1327: TSGetTimeStep(ts,&dt);
1328: TSGetTime(ts,&tnow);
1329: if (dt > 0.5/ctx->cfl_idt) {
1330: if (1) {
1331: PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt));
1332: } else SETERRQ2(PETSC_COMM_SELF,1,"Stability constraint exceeded, %g > %g",(double)dt,(double)(ctx->cfl/ctx->cfl_idt));
1333: }
1334: }
1335: return(0);
1336: }
1340: static PetscErrorCode SmallMatMultADB(PetscScalar *C,PetscInt bs,const PetscScalar *A,const PetscReal *D,const PetscScalar *B)
1341: {
1342: PetscInt i,j,k;
1345: for (i=0; i<bs; i++) {
1346: for (j=0; j<bs; j++) {
1347: PetscScalar tmp = 0;
1348: for (k=0; k<bs; k++) tmp += A[i*bs+k] * D[k] * B[k*bs+j];
1349: C[i*bs+j] = tmp;
1350: }
1351: }
1352: return(0);
1353: }
1358: static PetscErrorCode FVIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *vctx)
1359: {
1360: FVCtx *ctx = (FVCtx*)vctx;
1362: PetscInt i,j,dof = ctx->physics.dof;
1363: PetscScalar *x,*J;
1364: PetscReal hx;
1365: DM da;
1366: DMDALocalInfo dainfo;
1369: TSGetDM(ts,&da);
1370: DMDAVecGetArray(da,X,&x);
1371: DMDAGetLocalInfo(da,&dainfo);
1372: hx = (ctx->xmax - ctx->xmin)/dainfo.mx;
1373: PetscMalloc1(dof*dof,&J);
1374: for (i=dainfo.xs; i<dainfo.xs+dainfo.xm; i++) {
1375: (*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds);
1376: for (j=0; j<dof; j++) ctx->speeds[j] = PetscAbs(ctx->speeds[j]);
1377: SmallMatMultADB(J,dof,ctx->R,ctx->speeds,ctx->Rinv);
1378: for (j=0; j<dof*dof; j++) J[j] = J[j]/hx + shift*(j/dof == j%dof);
1379: MatSetValuesBlocked(*B,1,&i,1,&i,J,INSERT_VALUES);
1380: }
1381: PetscFree(J);
1382: DMDAVecRestoreArray(da,X,&x);
1384: MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1385: MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1386: if (*A != *B) {
1387: MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);
1388: MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);
1389: }
1390: return(0);
1391: }
1395: static PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U)
1396: {
1398: PetscScalar *u,*uj;
1399: PetscInt i,j,k,dof,xs,xm,Mx;
1402: if (!ctx->physics.sample) SETERRQ(PETSC_COMM_SELF,1,"Physics has not provided a sampling function");
1403: DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1404: DMDAGetCorners(da,&xs,0,0,&xm,0,0);
1405: DMDAVecGetArray(da,U,&u);
1406: PetscMalloc1(dof,&uj);
1407: for (i=xs; i<xs+xm; i++) {
1408: const PetscReal h = (ctx->xmax-ctx->xmin)/Mx,xi = ctx->xmin+h/2+i*h;
1409: const PetscInt N = 200;
1410: /* Integrate over cell i using trapezoid rule with N points. */
1411: for (k=0; k<dof; k++) u[i*dof+k] = 0;
1412: for (j=0; j<N+1; j++) {
1413: PetscScalar xj = xi+h*(j-N/2)/(PetscReal)N;
1414: (*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj);
1415: for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
1416: }
1417: }
1418: DMDAVecRestoreArray(da,U,&u);
1419: PetscFree(uj);
1420: return(0);
1421: }
1425: static PetscErrorCode SolutionStatsView(DM da,Vec X,PetscViewer viewer)
1426: {
1428: PetscReal xmin,xmax;
1429: PetscScalar sum,*x,tvsum,tvgsum;
1430: PetscInt imin,imax,Mx,i,j,xs,xm,dof;
1431: Vec Xloc;
1432: PetscBool iascii;
1435: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1436: if (iascii) {
1437: /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */
1438: DMGetLocalVector(da,&Xloc);
1439: DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
1440: DMGlobalToLocalEnd (da,X,INSERT_VALUES,Xloc);
1441: DMDAVecGetArray(da,Xloc,&x);
1442: DMDAGetCorners(da,&xs,0,0,&xm,0,0);
1443: DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1444: tvsum = 0;
1445: for (i=xs; i<xs+xm; i++) {
1446: for (j=0; j<dof; j++) tvsum += PetscAbsScalar(x[i*dof+j] - x[(i-1)*dof+j]);
1447: }
1448: MPI_Allreduce(&tvsum,&tvgsum,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da));
1449: DMDAVecRestoreArray(da,Xloc,&x);
1450: DMRestoreLocalVector(da,&Xloc);
1452: VecMin(X,&imin,&xmin);
1453: VecMax(X,&imax,&xmax);
1454: VecSum(X,&sum);
1455: PetscViewerASCIIPrintf(viewer,"Solution range [%8.5f,%8.5f] with extrema at %D and %D, mean %8.5f, ||x||_TV %8.5f\n",(double)xmin,(double)xmax,imin,imax,(double)(sum/Mx),(double)(tvgsum/Mx));
1456: } else SETERRQ(PETSC_COMM_SELF,1,"Viewer type not supported");
1457: return(0);
1458: }
1462: static PetscErrorCode SolutionErrorNorms(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1,PetscReal *nrmsup)
1463: {
1465: Vec Y;
1466: PetscInt Mx;
1469: VecGetSize(X,&Mx);
1470: VecDuplicate(X,&Y);
1471: FVSample(ctx,da,t,Y);
1472: VecAYPX(Y,-1,X);
1473: VecNorm(Y,NORM_1,nrm1);
1474: VecNorm(Y,NORM_INFINITY,nrmsup);
1475: *nrm1 /= Mx;
1476: VecDestroy(&Y);
1477: return(0);
1478: }
1482: int main(int argc,char *argv[])
1483: {
1484: char lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
1485: PetscFunctionList limiters = 0,physics = 0;
1486: MPI_Comm comm;
1487: TS ts;
1488: DM da;
1489: Vec X,X0,R;
1490: Mat B;
1491: FVCtx ctx;
1492: PetscInt i,dof,xs,xm,Mx,draw = 0;
1493: PetscBool view_final = PETSC_FALSE;
1494: PetscReal ptime;
1495: PetscErrorCode ierr;
1497: PetscInitialize(&argc,&argv,0,help);
1498: comm = PETSC_COMM_WORLD;
1499: PetscMemzero(&ctx,sizeof(ctx));
1501: /* Register limiters to be available on the command line */
1502: PetscFunctionListAdd(&limiters,"upwind" ,Limit_Upwind);
1503: PetscFunctionListAdd(&limiters,"lax-wendroff" ,Limit_LaxWendroff);
1504: PetscFunctionListAdd(&limiters,"beam-warming" ,Limit_BeamWarming);
1505: PetscFunctionListAdd(&limiters,"fromm" ,Limit_Fromm);
1506: PetscFunctionListAdd(&limiters,"minmod" ,Limit_Minmod);
1507: PetscFunctionListAdd(&limiters,"superbee" ,Limit_Superbee);
1508: PetscFunctionListAdd(&limiters,"mc" ,Limit_MC);
1509: PetscFunctionListAdd(&limiters,"vanleer" ,Limit_VanLeer);
1510: PetscFunctionListAdd(&limiters,"vanalbada" ,Limit_VanAlbada);
1511: PetscFunctionListAdd(&limiters,"vanalbadatvd" ,Limit_VanAlbadaTVD);
1512: PetscFunctionListAdd(&limiters,"koren" ,Limit_Koren);
1513: PetscFunctionListAdd(&limiters,"korensym" ,Limit_KorenSym);
1514: PetscFunctionListAdd(&limiters,"koren3" ,Limit_Koren3);
1515: PetscFunctionListAdd(&limiters,"cada-torrilhon2" ,Limit_CadaTorrilhon2);
1516: PetscFunctionListAdd(&limiters,"cada-torrilhon3-r0p1",Limit_CadaTorrilhon3R0p1);
1517: PetscFunctionListAdd(&limiters,"cada-torrilhon3-r1" ,Limit_CadaTorrilhon3R1);
1518: PetscFunctionListAdd(&limiters,"cada-torrilhon3-r10" ,Limit_CadaTorrilhon3R10);
1519: PetscFunctionListAdd(&limiters,"cada-torrilhon3-r100",Limit_CadaTorrilhon3R100);
1521: /* Register physical models to be available on the command line */
1522: PetscFunctionListAdd(&physics,"advect" ,PhysicsCreate_Advect);
1523: PetscFunctionListAdd(&physics,"burgers" ,PhysicsCreate_Burgers);
1524: PetscFunctionListAdd(&physics,"traffic" ,PhysicsCreate_Traffic);
1525: PetscFunctionListAdd(&physics,"acoustics" ,PhysicsCreate_Acoustics);
1526: PetscFunctionListAdd(&physics,"isogas" ,PhysicsCreate_IsoGas);
1527: PetscFunctionListAdd(&physics,"shallow" ,PhysicsCreate_Shallow);
1529: ctx.comm = comm;
1530: ctx.cfl = 0.9; ctx.bctype = FVBC_PERIODIC;
1531: ctx.xmin = -1; ctx.xmax = 1;
1532: PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");
1533: {
1534: PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL);
1535: PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL);
1536: PetscOptionsFList("-limit","Name of flux limiter to use","",limiters,lname,lname,sizeof(lname),NULL);
1537: PetscOptionsFList("-physics","Name of physics (Riemann solver and characteristics) to use","",physics,physname,physname,sizeof(physname),NULL);
1538: PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL);
1539: PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final);
1540: PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL);
1541: PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL);
1542: PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL);
1543: PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL);
1544: }
1545: PetscOptionsEnd();
1547: /* Choose the limiter from the list of registered limiters */
1548: PetscFunctionListFind(limiters,lname,&ctx.limit);
1549: if (!ctx.limit) SETERRQ1(PETSC_COMM_SELF,1,"Limiter '%s' not found",lname);
1551: /* Choose the physics from the list of registered models */
1552: {
1553: PetscErrorCode (*r)(FVCtx*);
1554: PetscFunctionListFind(physics,physname,&r);
1555: if (!r) SETERRQ1(PETSC_COMM_SELF,1,"Physics '%s' not found",physname);
1556: /* Create the physics, will set the number of fields and their names */
1557: (*r)(&ctx);
1558: }
1560: /* Create a DMDA to manage the parallel grid */
1561: DMDACreate1d(comm,DMDA_BOUNDARY_PERIODIC,-50,ctx.physics.dof,2,NULL,&da);
1562: /* Inform the DMDA of the field names provided by the physics. */
1563: /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
1564: for (i=0; i<ctx.physics.dof; i++) {
1565: DMDASetFieldName(da,i,ctx.physics.fieldname[i]);
1566: }
1567: DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0);
1568: DMDAGetCorners(da,&xs,0,0,&xm,0,0);
1570: /* Set coordinates of cell centers */
1571: DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0);
1573: /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
1574: PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope);
1575: PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds);
1577: /* Create a vector to store the solution and to save the initial state */
1578: DMCreateGlobalVector(da,&X);
1579: VecDuplicate(X,&X0);
1580: VecDuplicate(X,&R);
1582: DMCreateMatrix(da,&B);
1584: /* Create a time-stepping object */
1585: TSCreate(comm,&ts);
1586: TSSetDM(ts,da);
1587: TSSetRHSFunction(ts,R,FVRHSFunction,&ctx);
1588: TSSetIJacobian(ts,B,B,FVIJacobian,&ctx);
1589: TSSetType(ts,TSSSP);
1590: TSSetDuration(ts,1000,10);
1592: /* Compute initial conditions and starting time step */
1593: FVSample(&ctx,da,0,X0);
1594: FVRHSFunction(ts,0,X0,X,(void*)&ctx); /* Initial function evaluation, only used to determine max speed */
1595: VecCopy(X0,X); /* The function value was not used so we set X=X0 again */
1596: TSSetInitialTimeStep(ts,0,ctx.cfl/ctx.cfl_idt);
1598: TSSetFromOptions(ts); /* Take runtime options */
1600: SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);
1602: {
1603: PetscReal nrm1,nrmsup;
1604: PetscInt steps;
1606: TSSolve(ts,X);
1607: TSGetSolveTime(ts,&ptime);
1608: TSGetTimeStepNumber(ts,&steps);
1610: PetscPrintf(comm,"Final time %8.5f, steps %D\n",(double)ptime,steps);
1611: if (ctx.exact) {
1612: SolutionErrorNorms(&ctx,da,ptime,X,&nrm1,&nrmsup);
1613: PetscPrintf(comm,"Error ||x-x_e||_1 %8.4e ||x-x_e||_sup %8.4e\n",(double)nrm1,(double)nrmsup);
1614: }
1615: }
1617: SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD);
1619: if (draw & 0x1) {VecView(X0,PETSC_VIEWER_DRAW_WORLD);}
1620: if (draw & 0x2) {VecView(X,PETSC_VIEWER_DRAW_WORLD);}
1621: if (draw & 0x4) {
1622: Vec Y;
1623: VecDuplicate(X,&Y);
1624: FVSample(&ctx,da,ptime,Y);
1625: VecAYPX(Y,-1,X);
1626: VecView(Y,PETSC_VIEWER_DRAW_WORLD);
1627: VecDestroy(&Y);
1628: }
1630: if (view_final) {
1631: PetscViewer viewer;
1632: PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer);
1633: PetscViewerSetFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
1634: VecView(X,viewer);
1635: PetscViewerDestroy(&viewer);
1636: }
1638: /* Clean up */
1639: (*ctx.physics.destroy)(ctx.physics.user);
1640: for (i=0; i<ctx.physics.dof; i++) {PetscFree(ctx.physics.fieldname[i]);}
1641: PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope);
1642: PetscFree3(ctx.uLR,ctx.flux,ctx.speeds);
1643: VecDestroy(&X);
1644: VecDestroy(&X0);
1645: VecDestroy(&R);
1646: MatDestroy(&B);
1647: DMDestroy(&da);
1648: TSDestroy(&ts);
1649: PetscFunctionListDestroy(&limiters);
1650: PetscFunctionListDestroy(&physics);
1651: PetscFinalize();
1652: return 0;
1653: }