xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex6.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] ="Model Equations for Advection \n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*
5*c4762a1bSJed Brown     Modified from ex3.c
6*c4762a1bSJed Brown     Page 9, Section 1.2 Model Equations for Advection-Diffusion
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown           u_t + a u_x = 0, 0<= x <= 1.0
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown    The initial conditions used here different from the book.
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown    Example:
13*c4762a1bSJed Brown      ./ex6 -ts_monitor -ts_view_solution -ts_max_steps 100 -ts_monitor_solution draw -draw_pause .1
14*c4762a1bSJed Brown      ./ex6 -ts_monitor -ts_max_steps 100 -ts_monitor_lg_error -draw_pause .1
15*c4762a1bSJed Brown */
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown #include <petscts.h>
18*c4762a1bSJed Brown #include <petscdm.h>
19*c4762a1bSJed Brown #include <petscdmda.h>
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown /*
22*c4762a1bSJed Brown    User-defined application context - contains data needed by the
23*c4762a1bSJed Brown    application-provided call-back routines.
24*c4762a1bSJed Brown */
25*c4762a1bSJed Brown typedef struct {
26*c4762a1bSJed Brown   PetscReal a;   /* advection strength */
27*c4762a1bSJed Brown } AppCtx;
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown /* User-defined routines */
30*c4762a1bSJed Brown extern PetscErrorCode InitialConditions(TS,Vec,AppCtx*);
31*c4762a1bSJed Brown extern PetscErrorCode Solution(TS,PetscReal,Vec,AppCtx*);
32*c4762a1bSJed Brown extern PetscErrorCode IFunction_LaxFriedrichs(TS,PetscReal,Vec,Vec,Vec,void*);
33*c4762a1bSJed Brown extern PetscErrorCode IFunction_LaxWendroff(TS,PetscReal,Vec,Vec,Vec,void*);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown int main(int argc,char **argv)
36*c4762a1bSJed Brown {
37*c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
38*c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
39*c4762a1bSJed Brown   Vec            U;                      /* approximate solution vector */
40*c4762a1bSJed Brown   PetscErrorCode ierr;
41*c4762a1bSJed Brown   PetscReal      dt;
42*c4762a1bSJed Brown   DM             da;
43*c4762a1bSJed Brown   PetscInt       M;
44*c4762a1bSJed Brown   PetscMPIInt    rank;
45*c4762a1bSJed Brown   PetscBool      useLaxWendroff = PETSC_TRUE;
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   /* Initialize program and set problem parameters */
48*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
49*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   appctx.a  = -1.0;
52*c4762a1bSJed Brown   ierr      = PetscOptionsGetReal(NULL,NULL,"-a",&appctx.a,NULL);CHKERRQ(ierr);
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC, 60, 1, 1,NULL,&da);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   /* Create vector data structures for approximate and exact solutions */
59*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&U);CHKERRQ(ierr);
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown   /* Create timestepping solver context */
62*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   /* Function evaluation */
66*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-useLaxWendroff",&useLaxWendroff,NULL);CHKERRQ(ierr);
67*c4762a1bSJed Brown   if (useLaxWendroff) {
68*c4762a1bSJed Brown     if (!rank) {
69*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"... Use Lax-Wendroff finite volume\n");CHKERRQ(ierr);
70*c4762a1bSJed Brown     }
71*c4762a1bSJed Brown     ierr = TSSetIFunction(ts,NULL,IFunction_LaxWendroff,&appctx);CHKERRQ(ierr);
72*c4762a1bSJed Brown   } else {
73*c4762a1bSJed Brown     if (!rank) {
74*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"... Use Lax-LaxFriedrichs finite difference\n");CHKERRQ(ierr);
75*c4762a1bSJed Brown     }
76*c4762a1bSJed Brown     ierr = TSSetIFunction(ts,NULL,IFunction_LaxFriedrichs,&appctx);CHKERRQ(ierr);
77*c4762a1bSJed Brown   }
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   /* Customize timestepping solver */
80*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
81*c4762a1bSJed Brown   dt = 1.0/(PetscAbsReal(appctx.a)*M);
82*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,100);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,100.0);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown   /* Evaluate initial conditions */
90*c4762a1bSJed Brown   ierr = InitialConditions(ts,U,&appctx);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   /* For testing accuracy of TS with already known solution, e.g., '-ts_monitor_lg_error' */
93*c4762a1bSJed Brown   ierr = TSSetSolutionFunction(ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void*))Solution,&appctx);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   /* Run the timestepping solver */
96*c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   /* Free work space */
99*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
102*c4762a1bSJed Brown 
103*c4762a1bSJed Brown   ierr = PetscFinalize();
104*c4762a1bSJed Brown   return ierr;
105*c4762a1bSJed Brown }
106*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
107*c4762a1bSJed Brown /*
108*c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown    Input Parameter:
111*c4762a1bSJed Brown    u - uninitialized solution vector (global)
112*c4762a1bSJed Brown    appctx - user-defined application context
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown    Output Parameter:
115*c4762a1bSJed Brown    u - vector with solution at initial time (global)
116*c4762a1bSJed Brown */
117*c4762a1bSJed Brown PetscErrorCode InitialConditions(TS ts,Vec U,AppCtx *appctx)
118*c4762a1bSJed Brown {
119*c4762a1bSJed Brown   PetscScalar    *u;
120*c4762a1bSJed Brown   PetscErrorCode ierr;
121*c4762a1bSJed Brown   PetscInt       i,mstart,mend,um,M;
122*c4762a1bSJed Brown   DM             da;
123*c4762a1bSJed Brown   PetscReal      h;
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
126*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&um,0,0);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
128*c4762a1bSJed Brown   h    = 1.0/M;
129*c4762a1bSJed Brown   mend = mstart + um;
130*c4762a1bSJed Brown   /*
131*c4762a1bSJed Brown     Get a pointer to vector data.
132*c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
133*c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
134*c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
135*c4762a1bSJed Brown       the array.
136*c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
137*c4762a1bSJed Brown       C version.  See the users manual for details.
138*c4762a1bSJed Brown   */
139*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown   /*
142*c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
143*c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
144*c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
145*c4762a1bSJed Brown   */
146*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) u[i] = PetscSinReal(PETSC_PI*i*6.*h) + 3.*PetscSinReal(PETSC_PI*i*2.*h);
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown   /* Restore vector */
149*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
150*c4762a1bSJed Brown   return 0;
151*c4762a1bSJed Brown }
152*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
153*c4762a1bSJed Brown /*
154*c4762a1bSJed Brown    Solution - Computes the exact solution at a given time
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown    Input Parameters:
157*c4762a1bSJed Brown    t - current time
158*c4762a1bSJed Brown    solution - vector in which exact solution will be computed
159*c4762a1bSJed Brown    appctx - user-defined application context
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown    Output Parameter:
162*c4762a1bSJed Brown    solution - vector with the newly computed exact solution
163*c4762a1bSJed Brown               u(x,t) = sin(6*PI*(x - a*t)) + 3 * sin(2*PI*(x - a*t))
164*c4762a1bSJed Brown */
165*c4762a1bSJed Brown PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *appctx)
166*c4762a1bSJed Brown {
167*c4762a1bSJed Brown   PetscScalar    *u;
168*c4762a1bSJed Brown   PetscReal      a=appctx->a,h,PI6,PI2;
169*c4762a1bSJed Brown   PetscErrorCode ierr;
170*c4762a1bSJed Brown   PetscInt       i,mstart,mend,um,M;
171*c4762a1bSJed Brown   DM             da;
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&um,0,0);CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
176*c4762a1bSJed Brown   h    = 1.0/M;
177*c4762a1bSJed Brown   mend = mstart + um;
178*c4762a1bSJed Brown 
179*c4762a1bSJed Brown   /* Get a pointer to vector data. */
180*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   /* u[i] = sin(6*PI*(x[i] - a*t)) + 3 * sin(2*PI*(x[i] - a*t)) */
183*c4762a1bSJed Brown   PI6 = PETSC_PI*6.;
184*c4762a1bSJed Brown   PI2 = PETSC_PI*2.;
185*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
186*c4762a1bSJed Brown     u[i] = PetscSinReal(PI6*(i*h - a*t)) + 3.*PetscSinReal(PI2*(i*h - a*t));
187*c4762a1bSJed Brown   }
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown   /* Restore vector */
190*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
191*c4762a1bSJed Brown   return 0;
192*c4762a1bSJed Brown }
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
195*c4762a1bSJed Brown /*
196*c4762a1bSJed Brown  Use Lax-Friedrichs method to evaluate F(u,t) = du/dt + a *  du/dx
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown  See https://en.wikipedia.org/wiki/Lax%E2%80%93Friedrichs_method
199*c4762a1bSJed Brown  */
200*c4762a1bSJed Brown PetscErrorCode IFunction_LaxFriedrichs(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void* ctx)
201*c4762a1bSJed Brown {
202*c4762a1bSJed Brown   PetscErrorCode ierr;
203*c4762a1bSJed Brown   AppCtx         *appctx=(AppCtx*)ctx;
204*c4762a1bSJed Brown   PetscInt       mstart,mend,M,i,um;
205*c4762a1bSJed Brown   DM             da;
206*c4762a1bSJed Brown   Vec            Uold,localUold;
207*c4762a1bSJed Brown   PetscScalar    *uarray,*f,*uoldarray,h,uave,c;
208*c4762a1bSJed Brown   PetscReal      dt;
209*c4762a1bSJed Brown 
210*c4762a1bSJed Brown   PetscFunctionBegin;
211*c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
212*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&Uold);CHKERRQ(ierr);
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
215*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&um,0,0);CHKERRQ(ierr);
217*c4762a1bSJed Brown   h    = 1.0/M;
218*c4762a1bSJed Brown   mend = mstart + um;
219*c4762a1bSJed Brown   /* printf(" mstart %d, um %d\n",mstart,um); */
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localUold);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,Uold,INSERT_VALUES,localUold);CHKERRQ(ierr);
223*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,Uold,INSERT_VALUES,localUold);CHKERRQ(ierr);
224*c4762a1bSJed Brown 
225*c4762a1bSJed Brown   /* Get pointers to vector data */
226*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,U,&uarray);CHKERRQ(ierr);
227*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localUold,&uoldarray);CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
229*c4762a1bSJed Brown 
230*c4762a1bSJed Brown   /* advection */
231*c4762a1bSJed Brown   c = appctx->a*dt/h; /* Courant-Friedrichs-Lewy number (CFL number) */
232*c4762a1bSJed Brown 
233*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
234*c4762a1bSJed Brown     uave = 0.5*(uoldarray[i-1] + uoldarray[i+1]);
235*c4762a1bSJed Brown     f[i] = uarray[i] - uave + c*0.5*(uoldarray[i+1] - uoldarray[i-1]);
236*c4762a1bSJed Brown   }
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown   /* Restore vectors */
239*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,U,&uarray);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localUold,&uoldarray);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
242*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localUold);CHKERRQ(ierr);
243*c4762a1bSJed Brown   PetscFunctionReturn(0);
244*c4762a1bSJed Brown }
245*c4762a1bSJed Brown 
246*c4762a1bSJed Brown /*
247*c4762a1bSJed Brown  Use Lax-Wendroff method to evaluate F(u,t) = du/dt + a *  du/dx
248*c4762a1bSJed Brown */
249*c4762a1bSJed Brown PetscErrorCode IFunction_LaxWendroff(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void* ctx)
250*c4762a1bSJed Brown {
251*c4762a1bSJed Brown   PetscErrorCode ierr;
252*c4762a1bSJed Brown   AppCtx         *appctx=(AppCtx*)ctx;
253*c4762a1bSJed Brown   PetscInt       mstart,mend,M,i,um;
254*c4762a1bSJed Brown   DM             da;
255*c4762a1bSJed Brown   Vec            Uold,localUold;
256*c4762a1bSJed Brown   PetscScalar    *uarray,*f,*uoldarray,h,RFlux,LFlux,lambda;
257*c4762a1bSJed Brown   PetscReal      dt,a;
258*c4762a1bSJed Brown 
259*c4762a1bSJed Brown   PetscFunctionBegin;
260*c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
261*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&Uold);CHKERRQ(ierr);
262*c4762a1bSJed Brown 
263*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
264*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
265*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&mstart,0,0,&um,0,0);CHKERRQ(ierr);
266*c4762a1bSJed Brown   h    = 1.0/M;
267*c4762a1bSJed Brown   mend = mstart + um;
268*c4762a1bSJed Brown   /* printf(" mstart %d, um %d\n",mstart,um); */
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localUold);CHKERRQ(ierr);
271*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,Uold,INSERT_VALUES,localUold);CHKERRQ(ierr);
272*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,Uold,INSERT_VALUES,localUold);CHKERRQ(ierr);
273*c4762a1bSJed Brown 
274*c4762a1bSJed Brown   /* Get pointers to vector data */
275*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,U,&uarray);CHKERRQ(ierr);
276*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localUold,&uoldarray);CHKERRQ(ierr);
277*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
278*c4762a1bSJed Brown 
279*c4762a1bSJed Brown   /* advection -- finite volume (appctx->a < 0 -- can be relaxed?) */
280*c4762a1bSJed Brown   lambda = dt/h;
281*c4762a1bSJed Brown   a = appctx->a;
282*c4762a1bSJed Brown 
283*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
284*c4762a1bSJed Brown     RFlux = 0.5 * a * (uoldarray[i+1] + uoldarray[i]) - a*a*0.5*lambda * (uoldarray[i+1] - uoldarray[i]);
285*c4762a1bSJed Brown     LFlux = 0.5 * a * (uoldarray[i-1] + uoldarray[i]) - a*a*0.5*lambda * (uoldarray[i] - uoldarray[i-1]);
286*c4762a1bSJed Brown     f[i]  = uarray[i] - uoldarray[i] + lambda * (RFlux - LFlux);
287*c4762a1bSJed Brown   }
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown   /* Restore vectors */
290*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,U,&uarray);CHKERRQ(ierr);
291*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localUold,&uoldarray);CHKERRQ(ierr);
292*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
293*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localUold);CHKERRQ(ierr);
294*c4762a1bSJed Brown   PetscFunctionReturn(0);
295*c4762a1bSJed Brown }
296*c4762a1bSJed Brown 
297*c4762a1bSJed Brown 
298*c4762a1bSJed Brown /*TEST
299*c4762a1bSJed Brown 
300*c4762a1bSJed Brown    test:
301*c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor
302*c4762a1bSJed Brown 
303*c4762a1bSJed Brown    test:
304*c4762a1bSJed Brown       suffix: 2
305*c4762a1bSJed Brown       nsize: 3
306*c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor
307*c4762a1bSJed Brown       output_file: output/ex6_1.out
308*c4762a1bSJed Brown 
309*c4762a1bSJed Brown    test:
310*c4762a1bSJed Brown       suffix: 3
311*c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false
312*c4762a1bSJed Brown 
313*c4762a1bSJed Brown    test:
314*c4762a1bSJed Brown       suffix: 4
315*c4762a1bSJed Brown       nsize: 3
316*c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_monitor -useLaxWendroff false
317*c4762a1bSJed Brown       output_file: output/ex6_3.out
318*c4762a1bSJed Brown 
319*c4762a1bSJed Brown TEST*/
320