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