xref: /petsc/src/ts/tutorials/ex9.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
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";
32 
33 #include <petscts.h>
34 #include <petscdm.h>
35 #include <petscdmda.h>
36 #include <petscdraw.h>
37 
38 #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */
39 
40 static inline PetscReal Sgn(PetscReal a) { return (a<0) ? -1 : 1; }
41 static inline PetscReal Abs(PetscReal a) { return (a<0) ? 0 : a; }
42 static inline PetscReal Sqr(PetscReal a) { return a*a; }
43 static inline PetscReal MaxAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) > PetscAbs(b)) ? a : b; }
44 PETSC_UNUSED static inline PetscReal MinAbs(PetscReal a,PetscReal b) { return (PetscAbs(a) < PetscAbs(b)) ? a : b; }
45 static inline PetscReal MinMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMin(PetscAbs(a),PetscAbs(b)); }
46 static inline PetscReal MaxMod2(PetscReal a,PetscReal b) { return (a*b<0) ? 0 : Sgn(a)*PetscMax(PetscAbs(a),PetscAbs(b)); }
47 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 
49 static inline PetscReal RangeMod(PetscReal a,PetscReal xmin,PetscReal xmax) { PetscReal range = xmax-xmin; return xmin +PetscFmodReal(range+PetscFmodReal(a,range),range); }
50 
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 : (jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15);
105 }
106 static void Limit_Koren(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
107 { /* phi = (t + 2*t^2) / (2 - t + 2*t^2) */
108   PetscInt i;
109   for (i=0; i<info->m; i++) lmt[i] = ((jL[i]*Sqr(jR[i]) + 2*Sqr(jL[i])*jR[i])/(2*Sqr(jL[i]) - jL[i]*jR[i] + 2*Sqr(jR[i]) + 1e-15));
110 }
111 static void Limit_KorenSym(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt) /* differentiable */
112 { /* Symmetric version of above */
113   PetscInt i;
114   for (i=0; i<info->m; i++) lmt[i] = (1.5*(jL[i]*Sqr(jR[i]) + Sqr(jL[i])*jR[i])/(2*Sqr(jL[i]) - jL[i]*jR[i] + 2*Sqr(jR[i]) + 1e-15));
115 }
116 static void Limit_Koren3(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
117 { /* Eq 11 of Cada-Torrilhon 2009 */
118   PetscInt i;
119   for (i=0; i<info->m; i++) lmt[i] = MinMod3(2*jL[i],(jL[i]+2*jR[i])/3,2*jR[i]);
120 }
121 static PetscReal CadaTorrilhonPhiHatR_Eq13(PetscReal L,PetscReal R)
122 {
123   return PetscMax(0,PetscMin((L+2*R)/3,PetscMax(-0.5*L,PetscMin(2*L,PetscMin((L+2*R)/3,1.6*R)))));
124 }
125 static void Limit_CadaTorrilhon2(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
126 { /* Cada-Torrilhon 2009, Eq 13 */
127   PetscInt i;
128   for (i=0; i<info->m; i++) lmt[i] = CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]);
129 }
130 static void Limit_CadaTorrilhon3R(PetscReal r,LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
131 { /* Cada-Torrilhon 2009, Eq 22 */
132   /* They recommend 0.001 < r < 1, but larger values are more accurate in smooth regions */
133   const PetscReal eps = 1e-7,hx = info->hx;
134   PetscInt        i;
135   for (i=0; i<info->m; i++) {
136     const PetscReal eta = (Sqr(jL[i]) + Sqr(jR[i])) / Sqr(r*hx);
137     lmt[i] = ((eta < 1-eps) ? (jL[i] + 2*jR[i]) / 3 : ((eta > 1+eps) ? CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]) : 0.5*((1-(eta-1)/eps)*(jL[i]+2*jR[i])/3 + (1+(eta+1)/eps)*CadaTorrilhonPhiHatR_Eq13(jL[i],jR[i]))));
138   }
139 }
140 static void Limit_CadaTorrilhon3R0p1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
141 {
142   Limit_CadaTorrilhon3R(0.1,info,jL,jR,lmt);
143 }
144 static void Limit_CadaTorrilhon3R1(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
145 {
146   Limit_CadaTorrilhon3R(1,info,jL,jR,lmt);
147 }
148 static void Limit_CadaTorrilhon3R10(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
149 {
150   Limit_CadaTorrilhon3R(10,info,jL,jR,lmt);
151 }
152 static void Limit_CadaTorrilhon3R100(LimitInfo info,const PetscScalar *jL,const PetscScalar *jR,PetscScalar *lmt)
153 {
154   Limit_CadaTorrilhon3R(100,info,jL,jR,lmt);
155 }
156 
157 /* --------------------------------- Finite Volume data structures ----------------------------------- */
158 
159 typedef enum {FVBC_PERIODIC, FVBC_OUTFLOW} FVBCType;
160 static const char *FVBCTypes[] = {"PERIODIC","OUTFLOW","FVBCType","FVBC_",0};
161 typedef PetscErrorCode (*RiemannFunction)(void*,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscReal*);
162 typedef PetscErrorCode (*ReconstructFunction)(void*,PetscInt,const PetscScalar*,PetscScalar*,PetscScalar*,PetscReal*);
163 
164 typedef struct {
165   PetscErrorCode      (*sample)(void*,PetscInt,FVBCType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal*);
166   RiemannFunction     riemann;
167   ReconstructFunction characteristic;
168   PetscErrorCode      (*destroy)(void*);
169   void                *user;
170   PetscInt            dof;
171   char                *fieldname[16];
172 } PhysicsCtx;
173 
174 typedef struct {
175   void        (*limit)(LimitInfo,const PetscScalar*,const PetscScalar*,PetscScalar*);
176   PhysicsCtx  physics;
177   MPI_Comm    comm;
178   char        prefix[256];
179 
180   /* Local work arrays */
181   PetscScalar *R,*Rinv;         /* Characteristic basis, and it's inverse.  COLUMN-MAJOR */
182   PetscScalar *cjmpLR;          /* Jumps at left and right edge of cell, in characteristic basis, len=2*dof */
183   PetscScalar *cslope;          /* Limited slope, written in characteristic basis */
184   PetscScalar *uLR;             /* Solution at left and right of interface, conservative variables, len=2*dof */
185   PetscScalar *flux;            /* Flux across interface */
186   PetscReal   *speeds;          /* Speeds of each wave */
187 
188   PetscReal   cfl_idt;            /* Max allowable value of 1/Delta t */
189   PetscReal   cfl;
190   PetscReal   xmin,xmax;
191   PetscInt    initial;
192   PetscBool   exact;
193   FVBCType    bctype;
194 } FVCtx;
195 
196 PetscErrorCode RiemannListAdd(PetscFunctionList *flist,const char *name,RiemannFunction rsolve)
197 {
198   PetscFunctionBeginUser;
199   PetscCall(PetscFunctionListAdd(flist,name,rsolve));
200   PetscFunctionReturn(0);
201 }
202 
203 PetscErrorCode RiemannListFind(PetscFunctionList flist,const char *name,RiemannFunction *rsolve)
204 {
205   PetscFunctionBeginUser;
206   PetscCall(PetscFunctionListFind(flist,name,rsolve));
207   PetscCheck(*rsolve,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Riemann solver \"%s\" could not be found",name);
208   PetscFunctionReturn(0);
209 }
210 
211 PetscErrorCode ReconstructListAdd(PetscFunctionList *flist,const char *name,ReconstructFunction r)
212 {
213   PetscFunctionBeginUser;
214   PetscCall(PetscFunctionListAdd(flist,name,r));
215   PetscFunctionReturn(0);
216 }
217 
218 PetscErrorCode ReconstructListFind(PetscFunctionList flist,const char *name,ReconstructFunction *r)
219 {
220   PetscFunctionBeginUser;
221   PetscCall(PetscFunctionListFind(flist,name,r));
222   PetscCheck(*r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Reconstruction \"%s\" could not be found",name);
223   PetscFunctionReturn(0);
224 }
225 
226 /* --------------------------------- Physics ----------------------------------- */
227 /*
228   Each physical model consists of Riemann solver and a function to determine the basis to use for reconstruction.  These
229   are set with the PhysicsCreate_XXX function which allocates private storage and sets these methods as well as the
230   number of fields and their names, and a function to deallocate private storage.
231 */
232 
233 /* First a few functions useful to several different physics */
234 static PetscErrorCode PhysicsCharacteristic_Conservative(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
235 {
236   PetscInt i,j;
237 
238   PetscFunctionBeginUser;
239   for (i=0; i<m; i++) {
240     for (j=0; j<m; j++) Xi[i*m+j] = X[i*m+j] = (PetscScalar)(i==j);
241     speeds[i] = PETSC_MAX_REAL; /* Indicates invalid */
242   }
243   PetscFunctionReturn(0);
244 }
245 
246 static PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx)
247 {
248   PetscFunctionBeginUser;
249   PetscCall(PetscFree(vctx));
250   PetscFunctionReturn(0);
251 }
252 
253 /* --------------------------------- Advection ----------------------------------- */
254 
255 typedef struct {
256   PetscReal a;                  /* advective velocity */
257 } AdvectCtx;
258 
259 static PetscErrorCode PhysicsRiemann_Advect(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
260 {
261   AdvectCtx *ctx = (AdvectCtx*)vctx;
262   PetscReal speed;
263 
264   PetscFunctionBeginUser;
265   speed     = ctx->a;
266   flux[0]   = PetscMax(0,speed)*uL[0] + PetscMin(0,speed)*uR[0];
267   *maxspeed = speed;
268   PetscFunctionReturn(0);
269 }
270 
271 static PetscErrorCode PhysicsCharacteristic_Advect(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
272 {
273   AdvectCtx *ctx = (AdvectCtx*)vctx;
274 
275   PetscFunctionBeginUser;
276   X[0]      = 1.;
277   Xi[0]     = 1.;
278   speeds[0] = ctx->a;
279   PetscFunctionReturn(0);
280 }
281 
282 static PetscErrorCode PhysicsSample_Advect(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
283 {
284   AdvectCtx *ctx = (AdvectCtx*)vctx;
285   PetscReal a    = ctx->a,x0;
286 
287   PetscFunctionBeginUser;
288   switch (bctype) {
289     case FVBC_OUTFLOW:   x0 = x-a*t; break;
290     case FVBC_PERIODIC: x0 = RangeMod(x-a*t,xmin,xmax); break;
291     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType");
292   }
293   switch (initial) {
294     case 0: u[0] = (x0 < 0) ? 1 : -1; break;
295     case 1: u[0] = (x0 < 0) ? -1 : 1; break;
296     case 2: u[0] = (0 < x0 && x0 < 1) ? 1 : 0; break;
297     case 3: u[0] = PetscSinReal(2*PETSC_PI*x0); break;
298     case 4: u[0] = PetscAbs(x0); break;
299     case 5: u[0] = (x0 < 0 || x0 > 0.5) ? 0 : PetscSqr(PetscSinReal(2*PETSC_PI*x0)); break;
300     case 6: u[0] = (x0 < 0) ? 0 : ((x0 < 1) ? x0 : ((x0 < 2) ? 2-x0 : 0)); break;
301     case 7: u[0] = PetscPowReal(PetscSinReal(PETSC_PI*x0),10.0);break;
302     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
303   }
304   PetscFunctionReturn(0);
305 }
306 
307 static PetscErrorCode PhysicsCreate_Advect(FVCtx *ctx)
308 {
309   AdvectCtx      *user;
310 
311   PetscFunctionBeginUser;
312   PetscCall(PetscNew(&user));
313   ctx->physics.sample         = PhysicsSample_Advect;
314   ctx->physics.riemann        = PhysicsRiemann_Advect;
315   ctx->physics.characteristic = PhysicsCharacteristic_Advect;
316   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
317   ctx->physics.user           = user;
318   ctx->physics.dof            = 1;
319   PetscCall(PetscStrallocpy("u",&ctx->physics.fieldname[0]));
320   user->a = 1;
321   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");
322   {
323     PetscCall(PetscOptionsReal("-physics_advect_a","Speed","",user->a,&user->a,NULL));
324   }
325   PetscOptionsEnd();
326   PetscFunctionReturn(0);
327 }
328 
329 /* --------------------------------- Burgers ----------------------------------- */
330 
331 typedef struct {
332   PetscReal lxf_speed;
333 } BurgersCtx;
334 
335 static PetscErrorCode PhysicsSample_Burgers(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
336 {
337   PetscFunctionBeginUser;
338   PetscCheck(bctype != FVBC_PERIODIC || t <= 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exact solution not implemented for periodic");
339   switch (initial) {
340     case 0: u[0] = (x < 0) ? 1 : -1; break;
341     case 1:
342       if       (x < -t) u[0] = -1;
343       else if  (x < t)  u[0] = x/t;
344       else              u[0] = 1;
345       break;
346     case 2:
347       if      (x <= 0)      u[0] = 0;
348       else if (x < t)       u[0] = x/t;
349       else if (x < 1+0.5*t) u[0] = 1;
350       else                  u[0] = 0;
351       break;
352     case 3:
353       if       (x < 0.2*t) u[0] = 0.2;
354       else if  (x < t) u[0] = x/t;
355       else             u[0] = 1;
356       break;
357     case 4:
358       PetscCheck(t <= 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only initial condition available");
359       u[0] = 0.7 + 0.3*PetscSinReal(2*PETSC_PI*((x-xmin)/(xmax-xmin)));
360       break;
361     case 5:                     /* Pure shock solution */
362       if (x < 0.5*t) u[0] = 1;
363       else u[0] = 0;
364       break;
365     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
366   }
367   PetscFunctionReturn(0);
368 }
369 
370 static PetscErrorCode PhysicsRiemann_Burgers_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
371 {
372   PetscFunctionBeginUser;
373   if (uL[0] < uR[0]) {          /* rarefaction */
374     flux[0] = (uL[0]*uR[0] < 0)
375       ? 0                       /* sonic rarefaction */
376       : 0.5*PetscMin(PetscSqr(uL[0]),PetscSqr(uR[0]));
377   } else {                      /* shock */
378     flux[0] = 0.5*PetscMax(PetscSqr(uL[0]),PetscSqr(uR[0]));
379   }
380   *maxspeed = (PetscAbs(uL[0]) > PetscAbs(uR[0])) ? uL[0] : uR[0];
381   PetscFunctionReturn(0);
382 }
383 
384 static PetscErrorCode PhysicsRiemann_Burgers_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
385 {
386   PetscReal speed;
387 
388   PetscFunctionBeginUser;
389   speed   = 0.5*(uL[0] + uR[0]);
390   flux[0] = 0.25*(PetscSqr(uL[0]) + PetscSqr(uR[0])) - 0.5*PetscAbs(speed)*(uR[0]-uL[0]);
391   if (uL[0] <= 0 && 0 <= uR[0]) flux[0] = 0; /* Entropy fix for sonic rarefaction */
392   *maxspeed = speed;
393   PetscFunctionReturn(0);
394 }
395 
396 static PetscErrorCode PhysicsRiemann_Burgers_LxF(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
397 {
398   PetscReal   c;
399   PetscScalar fL,fR;
400 
401   PetscFunctionBeginUser;
402   c         = ((BurgersCtx*)vctx)->lxf_speed;
403   fL        = 0.5*PetscSqr(uL[0]);
404   fR        = 0.5*PetscSqr(uR[0]);
405   flux[0]   = 0.5*(fL + fR) - 0.5*c*(uR[0] - uL[0]);
406   *maxspeed = c;
407   PetscFunctionReturn(0);
408 }
409 
410 static PetscErrorCode PhysicsRiemann_Burgers_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
411 {
412   PetscReal   c;
413   PetscScalar fL,fR;
414 
415   PetscFunctionBeginUser;
416   c         = PetscMax(PetscAbs(uL[0]),PetscAbs(uR[0]));
417   fL        = 0.5*PetscSqr(uL[0]);
418   fR        = 0.5*PetscSqr(uR[0]);
419   flux[0]   = 0.5*(fL + fR) - 0.5*c*(uR[0] - uL[0]);
420   *maxspeed = c;
421   PetscFunctionReturn(0);
422 }
423 
424 static PetscErrorCode PhysicsCreate_Burgers(FVCtx *ctx)
425 {
426   BurgersCtx        *user;
427   RiemannFunction   r;
428   PetscFunctionList rlist      = 0;
429   char              rname[256] = "exact";
430 
431   PetscFunctionBeginUser;
432   PetscCall(PetscNew(&user));
433 
434   ctx->physics.sample         = PhysicsSample_Burgers;
435   ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
436   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
437   ctx->physics.user           = user;
438   ctx->physics.dof            = 1;
439 
440   PetscCall(PetscStrallocpy("u",&ctx->physics.fieldname[0]));
441   PetscCall(RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Burgers_Exact));
442   PetscCall(RiemannListAdd(&rlist,"roe",    PhysicsRiemann_Burgers_Roe));
443   PetscCall(RiemannListAdd(&rlist,"lxf",    PhysicsRiemann_Burgers_LxF));
444   PetscCall(RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Burgers_Rusanov));
445   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for advection","");
446   {
447     PetscCall(PetscOptionsFList("-physics_burgers_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL));
448   }
449   PetscOptionsEnd();
450   PetscCall(RiemannListFind(rlist,rname,&r));
451   PetscCall(PetscFunctionListDestroy(&rlist));
452   ctx->physics.riemann = r;
453 
454   /* *
455   * Hack to deal with LxF in semi-discrete form
456   * max speed is 1 for the basic initial conditions (where |u| <= 1)
457   * */
458   if (r == PhysicsRiemann_Burgers_LxF) user->lxf_speed = 1;
459   PetscFunctionReturn(0);
460 }
461 
462 /* --------------------------------- Traffic ----------------------------------- */
463 
464 typedef struct {
465   PetscReal lxf_speed;
466   PetscReal a;
467 } TrafficCtx;
468 
469 static inline PetscScalar TrafficFlux(PetscScalar a,PetscScalar u) { return a*u*(1-u); }
470 
471 static PetscErrorCode PhysicsSample_Traffic(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
472 {
473   PetscReal a = ((TrafficCtx*)vctx)->a;
474 
475   PetscFunctionBeginUser;
476   PetscCheck(bctype != FVBC_PERIODIC || t <= 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exact solution not implemented for periodic");
477   switch (initial) {
478     case 0:
479       u[0] = (-a*t < x) ? 2 : 0; break;
480     case 1:
481       if      (x < PetscMin(2*a*t,0.5+a*t)) u[0] = -1;
482       else if (x < 1)                       u[0] = 0;
483       else                                  u[0] = 1;
484       break;
485     case 2:
486       PetscCheck(t <= 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only initial condition available");
487       u[0] = 0.7 + 0.3*PetscSinReal(2*PETSC_PI*((x-xmin)/(xmax-xmin)));
488       break;
489     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
490   }
491   PetscFunctionReturn(0);
492 }
493 
494 static PetscErrorCode PhysicsRiemann_Traffic_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
495 {
496   PetscReal a = ((TrafficCtx*)vctx)->a;
497 
498   PetscFunctionBeginUser;
499   if (uL[0] < uR[0]) {
500     flux[0] = PetscMin(TrafficFlux(a,uL[0]),TrafficFlux(a,uR[0]));
501   } else {
502     flux[0] = (uR[0] < 0.5 && 0.5 < uL[0]) ? TrafficFlux(a,0.5) : PetscMax(TrafficFlux(a,uL[0]),TrafficFlux(a,uR[0]));
503   }
504   *maxspeed = a*MaxAbs(1-2*uL[0],1-2*uR[0]);
505   PetscFunctionReturn(0);
506 }
507 
508 static PetscErrorCode PhysicsRiemann_Traffic_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
509 {
510   PetscReal a = ((TrafficCtx*)vctx)->a;
511   PetscReal speed;
512 
513   PetscFunctionBeginUser;
514   speed = a*(1 - (uL[0] + uR[0]));
515   flux[0] = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*PetscAbs(speed)*(uR[0]-uL[0]);
516   *maxspeed = speed;
517   PetscFunctionReturn(0);
518 }
519 
520 static PetscErrorCode PhysicsRiemann_Traffic_LxF(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
521 {
522   TrafficCtx *phys = (TrafficCtx*)vctx;
523   PetscReal  a     = phys->a;
524   PetscReal  speed;
525 
526   PetscFunctionBeginUser;
527   speed     = a*(1 - (uL[0] + uR[0]));
528   flux[0]   = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*phys->lxf_speed*(uR[0]-uL[0]);
529   *maxspeed = speed;
530   PetscFunctionReturn(0);
531 }
532 
533 static PetscErrorCode PhysicsRiemann_Traffic_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
534 {
535   PetscReal a = ((TrafficCtx*)vctx)->a;
536   PetscReal speed;
537 
538   PetscFunctionBeginUser;
539   speed     = a*PetscMax(PetscAbs(1-2*uL[0]),PetscAbs(1-2*uR[0]));
540   flux[0]   = 0.5*(TrafficFlux(a,uL[0]) + TrafficFlux(a,uR[0])) - 0.5*speed*(uR[0]-uL[0]);
541   *maxspeed = speed;
542   PetscFunctionReturn(0);
543 }
544 
545 static PetscErrorCode PhysicsCreate_Traffic(FVCtx *ctx)
546 {
547   TrafficCtx        *user;
548   RiemannFunction   r;
549   PetscFunctionList rlist      = 0;
550   char              rname[256] = "exact";
551 
552   PetscFunctionBeginUser;
553   PetscCall(PetscNew(&user));
554   ctx->physics.sample         = PhysicsSample_Traffic;
555   ctx->physics.characteristic = PhysicsCharacteristic_Conservative;
556   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
557   ctx->physics.user           = user;
558   ctx->physics.dof            = 1;
559 
560   PetscCall(PetscStrallocpy("density",&ctx->physics.fieldname[0]));
561   user->a = 0.5;
562   PetscCall(RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Traffic_Exact));
563   PetscCall(RiemannListAdd(&rlist,"roe",    PhysicsRiemann_Traffic_Roe));
564   PetscCall(RiemannListAdd(&rlist,"lxf",    PhysicsRiemann_Traffic_LxF));
565   PetscCall(RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Traffic_Rusanov));
566   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Traffic","");
567     PetscCall(PetscOptionsReal("-physics_traffic_a","Flux = a*u*(1-u)","",user->a,&user->a,NULL));
568     PetscCall(PetscOptionsFList("-physics_traffic_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL));
569   PetscOptionsEnd();
570 
571   PetscCall(RiemannListFind(rlist,rname,&r));
572   PetscCall(PetscFunctionListDestroy(&rlist));
573 
574   ctx->physics.riemann = r;
575 
576   /* *
577   * Hack to deal with LxF in semi-discrete form
578   * max speed is 3*a for the basic initial conditions (-1 <= u <= 2)
579   * */
580   if (r == PhysicsRiemann_Traffic_LxF) user->lxf_speed = 3*user->a;
581   PetscFunctionReturn(0);
582 }
583 
584 /* --------------------------------- Linear Acoustics ----------------------------------- */
585 
586 /* Flux: u_t + (A u)_x
587  * z = sqrt(rho*bulk), c = sqrt(rho/bulk)
588  * Spectral decomposition: A = R * D * Rinv
589  * [    cz] = [-z   z] [-c    ] [-1/2z  1/2]
590  * [c/z   ] = [ 1   1] [     c] [ 1/2z  1/2]
591  *
592  * We decompose this into the left-traveling waves Al = R * D^- Rinv
593  * and the right-traveling waves Ar = R * D^+ * Rinv
594  * Multiplying out these expressions produces the following two matrices
595  */
596 
597 typedef struct {
598   PetscReal c;                  /* speed of sound: c = sqrt(bulk/rho) */
599   PetscReal z;                  /* impedence: z = sqrt(rho*bulk) */
600 } AcousticsCtx;
601 
602 PETSC_UNUSED static inline void AcousticsFlux(AcousticsCtx *ctx,const PetscScalar *u,PetscScalar *f)
603 {
604   f[0] = ctx->c*ctx->z*u[1];
605   f[1] = ctx->c/ctx->z*u[0];
606 }
607 
608 static PetscErrorCode PhysicsCharacteristic_Acoustics(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
609 {
610   AcousticsCtx *phys = (AcousticsCtx*)vctx;
611   PetscReal    z     = phys->z,c = phys->c;
612 
613   PetscFunctionBeginUser;
614   X[0*2+0]  = -z;
615   X[0*2+1]  = z;
616   X[1*2+0]  = 1;
617   X[1*2+1]  = 1;
618   Xi[0*2+0] = -1./(2*z);
619   Xi[0*2+1] = 1./2;
620   Xi[1*2+0] = 1./(2*z);
621   Xi[1*2+1] = 1./2;
622   speeds[0] = -c;
623   speeds[1] = c;
624   PetscFunctionReturn(0);
625 }
626 
627 static PetscErrorCode PhysicsSample_Acoustics_Initial(AcousticsCtx *phys,PetscInt initial,PetscReal xmin,PetscReal xmax,PetscReal x,PetscReal *u)
628 {
629   PetscFunctionBeginUser;
630   switch (initial) {
631   case 0:
632     u[0] = (PetscAbs((x - xmin)/(xmax - xmin) - 0.2) < 0.1) ? 1 : 0.5;
633     u[1] = (PetscAbs((x - xmin)/(xmax - xmin) - 0.7) < 0.1) ? 1 : -0.5;
634     break;
635   case 1:
636     u[0] = PetscCosReal(3 * 2*PETSC_PI*x/(xmax-xmin));
637     u[1] = PetscExpReal(-PetscSqr(x - (xmax + xmin)/2) / (2*PetscSqr(0.2*(xmax - xmin)))) - 0.5;
638     break;
639   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
640   }
641   PetscFunctionReturn(0);
642 }
643 
644 static PetscErrorCode PhysicsSample_Acoustics(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
645 {
646   AcousticsCtx   *phys = (AcousticsCtx*)vctx;
647   PetscReal      c     = phys->c;
648   PetscReal      x0a,x0b,u0a[2],u0b[2],tmp[2];
649   PetscReal      X[2][2],Xi[2][2],dummy[2];
650 
651   PetscFunctionBeginUser;
652   switch (bctype) {
653   case FVBC_OUTFLOW:
654     x0a = x+c*t;
655     x0b = x-c*t;
656     break;
657   case FVBC_PERIODIC:
658     x0a = RangeMod(x+c*t,xmin,xmax);
659     x0b = RangeMod(x-c*t,xmin,xmax);
660     break;
661   default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown BCType");
662   }
663   PetscCall(PhysicsSample_Acoustics_Initial(phys,initial,xmin,xmax,x0a,u0a));
664   PetscCall(PhysicsSample_Acoustics_Initial(phys,initial,xmin,xmax,x0b,u0b));
665   PetscCall(PhysicsCharacteristic_Acoustics(vctx,2,u,&X[0][0],&Xi[0][0],dummy));
666   tmp[0] = Xi[0][0]*u0a[0] + Xi[0][1]*u0a[1];
667   tmp[1] = Xi[1][0]*u0b[0] + Xi[1][1]*u0b[1];
668   u[0]   = X[0][0]*tmp[0] + X[0][1]*tmp[1];
669   u[1]   = X[1][0]*tmp[0] + X[1][1]*tmp[1];
670   PetscFunctionReturn(0);
671 }
672 
673 static PetscErrorCode PhysicsRiemann_Acoustics_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
674 {
675   AcousticsCtx *phys = (AcousticsCtx*)vctx;
676   PetscReal    c     = phys->c,z = phys->z;
677   PetscReal
678     Al[2][2] = {{-c/2     , c*z/2  },
679                 {c/(2*z)  , -c/2   }}, /* Left traveling waves */
680     Ar[2][2] = {{c/2      , c*z/2  },
681                 {c/(2*z)  , c/2    }}; /* Right traveling waves */
682 
683   PetscFunctionBeginUser;
684   flux[0]   = Al[0][0]*uR[0] + Al[0][1]*uR[1] + Ar[0][0]*uL[0] + Ar[0][1]*uL[1];
685   flux[1]   = Al[1][0]*uR[0] + Al[1][1]*uR[1] + Ar[1][0]*uL[0] + Ar[1][1]*uL[1];
686   *maxspeed = c;
687   PetscFunctionReturn(0);
688 }
689 
690 static PetscErrorCode PhysicsCreate_Acoustics(FVCtx *ctx)
691 {
692   AcousticsCtx      *user;
693   PetscFunctionList rlist      = 0,rclist = 0;
694   char              rname[256] = "exact",rcname[256] = "characteristic";
695 
696   PetscFunctionBeginUser;
697   PetscCall(PetscNew(&user));
698   ctx->physics.sample         = PhysicsSample_Acoustics;
699   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
700   ctx->physics.user           = user;
701   ctx->physics.dof            = 2;
702 
703   PetscCall(PetscStrallocpy("u",&ctx->physics.fieldname[0]));
704   PetscCall(PetscStrallocpy("v",&ctx->physics.fieldname[1]));
705 
706   user->c = 1;
707   user->z = 1;
708 
709   PetscCall(RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Acoustics_Exact));
710   PetscCall(ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_Acoustics));
711   PetscCall(ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative));
712   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for linear Acoustics","");
713   {
714     PetscCall(PetscOptionsReal("-physics_acoustics_c","c = sqrt(bulk/rho)","",user->c,&user->c,NULL));
715     PetscCall(PetscOptionsReal("-physics_acoustics_z","z = sqrt(bulk*rho)","",user->z,&user->z,NULL));
716     PetscCall(PetscOptionsFList("-physics_acoustics_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL));
717     PetscCall(PetscOptionsFList("-physics_acoustics_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL));
718   }
719   PetscOptionsEnd();
720   PetscCall(RiemannListFind(rlist,rname,&ctx->physics.riemann));
721   PetscCall(ReconstructListFind(rclist,rcname,&ctx->physics.characteristic));
722   PetscCall(PetscFunctionListDestroy(&rlist));
723   PetscCall(PetscFunctionListDestroy(&rclist));
724   PetscFunctionReturn(0);
725 }
726 
727 /* --------------------------------- Isothermal Gas Dynamics ----------------------------------- */
728 
729 typedef struct {
730   PetscReal acoustic_speed;
731 } IsoGasCtx;
732 
733 static inline void IsoGasFlux(PetscReal c,const PetscScalar *u,PetscScalar *f)
734 {
735   f[0] = u[1];
736   f[1] = PetscSqr(u[1])/u[0] + c*c*u[0];
737 }
738 
739 static PetscErrorCode PhysicsSample_IsoGas(void *vctx,PetscInt initial,FVBCType bctype,PetscReal xmin,PetscReal xmax,PetscReal t,PetscReal x,PetscReal *u)
740 {
741   PetscFunctionBeginUser;
742   PetscCheck(t <= 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exact solutions not implemented for t > 0");
743   switch (initial) {
744     case 0:
745       u[0] = (x < 0) ? 1 : 0.5;
746       u[1] = (x < 0) ? 1 : 0.7;
747       break;
748     case 1:
749       u[0] = 1+0.5*PetscSinReal(2*PETSC_PI*x);
750       u[1] = 1*u[0];
751       break;
752     default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"unknown initial condition");
753   }
754   PetscFunctionReturn(0);
755 }
756 
757 static PetscErrorCode PhysicsRiemann_IsoGas_Roe(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
758 {
759   IsoGasCtx   *phys = (IsoGasCtx*)vctx;
760   PetscReal   c     = phys->acoustic_speed;
761   PetscScalar ubar,du[2],a[2],fL[2],fR[2],lam[2],ustar[2],R[2][2];
762   PetscInt    i;
763 
764   PetscFunctionBeginUser;
765   ubar = (uL[1]/PetscSqrtScalar(uL[0]) + uR[1]/PetscSqrtScalar(uR[0])) / (PetscSqrtScalar(uL[0]) + PetscSqrtScalar(uR[0]));
766   /* write fluxuations in characteristic basis */
767   du[0] = uR[0] - uL[0];
768   du[1] = uR[1] - uL[1];
769   a[0]  = (1/(2*c)) * ((ubar + c)*du[0] - du[1]);
770   a[1]  = (1/(2*c)) * ((-ubar + c)*du[0] + du[1]);
771   /* wave speeds */
772   lam[0] = ubar - c;
773   lam[1] = ubar + c;
774   /* Right eigenvectors */
775   R[0][0] = 1; R[0][1] = ubar-c;
776   R[1][0] = 1; R[1][1] = ubar+c;
777   /* Compute state in star region (between the 1-wave and 2-wave) */
778   for (i=0; i<2; i++) ustar[i] = uL[i] + a[0]*R[0][i];
779   if (uL[1]/uL[0] < c && c < ustar[1]/ustar[0]) { /* 1-wave is sonic rarefaction */
780     PetscScalar ufan[2];
781     ufan[0] = uL[0]*PetscExpScalar(uL[1]/(uL[0]*c) - 1);
782     ufan[1] = c*ufan[0];
783     IsoGasFlux(c,ufan,flux);
784   } else if (ustar[1]/ustar[0] < -c && -c < uR[1]/uR[0]) { /* 2-wave is sonic rarefaction */
785     PetscScalar ufan[2];
786     ufan[0] = uR[0]*PetscExpScalar(-uR[1]/(uR[0]*c) - 1);
787     ufan[1] = -c*ufan[0];
788     IsoGasFlux(c,ufan,flux);
789   } else {                      /* Centered form */
790     IsoGasFlux(c,uL,fL);
791     IsoGasFlux(c,uR,fR);
792     for (i=0; i<2; i++) {
793       PetscScalar absdu = PetscAbsScalar(lam[0])*a[0]*R[0][i] + PetscAbsScalar(lam[1])*a[1]*R[1][i];
794       flux[i] = 0.5*(fL[i]+fR[i]) - 0.5*absdu;
795     }
796   }
797   *maxspeed = MaxAbs(lam[0],lam[1]);
798   PetscFunctionReturn(0);
799 }
800 
801 static PetscErrorCode PhysicsRiemann_IsoGas_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
802 {
803   IsoGasCtx                   *phys = (IsoGasCtx*)vctx;
804   PetscReal                   c     = phys->acoustic_speed;
805   PetscScalar                 ustar[2];
806   struct {PetscScalar rho,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
807   PetscInt                    i;
808 
809   PetscFunctionBeginUser;
810   PetscCheck((L.rho > 0 && R.rho > 0),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed density is negative");
811   {
812     /* Solve for star state */
813     PetscScalar res,tmp,rho = 0.5*(L.rho + R.rho); /* initial guess */
814     for (i=0; i<20; i++) {
815       PetscScalar fr,fl,dfr,dfl;
816       fl = (L.rho < rho)
817         ? (rho-L.rho)/PetscSqrtScalar(L.rho*rho)       /* shock */
818         : PetscLogScalar(rho) - PetscLogScalar(L.rho); /* rarefaction */
819       fr = (R.rho < rho)
820         ? (rho-R.rho)/PetscSqrtScalar(R.rho*rho)       /* shock */
821         : PetscLogScalar(rho) - PetscLogScalar(R.rho); /* rarefaction */
822       res = R.u-L.u + c*(fr+fl);
823       PetscCheck(!PetscIsInfOrNanScalar(res),PETSC_COMM_SELF,PETSC_ERR_FP,"Infinity or Not-a-Number generated in computation");
824       if (PetscAbsScalar(res) < 1e-10) {
825         star.rho = rho;
826         star.u   = L.u - c*fl;
827         goto converged;
828       }
829       dfl = (L.rho < rho) ? 1/PetscSqrtScalar(L.rho*rho)*(1 - 0.5*(rho-L.rho)/rho) : 1/rho;
830       dfr = (R.rho < rho) ? 1/PetscSqrtScalar(R.rho*rho)*(1 - 0.5*(rho-R.rho)/rho) : 1/rho;
831       tmp = rho - res/(c*(dfr+dfl));
832       if (tmp <= 0) rho /= 2;   /* Guard against Newton shooting off to a negative density */
833       else rho = tmp;
834       PetscCheck(((rho > 0) && PetscIsNormalScalar(rho)),PETSC_COMM_SELF,PETSC_ERR_FP,"non-normal iterate rho=%g",(double)PetscRealPart(rho));
835     }
836     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Newton iteration for star.rho diverged after %" PetscInt_FMT " iterations",i);
837   }
838 converged:
839   if (L.u-c < 0 && 0 < star.u-c) { /* 1-wave is sonic rarefaction */
840     PetscScalar ufan[2];
841     ufan[0] = L.rho*PetscExpScalar(L.u/c - 1);
842     ufan[1] = c*ufan[0];
843     IsoGasFlux(c,ufan,flux);
844   } else if (star.u+c < 0 && 0 < R.u+c) { /* 2-wave is sonic rarefaction */
845     PetscScalar ufan[2];
846     ufan[0] = R.rho*PetscExpScalar(-R.u/c - 1);
847     ufan[1] = -c*ufan[0];
848     IsoGasFlux(c,ufan,flux);
849   } else if ((L.rho >= star.rho && L.u-c >= 0) || (L.rho < star.rho && (star.rho*star.u-L.rho*L.u)/(star.rho-L.rho) > 0)) {
850     /* 1-wave is supersonic rarefaction, or supersonic shock */
851     IsoGasFlux(c,uL,flux);
852   } else if ((star.rho <= R.rho && R.u+c <= 0) || (star.rho > R.rho && (R.rho*R.u-star.rho*star.u)/(R.rho-star.rho) < 0)) {
853     /* 2-wave is supersonic rarefaction or supersonic shock */
854     IsoGasFlux(c,uR,flux);
855   } else {
856     ustar[0] = star.rho;
857     ustar[1] = star.rho*star.u;
858     IsoGasFlux(c,ustar,flux);
859   }
860   *maxspeed = MaxAbs(MaxAbs(star.u-c,star.u+c),MaxAbs(L.u-c,R.u+c));
861   PetscFunctionReturn(0);
862 }
863 
864 static PetscErrorCode PhysicsRiemann_IsoGas_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
865 {
866   IsoGasCtx                   *phys = (IsoGasCtx*)vctx;
867   PetscScalar                 c = phys->acoustic_speed,fL[2],fR[2],s;
868   struct {PetscScalar rho,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};
869 
870   PetscFunctionBeginUser;
871   PetscCheck((L.rho > 0 && R.rho > 0),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed density is negative");
872   IsoGasFlux(c,uL,fL);
873   IsoGasFlux(c,uR,fR);
874   s         = PetscMax(PetscAbs(L.u),PetscAbs(R.u))+c;
875   flux[0]   = 0.5*(fL[0] + fR[0]) + 0.5*s*(uL[0] - uR[0]);
876   flux[1]   = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
877   *maxspeed = s;
878   PetscFunctionReturn(0);
879 }
880 
881 static PetscErrorCode PhysicsCharacteristic_IsoGas(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
882 {
883   IsoGasCtx      *phys = (IsoGasCtx*)vctx;
884   PetscReal      c     = phys->acoustic_speed;
885 
886   PetscFunctionBeginUser;
887   speeds[0] = u[1]/u[0] - c;
888   speeds[1] = u[1]/u[0] + c;
889   X[0*2+0]  = 1;
890   X[0*2+1]  = speeds[0];
891   X[1*2+0]  = 1;
892   X[1*2+1]  = speeds[1];
893   PetscCall(PetscArraycpy(Xi,X,4));
894   PetscCall(PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL));
895   PetscFunctionReturn(0);
896 }
897 
898 static PetscErrorCode PhysicsCreate_IsoGas(FVCtx *ctx)
899 {
900   IsoGasCtx         *user;
901   PetscFunctionList rlist = 0,rclist = 0;
902   char              rname[256] = "exact",rcname[256] = "characteristic";
903 
904   PetscFunctionBeginUser;
905   PetscCall(PetscNew(&user));
906   ctx->physics.sample         = PhysicsSample_IsoGas;
907   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
908   ctx->physics.user           = user;
909   ctx->physics.dof            = 2;
910 
911   PetscCall(PetscStrallocpy("density",&ctx->physics.fieldname[0]));
912   PetscCall(PetscStrallocpy("momentum",&ctx->physics.fieldname[1]));
913 
914   user->acoustic_speed = 1;
915 
916   PetscCall(RiemannListAdd(&rlist,"exact",  PhysicsRiemann_IsoGas_Exact));
917   PetscCall(RiemannListAdd(&rlist,"roe",    PhysicsRiemann_IsoGas_Roe));
918   PetscCall(RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_IsoGas_Rusanov));
919   PetscCall(ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_IsoGas));
920   PetscCall(ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative));
921   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for IsoGas","");
922     PetscCall(PetscOptionsReal("-physics_isogas_acoustic_speed","Acoustic speed","",user->acoustic_speed,&user->acoustic_speed,NULL));
923     PetscCall(PetscOptionsFList("-physics_isogas_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL));
924     PetscCall(PetscOptionsFList("-physics_isogas_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL));
925   PetscOptionsEnd();
926   PetscCall(RiemannListFind(rlist,rname,&ctx->physics.riemann));
927   PetscCall(ReconstructListFind(rclist,rcname,&ctx->physics.characteristic));
928   PetscCall(PetscFunctionListDestroy(&rlist));
929   PetscCall(PetscFunctionListDestroy(&rclist));
930   PetscFunctionReturn(0);
931 }
932 
933 /* --------------------------------- Shallow Water ----------------------------------- */
934 typedef struct {
935   PetscReal gravity;
936 } ShallowCtx;
937 
938 static inline void ShallowFlux(ShallowCtx *phys,const PetscScalar *u,PetscScalar *f)
939 {
940   f[0] = u[1];
941   f[1] = PetscSqr(u[1])/u[0] + 0.5*phys->gravity*PetscSqr(u[0]);
942 }
943 
944 static PetscErrorCode PhysicsRiemann_Shallow_Exact(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
945 {
946   ShallowCtx                *phys = (ShallowCtx*)vctx;
947   PetscScalar               g    = phys->gravity,ustar[2],cL,cR,c,cstar;
948   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]},star;
949   PetscInt                  i;
950 
951   PetscFunctionBeginUser;
952   PetscCheck((L.h > 0 && R.h > 0),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative");
953   cL = PetscSqrtScalar(g*L.h);
954   cR = PetscSqrtScalar(g*R.h);
955   c  = PetscMax(cL,cR);
956   {
957     /* Solve for star state */
958     const PetscInt maxits = 50;
959     PetscScalar tmp,res,res0=0,h0,h = 0.5*(L.h + R.h); /* initial guess */
960     h0 = h;
961     for (i=0; i<maxits; i++) {
962       PetscScalar fr,fl,dfr,dfl;
963       fl = (L.h < h)
964         ? PetscSqrtScalar(0.5*g*(h*h - L.h*L.h)*(1/L.h - 1/h)) /* shock */
965         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*L.h);   /* rarefaction */
966       fr = (R.h < h)
967         ? PetscSqrtScalar(0.5*g*(h*h - R.h*R.h)*(1/R.h - 1/h)) /* shock */
968         : 2*PetscSqrtScalar(g*h) - 2*PetscSqrtScalar(g*R.h);   /* rarefaction */
969       res = R.u - L.u + fr + fl;
970       PetscCheck(!PetscIsInfOrNanScalar(res),PETSC_COMM_SELF,PETSC_ERR_FP,"Infinity or Not-a-Number generated in computation");
971       if (PetscAbsScalar(res) < 1e-8 || (i > 0 && PetscAbsScalar(h-h0) < 1e-8)) {
972         star.h = h;
973         star.u = L.u - fl;
974         goto converged;
975       } else if (i > 0 && PetscAbsScalar(res) >= PetscAbsScalar(res0)) {        /* Line search */
976         h = 0.8*h0 + 0.2*h;
977         continue;
978       }
979       /* Accept the last step and take another */
980       res0 = res;
981       h0 = h;
982       dfl = (L.h < h) ? 0.5/fl*0.5*g*(-L.h*L.h/(h*h) - 1 + 2*h/L.h) : PetscSqrtScalar(g/h);
983       dfr = (R.h < h) ? 0.5/fr*0.5*g*(-R.h*R.h/(h*h) - 1 + 2*h/R.h) : PetscSqrtScalar(g/h);
984       tmp = h - res/(dfr+dfl);
985       if (tmp <= 0) h /= 2;   /* Guard against Newton shooting off to a negative thickness */
986       else h = tmp;
987       PetscCheck(((h > 0) && PetscIsNormalScalar(h)),PETSC_COMM_SELF,PETSC_ERR_FP,"non-normal iterate h=%g",(double)h);
988     }
989     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Newton iteration for star.h diverged after %" PetscInt_FMT " iterations",i);
990   }
991 converged:
992   cstar = PetscSqrtScalar(g*star.h);
993   if (L.u-cL < 0 && 0 < star.u-cstar) { /* 1-wave is sonic rarefaction */
994     PetscScalar ufan[2];
995     ufan[0] = 1/g*PetscSqr(L.u/3 + 2./3*cL);
996     ufan[1] = PetscSqrtScalar(g*ufan[0])*ufan[0];
997     ShallowFlux(phys,ufan,flux);
998   } else if (star.u+cstar < 0 && 0 < R.u+cR) { /* 2-wave is sonic rarefaction */
999     PetscScalar ufan[2];
1000     ufan[0] = 1/g*PetscSqr(R.u/3 - 2./3*cR);
1001     ufan[1] = -PetscSqrtScalar(g*ufan[0])*ufan[0];
1002     ShallowFlux(phys,ufan,flux);
1003   } else if ((L.h >= star.h && L.u-c >= 0) || (L.h<star.h && (star.h*star.u-L.h*L.u)/(star.h-L.h) > 0)) {
1004     /* 1-wave is right-travelling shock (supersonic) */
1005     ShallowFlux(phys,uL,flux);
1006   } else if ((star.h <= R.h && R.u+c <= 0) || (star.h>R.h && (R.h*R.u-star.h*star.h)/(R.h-star.h) < 0)) {
1007     /* 2-wave is left-travelling shock (supersonic) */
1008     ShallowFlux(phys,uR,flux);
1009   } else {
1010     ustar[0] = star.h;
1011     ustar[1] = star.h*star.u;
1012     ShallowFlux(phys,ustar,flux);
1013   }
1014   *maxspeed = MaxAbs(MaxAbs(star.u-cstar,star.u+cstar),MaxAbs(L.u-cL,R.u+cR));
1015   PetscFunctionReturn(0);
1016 }
1017 
1018 static PetscErrorCode PhysicsRiemann_Shallow_Rusanov(void *vctx,PetscInt m,const PetscScalar *uL,const PetscScalar *uR,PetscScalar *flux,PetscReal *maxspeed)
1019 {
1020   ShallowCtx                *phys = (ShallowCtx*)vctx;
1021   PetscScalar               g = phys->gravity,fL[2],fR[2],s;
1022   struct {PetscScalar h,u;} L = {uL[0],uL[1]/uL[0]},R = {uR[0],uR[1]/uR[0]};
1023 
1024   PetscFunctionBeginUser;
1025   PetscCheck((L.h > 0 && R.h > 0),PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reconstructed thickness is negative");
1026   ShallowFlux(phys,uL,fL);
1027   ShallowFlux(phys,uR,fR);
1028   s         = PetscMax(PetscAbs(L.u)+PetscSqrtScalar(g*L.h),PetscAbs(R.u)+PetscSqrtScalar(g*R.h));
1029   flux[0]   = 0.5*(fL[0] + fR[0]) + 0.5*s*(uL[0] - uR[0]);
1030   flux[1]   = 0.5*(fL[1] + fR[1]) + 0.5*s*(uL[1] - uR[1]);
1031   *maxspeed = s;
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 static PetscErrorCode PhysicsCharacteristic_Shallow(void *vctx,PetscInt m,const PetscScalar *u,PetscScalar *X,PetscScalar *Xi,PetscReal *speeds)
1036 {
1037   ShallowCtx     *phys = (ShallowCtx*)vctx;
1038   PetscReal      c;
1039 
1040   PetscFunctionBeginUser;
1041   c         = PetscSqrtScalar(u[0]*phys->gravity);
1042   speeds[0] = u[1]/u[0] - c;
1043   speeds[1] = u[1]/u[0] + c;
1044   X[0*2+0]  = 1;
1045   X[0*2+1]  = speeds[0];
1046   X[1*2+0]  = 1;
1047   X[1*2+1]  = speeds[1];
1048   PetscCall(PetscArraycpy(Xi,X,4));
1049   PetscCall(PetscKernel_A_gets_inverse_A_2(Xi,0,PETSC_FALSE,NULL));
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 static PetscErrorCode PhysicsCreate_Shallow(FVCtx *ctx)
1054 {
1055   ShallowCtx        *user;
1056   PetscFunctionList rlist = 0,rclist = 0;
1057   char              rname[256] = "exact",rcname[256] = "characteristic";
1058 
1059   PetscFunctionBeginUser;
1060   PetscCall(PetscNew(&user));
1061   /* Shallow water and Isothermal Gas dynamics are similar so we reuse initial conditions for now */
1062   ctx->physics.sample         = PhysicsSample_IsoGas;
1063   ctx->physics.destroy        = PhysicsDestroy_SimpleFree;
1064   ctx->physics.user           = user;
1065   ctx->physics.dof            = 2;
1066 
1067   PetscCall(PetscStrallocpy("density",&ctx->physics.fieldname[0]));
1068   PetscCall(PetscStrallocpy("momentum",&ctx->physics.fieldname[1]));
1069 
1070   user->gravity = 1;
1071 
1072   PetscCall(RiemannListAdd(&rlist,"exact",  PhysicsRiemann_Shallow_Exact));
1073   PetscCall(RiemannListAdd(&rlist,"rusanov",PhysicsRiemann_Shallow_Rusanov));
1074   PetscCall(ReconstructListAdd(&rclist,"characteristic",PhysicsCharacteristic_Shallow));
1075   PetscCall(ReconstructListAdd(&rclist,"conservative",PhysicsCharacteristic_Conservative));
1076   PetscOptionsBegin(ctx->comm,ctx->prefix,"Options for Shallow","");
1077     PetscCall(PetscOptionsReal("-physics_shallow_gravity","Gravity","",user->gravity,&user->gravity,NULL));
1078     PetscCall(PetscOptionsFList("-physics_shallow_riemann","Riemann solver","",rlist,rname,rname,sizeof(rname),NULL));
1079     PetscCall(PetscOptionsFList("-physics_shallow_reconstruct","Reconstruction","",rclist,rcname,rcname,sizeof(rcname),NULL));
1080   PetscOptionsEnd();
1081   PetscCall(RiemannListFind(rlist,rname,&ctx->physics.riemann));
1082   PetscCall(ReconstructListFind(rclist,rcname,&ctx->physics.characteristic));
1083   PetscCall(PetscFunctionListDestroy(&rlist));
1084   PetscCall(PetscFunctionListDestroy(&rclist));
1085   PetscFunctionReturn(0);
1086 }
1087 
1088 /* --------------------------------- Finite Volume Solver ----------------------------------- */
1089 
1090 static PetscErrorCode FVRHSFunction(TS ts,PetscReal time,Vec X,Vec F,void *vctx)
1091 {
1092   FVCtx             *ctx = (FVCtx*)vctx;
1093   PetscInt          i,j,k,Mx,dof,xs,xm;
1094   PetscReal         hx,cfl_idt = 0;
1095   PetscScalar       *x,*f,*slope;
1096   Vec               Xloc;
1097   DM                da;
1098 
1099   PetscFunctionBeginUser;
1100   PetscCall(TSGetDM(ts,&da));
1101   PetscCall(DMGetLocalVector(da,&Xloc));
1102   PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
1103   hx   = (ctx->xmax - ctx->xmin)/Mx;
1104   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
1105   PetscCall(DMGlobalToLocalEnd  (da,X,INSERT_VALUES,Xloc));
1106 
1107   PetscCall(VecZeroEntries(F));
1108 
1109   PetscCall(DMDAVecGetArray(da,Xloc,&x));
1110   PetscCall(DMDAVecGetArray(da,F,&f));
1111   PetscCall(DMDAGetArray(da,PETSC_TRUE,&slope));
1112 
1113   PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1114 
1115   if (ctx->bctype == FVBC_OUTFLOW) {
1116     for (i=xs-2; i<0; i++) {
1117       for (j=0; j<dof; j++) x[i*dof+j] = x[j];
1118     }
1119     for (i=Mx; i<xs+xm+2; i++) {
1120       for (j=0; j<dof; j++) x[i*dof+j] = x[(xs+xm-1)*dof+j];
1121     }
1122   }
1123   for (i=xs-1; i<xs+xm+1; i++) {
1124     struct _LimitInfo info;
1125     PetscScalar       *cjmpL,*cjmpR;
1126     /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */
1127     PetscCall((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
1128     /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */
1129     PetscCall(PetscArrayzero(ctx->cjmpLR,2*dof));
1130     cjmpL = &ctx->cjmpLR[0];
1131     cjmpR = &ctx->cjmpLR[dof];
1132     for (j=0; j<dof; j++) {
1133       PetscScalar jmpL,jmpR;
1134       jmpL = x[(i+0)*dof+j] - x[(i-1)*dof+j];
1135       jmpR = x[(i+1)*dof+j] - x[(i+0)*dof+j];
1136       for (k=0; k<dof; k++) {
1137         cjmpL[k] += ctx->Rinv[k+j*dof] * jmpL;
1138         cjmpR[k] += ctx->Rinv[k+j*dof] * jmpR;
1139       }
1140     }
1141     /* Apply limiter to the left and right characteristic jumps */
1142     info.m  = dof;
1143     info.hx = hx;
1144     (*ctx->limit)(&info,cjmpL,cjmpR,ctx->cslope);
1145     for (j=0; j<dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */
1146     for (j=0; j<dof; j++) {
1147       PetscScalar tmp = 0;
1148       for (k=0; k<dof; k++) tmp += ctx->R[j+k*dof] * ctx->cslope[k];
1149       slope[i*dof+j] = tmp;
1150     }
1151   }
1152 
1153   for (i=xs; i<xs+xm+1; i++) {
1154     PetscReal   maxspeed;
1155     PetscScalar *uL,*uR;
1156     uL = &ctx->uLR[0];
1157     uR = &ctx->uLR[dof];
1158     for (j=0; j<dof; j++) {
1159       uL[j] = x[(i-1)*dof+j] + slope[(i-1)*dof+j]*hx/2;
1160       uR[j] = x[(i-0)*dof+j] - slope[(i-0)*dof+j]*hx/2;
1161     }
1162     PetscCall((*ctx->physics.riemann)(ctx->physics.user,dof,uL,uR,ctx->flux,&maxspeed));
1163     cfl_idt = PetscMax(cfl_idt,PetscAbsScalar(maxspeed/hx)); /* Max allowable value of 1/Delta t */
1164 
1165     if (i > xs) {
1166       for (j=0; j<dof; j++) f[(i-1)*dof+j] -= ctx->flux[j]/hx;
1167     }
1168     if (i < xs+xm) {
1169       for (j=0; j<dof; j++) f[i*dof+j] += ctx->flux[j]/hx;
1170     }
1171   }
1172 
1173   PetscCall(DMDAVecRestoreArray(da,Xloc,&x));
1174   PetscCall(DMDAVecRestoreArray(da,F,&f));
1175   PetscCall(DMDARestoreArray(da,PETSC_TRUE,&slope));
1176   PetscCall(DMRestoreLocalVector(da,&Xloc));
1177 
1178   PetscCallMPI(MPI_Allreduce(&cfl_idt,&ctx->cfl_idt,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)da)));
1179   if (0) {
1180     /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */
1181     PetscReal dt,tnow;
1182     PetscCall(TSGetTimeStep(ts,&dt));
1183     PetscCall(TSGetTime(ts,&tnow));
1184     if (dt > 0.5/ctx->cfl_idt) {
1185       PetscCall(PetscPrintf(ctx->comm,"Stability constraint exceeded at t=%g, dt %g > %g\n",(double)tnow,(double)dt,(double)(0.5/ctx->cfl_idt)));
1186     }
1187   }
1188   PetscFunctionReturn(0);
1189 }
1190 
1191 static PetscErrorCode SmallMatMultADB(PetscScalar *C,PetscInt bs,const PetscScalar *A,const PetscReal *D,const PetscScalar *B)
1192 {
1193   PetscInt i,j,k;
1194 
1195   PetscFunctionBeginUser;
1196   for (i=0; i<bs; i++) {
1197     for (j=0; j<bs; j++) {
1198       PetscScalar tmp = 0;
1199       for (k=0; k<bs; k++) tmp += A[i*bs+k] * D[k] * B[k*bs+j];
1200       C[i*bs+j] = tmp;
1201     }
1202   }
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 static PetscErrorCode FVIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat A,Mat B,void *vctx)
1207 {
1208   FVCtx             *ctx = (FVCtx*)vctx;
1209   PetscInt          i,j,dof = ctx->physics.dof;
1210   PetscScalar       *J;
1211   const PetscScalar *x;
1212   PetscReal         hx;
1213   DM                da;
1214   DMDALocalInfo     dainfo;
1215 
1216   PetscFunctionBeginUser;
1217   PetscCall(TSGetDM(ts,&da));
1218   PetscCall(DMDAVecGetArrayRead(da,X,(void*)&x));
1219   PetscCall(DMDAGetLocalInfo(da,&dainfo));
1220   hx   = (ctx->xmax - ctx->xmin)/dainfo.mx;
1221   PetscCall(PetscMalloc1(dof*dof,&J));
1222   for (i=dainfo.xs; i<dainfo.xs+dainfo.xm; i++) {
1223     PetscCall((*ctx->physics.characteristic)(ctx->physics.user,dof,&x[i*dof],ctx->R,ctx->Rinv,ctx->speeds));
1224     for (j=0; j<dof; j++) ctx->speeds[j] = PetscAbs(ctx->speeds[j]);
1225     PetscCall(SmallMatMultADB(J,dof,ctx->R,ctx->speeds,ctx->Rinv));
1226     for (j=0; j<dof*dof; j++) J[j] = J[j]/hx + shift*(j/dof == j%dof);
1227     PetscCall(MatSetValuesBlocked(B,1,&i,1,&i,J,INSERT_VALUES));
1228   }
1229   PetscCall(PetscFree(J));
1230   PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&x));
1231 
1232   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1233   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1234   if (A != B) {
1235     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1236     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1237   }
1238   PetscFunctionReturn(0);
1239 }
1240 
1241 static PetscErrorCode FVSample(FVCtx *ctx,DM da,PetscReal time,Vec U)
1242 {
1243   PetscScalar    *u,*uj;
1244   PetscInt       i,j,k,dof,xs,xm,Mx;
1245 
1246   PetscFunctionBeginUser;
1247   PetscCheck(ctx->physics.sample,PETSC_COMM_SELF,PETSC_ERR_SUP,"Physics has not provided a sampling function");
1248   PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
1249   PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1250   PetscCall(DMDAVecGetArray(da,U,&u));
1251   PetscCall(PetscMalloc1(dof,&uj));
1252   for (i=xs; i<xs+xm; i++) {
1253     const PetscReal h = (ctx->xmax-ctx->xmin)/Mx,xi = ctx->xmin+h/2+i*h;
1254     const PetscInt  N = 200;
1255     /* Integrate over cell i using trapezoid rule with N points. */
1256     for (k=0; k<dof; k++) u[i*dof+k] = 0;
1257     for (j=0; j<N+1; j++) {
1258       PetscScalar xj = xi+h*(j-N/2)/(PetscReal)N;
1259       PetscCall((*ctx->physics.sample)(ctx->physics.user,ctx->initial,ctx->bctype,ctx->xmin,ctx->xmax,time,xj,uj));
1260       for (k=0; k<dof; k++) u[i*dof+k] += ((j==0 || j==N) ? 0.5 : 1.0)*uj[k]/N;
1261     }
1262   }
1263   PetscCall(DMDAVecRestoreArray(da,U,&u));
1264   PetscCall(PetscFree(uj));
1265   PetscFunctionReturn(0);
1266 }
1267 
1268 static PetscErrorCode SolutionStatsView(DM da,Vec X,PetscViewer viewer)
1269 {
1270   PetscReal         xmin,xmax;
1271   PetscScalar       sum,tvsum,tvgsum;
1272   const PetscScalar *x;
1273   PetscInt          imin,imax,Mx,i,j,xs,xm,dof;
1274   Vec               Xloc;
1275   PetscBool         iascii;
1276 
1277   PetscFunctionBeginUser;
1278   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1279   if (iascii) {
1280     /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */
1281     PetscCall(DMGetLocalVector(da,&Xloc));
1282     PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
1283     PetscCall(DMGlobalToLocalEnd  (da,X,INSERT_VALUES,Xloc));
1284     PetscCall(DMDAVecGetArrayRead(da,Xloc,(void*)&x));
1285     PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1286     PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
1287     tvsum = 0;
1288     for (i=xs; i<xs+xm; i++) {
1289       for (j=0; j<dof; j++) tvsum += PetscAbsScalar(x[i*dof+j] - x[(i-1)*dof+j]);
1290     }
1291     PetscCallMPI(MPI_Allreduce(&tvsum,&tvgsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)da)));
1292     PetscCall(DMDAVecRestoreArrayRead(da,Xloc,(void*)&x));
1293     PetscCall(DMRestoreLocalVector(da,&Xloc));
1294 
1295     PetscCall(VecMin(X,&imin,&xmin));
1296     PetscCall(VecMax(X,&imax,&xmax));
1297     PetscCall(VecSum(X,&sum));
1298     PetscCall(PetscViewerASCIIPrintf(viewer,"Solution range [%8.5f,%8.5f] with extrema at %" PetscInt_FMT " and %" PetscInt_FMT ", mean %8.5f, ||x||_TV %8.5f\n",(double)xmin,(double)xmax,imin,imax,(double)(sum/Mx),(double)(tvgsum/Mx)));
1299   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type not supported");
1300   PetscFunctionReturn(0);
1301 }
1302 
1303 static PetscErrorCode SolutionErrorNorms(FVCtx *ctx,DM da,PetscReal t,Vec X,PetscReal *nrm1,PetscReal *nrmsup)
1304 {
1305   Vec            Y;
1306   PetscInt       Mx;
1307 
1308   PetscFunctionBeginUser;
1309   PetscCall(VecGetSize(X,&Mx));
1310   PetscCall(VecDuplicate(X,&Y));
1311   PetscCall(FVSample(ctx,da,t,Y));
1312   PetscCall(VecAYPX(Y,-1,X));
1313   PetscCall(VecNorm(Y,NORM_1,nrm1));
1314   PetscCall(VecNorm(Y,NORM_INFINITY,nrmsup));
1315   *nrm1 /= Mx;
1316   PetscCall(VecDestroy(&Y));
1317   PetscFunctionReturn(0);
1318 }
1319 
1320 int main(int argc,char *argv[])
1321 {
1322   char              lname[256] = "mc",physname[256] = "advect",final_fname[256] = "solution.m";
1323   PetscFunctionList limiters   = 0,physics = 0;
1324   MPI_Comm          comm;
1325   TS                ts;
1326   DM                da;
1327   Vec               X,X0,R;
1328   Mat               B;
1329   FVCtx             ctx;
1330   PetscInt          i,dof,xs,xm,Mx,draw = 0;
1331   PetscBool         view_final = PETSC_FALSE;
1332   PetscReal         ptime;
1333 
1334   PetscCall(PetscInitialize(&argc,&argv,0,help));
1335   comm = PETSC_COMM_WORLD;
1336   PetscCall(PetscMemzero(&ctx,sizeof(ctx)));
1337 
1338   /* Register limiters to be available on the command line */
1339   PetscCall(PetscFunctionListAdd(&limiters,"upwind"              ,Limit_Upwind));
1340   PetscCall(PetscFunctionListAdd(&limiters,"lax-wendroff"        ,Limit_LaxWendroff));
1341   PetscCall(PetscFunctionListAdd(&limiters,"beam-warming"        ,Limit_BeamWarming));
1342   PetscCall(PetscFunctionListAdd(&limiters,"fromm"               ,Limit_Fromm));
1343   PetscCall(PetscFunctionListAdd(&limiters,"minmod"              ,Limit_Minmod));
1344   PetscCall(PetscFunctionListAdd(&limiters,"superbee"            ,Limit_Superbee));
1345   PetscCall(PetscFunctionListAdd(&limiters,"mc"                  ,Limit_MC));
1346   PetscCall(PetscFunctionListAdd(&limiters,"vanleer"             ,Limit_VanLeer));
1347   PetscCall(PetscFunctionListAdd(&limiters,"vanalbada"           ,Limit_VanAlbada));
1348   PetscCall(PetscFunctionListAdd(&limiters,"vanalbadatvd"        ,Limit_VanAlbadaTVD));
1349   PetscCall(PetscFunctionListAdd(&limiters,"koren"               ,Limit_Koren));
1350   PetscCall(PetscFunctionListAdd(&limiters,"korensym"            ,Limit_KorenSym));
1351   PetscCall(PetscFunctionListAdd(&limiters,"koren3"              ,Limit_Koren3));
1352   PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon2"     ,Limit_CadaTorrilhon2));
1353   PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r0p1",Limit_CadaTorrilhon3R0p1));
1354   PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r1"  ,Limit_CadaTorrilhon3R1));
1355   PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r10" ,Limit_CadaTorrilhon3R10));
1356   PetscCall(PetscFunctionListAdd(&limiters,"cada-torrilhon3-r100",Limit_CadaTorrilhon3R100));
1357 
1358   /* Register physical models to be available on the command line */
1359   PetscCall(PetscFunctionListAdd(&physics,"advect"          ,PhysicsCreate_Advect));
1360   PetscCall(PetscFunctionListAdd(&physics,"burgers"         ,PhysicsCreate_Burgers));
1361   PetscCall(PetscFunctionListAdd(&physics,"traffic"         ,PhysicsCreate_Traffic));
1362   PetscCall(PetscFunctionListAdd(&physics,"acoustics"       ,PhysicsCreate_Acoustics));
1363   PetscCall(PetscFunctionListAdd(&physics,"isogas"          ,PhysicsCreate_IsoGas));
1364   PetscCall(PetscFunctionListAdd(&physics,"shallow"         ,PhysicsCreate_Shallow));
1365 
1366   ctx.comm = comm;
1367   ctx.cfl  = 0.9; ctx.bctype = FVBC_PERIODIC;
1368   ctx.xmin = -1; ctx.xmax = 1;
1369   PetscOptionsBegin(comm,NULL,"Finite Volume solver options","");
1370     PetscCall(PetscOptionsReal("-xmin","X min","",ctx.xmin,&ctx.xmin,NULL));
1371     PetscCall(PetscOptionsReal("-xmax","X max","",ctx.xmax,&ctx.xmax,NULL));
1372     PetscCall(PetscOptionsFList("-limit","Name of flux limiter to use","",limiters,lname,lname,sizeof(lname),NULL));
1373     PetscCall(PetscOptionsFList("-physics","Name of physics (Riemann solver and characteristics) to use","",physics,physname,physname,sizeof(physname),NULL));
1374     PetscCall(PetscOptionsInt("-draw","Draw solution vector, bitwise OR of (1=initial,2=final,4=final error)","",draw,&draw,NULL));
1375     PetscCall(PetscOptionsString("-view_final","Write final solution in ASCII MATLAB format to given file name","",final_fname,final_fname,sizeof(final_fname),&view_final));
1376     PetscCall(PetscOptionsInt("-initial","Initial condition (depends on the physics)","",ctx.initial,&ctx.initial,NULL));
1377     PetscCall(PetscOptionsBool("-exact","Compare errors with exact solution","",ctx.exact,&ctx.exact,NULL));
1378     PetscCall(PetscOptionsReal("-cfl","CFL number to time step at","",ctx.cfl,&ctx.cfl,NULL));
1379     PetscCall(PetscOptionsEnum("-bc_type","Boundary condition","",FVBCTypes,(PetscEnum)ctx.bctype,(PetscEnum*)&ctx.bctype,NULL));
1380   PetscOptionsEnd();
1381 
1382   /* Choose the limiter from the list of registered limiters */
1383   PetscCall(PetscFunctionListFind(limiters,lname,&ctx.limit));
1384   PetscCheck(ctx.limit,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Limiter '%s' not found",lname);
1385 
1386   /* Choose the physics from the list of registered models */
1387   {
1388     PetscErrorCode (*r)(FVCtx*);
1389     PetscCall(PetscFunctionListFind(physics,physname,&r));
1390     PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Physics '%s' not found",physname);
1391     /* Create the physics, will set the number of fields and their names */
1392     PetscCall((*r)(&ctx));
1393   }
1394 
1395   /* Create a DMDA to manage the parallel grid */
1396   PetscCall(DMDACreate1d(comm,DM_BOUNDARY_PERIODIC,50,ctx.physics.dof,2,NULL,&da));
1397   PetscCall(DMSetFromOptions(da));
1398   PetscCall(DMSetUp(da));
1399   /* Inform the DMDA of the field names provided by the physics. */
1400   /* The names will be shown in the title bars when run with -ts_monitor_draw_solution */
1401   for (i=0; i<ctx.physics.dof; i++) {
1402     PetscCall(DMDASetFieldName(da,i,ctx.physics.fieldname[i]));
1403   }
1404   PetscCall(DMDAGetInfo(da,0, &Mx,0,0, 0,0,0, &dof,0,0,0,0,0));
1405   PetscCall(DMDAGetCorners(da,&xs,0,0,&xm,0,0));
1406 
1407   /* Set coordinates of cell centers */
1408   PetscCall(DMDASetUniformCoordinates(da,ctx.xmin+0.5*(ctx.xmax-ctx.xmin)/Mx,ctx.xmax+0.5*(ctx.xmax-ctx.xmin)/Mx,0,0,0,0));
1409 
1410   /* Allocate work space for the Finite Volume solver (so it doesn't have to be reallocated on each function evaluation) */
1411   PetscCall(PetscMalloc4(dof*dof,&ctx.R,dof*dof,&ctx.Rinv,2*dof,&ctx.cjmpLR,1*dof,&ctx.cslope));
1412   PetscCall(PetscMalloc3(2*dof,&ctx.uLR,dof,&ctx.flux,dof,&ctx.speeds));
1413 
1414   /* Create a vector to store the solution and to save the initial state */
1415   PetscCall(DMCreateGlobalVector(da,&X));
1416   PetscCall(VecDuplicate(X,&X0));
1417   PetscCall(VecDuplicate(X,&R));
1418 
1419   PetscCall(DMCreateMatrix(da,&B));
1420 
1421   /* Create a time-stepping object */
1422   PetscCall(TSCreate(comm,&ts));
1423   PetscCall(TSSetDM(ts,da));
1424   PetscCall(TSSetRHSFunction(ts,R,FVRHSFunction,&ctx));
1425   PetscCall(TSSetIJacobian(ts,B,B,FVIJacobian,&ctx));
1426   PetscCall(TSSetType(ts,TSSSP));
1427   PetscCall(TSSetMaxTime(ts,10));
1428   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1429 
1430   /* Compute initial conditions and starting time step */
1431   PetscCall(FVSample(&ctx,da,0,X0));
1432   PetscCall(FVRHSFunction(ts,0,X0,X,(void*)&ctx)); /* Initial function evaluation, only used to determine max speed */
1433   PetscCall(VecCopy(X0,X));                        /* The function value was not used so we set X=X0 again */
1434   PetscCall(TSSetTimeStep(ts,ctx.cfl/ctx.cfl_idt));
1435   PetscCall(TSSetFromOptions(ts)); /* Take runtime options */
1436   PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
1437   {
1438     PetscReal nrm1,nrmsup;
1439     PetscInt  steps;
1440 
1441     PetscCall(TSSolve(ts,X));
1442     PetscCall(TSGetSolveTime(ts,&ptime));
1443     PetscCall(TSGetStepNumber(ts,&steps));
1444 
1445     PetscCall(PetscPrintf(comm,"Final time %8.5f, steps %" PetscInt_FMT "\n",(double)ptime,steps));
1446     if (ctx.exact) {
1447       PetscCall(SolutionErrorNorms(&ctx,da,ptime,X,&nrm1,&nrmsup));
1448       PetscCall(PetscPrintf(comm,"Error ||x-x_e||_1 %8.4e  ||x-x_e||_sup %8.4e\n",(double)nrm1,(double)nrmsup));
1449     }
1450   }
1451 
1452   PetscCall(SolutionStatsView(da,X,PETSC_VIEWER_STDOUT_WORLD));
1453   if (draw & 0x1) PetscCall(VecView(X0,PETSC_VIEWER_DRAW_WORLD));
1454   if (draw & 0x2) PetscCall(VecView(X,PETSC_VIEWER_DRAW_WORLD));
1455   if (draw & 0x4) {
1456     Vec Y;
1457     PetscCall(VecDuplicate(X,&Y));
1458     PetscCall(FVSample(&ctx,da,ptime,Y));
1459     PetscCall(VecAYPX(Y,-1,X));
1460     PetscCall(VecView(Y,PETSC_VIEWER_DRAW_WORLD));
1461     PetscCall(VecDestroy(&Y));
1462   }
1463 
1464   if (view_final) {
1465     PetscViewer viewer;
1466     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,final_fname,&viewer));
1467     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
1468     PetscCall(VecView(X,viewer));
1469     PetscCall(PetscViewerPopFormat(viewer));
1470     PetscCall(PetscViewerDestroy(&viewer));
1471   }
1472 
1473   /* Clean up */
1474   PetscCall((*ctx.physics.destroy)(ctx.physics.user));
1475   for (i=0; i<ctx.physics.dof; i++) PetscCall(PetscFree(ctx.physics.fieldname[i]));
1476   PetscCall(PetscFree4(ctx.R,ctx.Rinv,ctx.cjmpLR,ctx.cslope));
1477   PetscCall(PetscFree3(ctx.uLR,ctx.flux,ctx.speeds));
1478   PetscCall(VecDestroy(&X));
1479   PetscCall(VecDestroy(&X0));
1480   PetscCall(VecDestroy(&R));
1481   PetscCall(MatDestroy(&B));
1482   PetscCall(DMDestroy(&da));
1483   PetscCall(TSDestroy(&ts));
1484   PetscCall(PetscFunctionListDestroy(&limiters));
1485   PetscCall(PetscFunctionListDestroy(&physics));
1486   PetscCall(PetscFinalize());
1487   return 0;
1488 }
1489 
1490 /*TEST
1491 
1492     build:
1493       requires: !complex
1494 
1495     test:
1496       args: -da_grid_x 100 -initial 1 -xmin -2 -xmax 5 -exact -limit mc
1497       requires: !complex !single
1498 
1499     test:
1500       suffix: 2
1501       args: -da_grid_x 100 -initial 2 -xmin -2 -xmax 2 -exact -limit mc -physics burgers -bc_type outflow -ts_max_time 1
1502       filter:  sed "s/at 48/at 0/g"
1503       requires: !complex !single
1504 
1505     test:
1506       suffix: 3
1507       args: -da_grid_x 100 -initial 2 -xmin -2 -xmax 2 -exact -limit mc -physics burgers -bc_type outflow -ts_max_time 1
1508       nsize: 3
1509       filter:  sed "s/at 48/at 0/g"
1510       requires: !complex !single
1511 
1512 TEST*/
1513