xref: /petsc/src/ts/tutorials/ex22.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown static const char help[] = "Time-dependent advection-reaction PDE in 1d, demonstrates IMEX methods.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    u_t + a1*u_x = -k1*u + k2*v + s1
4c4762a1bSJed Brown    v_t + a2*v_x = k1*u - k2*v + s2
5c4762a1bSJed Brown    0 < x < 1;
6c4762a1bSJed Brown    a1 = 1, k1 = 10^6, s1 = 0,
7c4762a1bSJed Brown    a2 = 0, k2 = 2*k1, s2 = 1
8c4762a1bSJed Brown 
9c4762a1bSJed Brown    Initial conditions:
10c4762a1bSJed Brown    u(x,0) = 1 + s2*x
11c4762a1bSJed Brown    v(x,0) = k0/k1*u(x,0) + s1/k1
12c4762a1bSJed Brown 
13c4762a1bSJed Brown    Upstream boundary conditions:
14c4762a1bSJed Brown    u(0,t) = 1-sin(12*t)^4
15c4762a1bSJed Brown 
16c4762a1bSJed Brown    Useful command-line parameters:
17c4762a1bSJed Brown    -ts_arkimex_fully_implicit - treats advection implicitly, only relevant with -ts_type arkimex (default)
18c4762a1bSJed Brown    -ts_type rosw - use Rosenbrock-W method (linearly implicit IMEX, amortizes assembly and preconditioner setup)
19c4762a1bSJed Brown */
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #include <petscdm.h>
22c4762a1bSJed Brown #include <petscdmda.h>
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown 
25c4762a1bSJed Brown typedef PetscScalar Field[2];
26c4762a1bSJed Brown 
27c4762a1bSJed Brown typedef struct _User *User;
28c4762a1bSJed Brown struct _User {
29c4762a1bSJed Brown   PetscReal a[2];              /* Advection speeds */
30c4762a1bSJed Brown   PetscReal k[2];              /* Reaction coefficients */
31c4762a1bSJed Brown   PetscReal s[2];              /* Source terms */
32c4762a1bSJed Brown };
33c4762a1bSJed Brown 
34c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
35c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
36c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
37c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*);
38c4762a1bSJed Brown 
39c4762a1bSJed Brown int main(int argc,char **argv)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   TS                ts;         /* time integrator */
42c4762a1bSJed Brown   SNES              snes;       /* nonlinear solver */
43c4762a1bSJed Brown   SNESLineSearch    linesearch; /* line search */
44c4762a1bSJed Brown   Vec               X;          /* solution, residual vectors */
45c4762a1bSJed Brown   Mat               J;          /* Jacobian matrix */
46c4762a1bSJed Brown   PetscInt          steps,mx;
47c4762a1bSJed Brown   DM                da;
48c4762a1bSJed Brown   PetscReal         ftime,dt;
49c4762a1bSJed Brown   struct _User      user;       /* user-defined work context */
50c4762a1bSJed Brown   TSConvergedReason reason;
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
56c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
579566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,2,2,NULL,&da));
589566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
599566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62c4762a1bSJed Brown      Extract global vectors from DMDA;
63c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
649566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&X));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* Initialize user application context */
67*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
68c4762a1bSJed Brown   {
699566063dSJacob Faibussowitsch     user.a[0] = 1;           PetscCall(PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL));
709566063dSJacob Faibussowitsch     user.a[1] = 0;           PetscCall(PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL));
719566063dSJacob Faibussowitsch     user.k[0] = 1e6;         PetscCall(PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL));
729566063dSJacob Faibussowitsch     user.k[1] = 2*user.k[0]; PetscCall(PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL));
739566063dSJacob Faibussowitsch     user.s[0] = 0;           PetscCall(PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL));
749566063dSJacob Faibussowitsch     user.s[1] = 1;           PetscCall(PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL));
75c4762a1bSJed Brown   }
76*d0609cedSBarry Smith   PetscOptionsEnd();
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79c4762a1bSJed Brown      Create timestepping solver context
80c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
819566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
829566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,da));
839566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSARKIMEX));
849566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,FormRHSFunction,&user));
859566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,FormIFunction,&user));
869566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da,MATAIJ));
879566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da,&J));
889566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,&user));
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since
91c4762a1bSJed Brown    * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use
92c4762a1bSJed Brown    * SNESSetType(snes,SNESKSPONLY). */
939566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&snes));
949566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes,&linesearch));
959566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   ftime = .1;
989566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
999566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102c4762a1bSJed Brown      Set initial conditions
103c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1049566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts,X,&user));
1059566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,X));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X,&mx));
107c4762a1bSJed Brown   dt   = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
1089566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4762a1bSJed Brown      Set runtime options
112c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1139566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown      Solve nonlinear system
117c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1189566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,X));
1199566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
1209566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
1219566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts,&reason));
1229566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125c4762a1bSJed Brown      Free work space.
126c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1299566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1309566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1319566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
132b122ec5aSJacob Faibussowitsch   return 0;
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   User           user = (User)ptr;
138c4762a1bSJed Brown   DM             da;
139c4762a1bSJed Brown   DMDALocalInfo  info;
140c4762a1bSJed Brown   PetscInt       i;
141c4762a1bSJed Brown   Field          *f;
142c4762a1bSJed Brown   const Field    *x,*xdot;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   PetscFunctionBeginUser;
1459566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
1469566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Get pointers to vector data */
1499566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,X,(void*)&x));
1509566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Xdot,(void*)&xdot));
1519566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
154c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
155c4762a1bSJed Brown     f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0];
156c4762a1bSJed Brown     f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1];
157c4762a1bSJed Brown   }
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Restore vectors */
1609566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&x));
1619566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot));
1629566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
163c4762a1bSJed Brown   PetscFunctionReturn(0);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
167c4762a1bSJed Brown {
168c4762a1bSJed Brown   User           user = (User)ptr;
169c4762a1bSJed Brown   DM             da;
170c4762a1bSJed Brown   Vec            Xloc;
171c4762a1bSJed Brown   DMDALocalInfo  info;
172c4762a1bSJed Brown   PetscInt       i,j;
173c4762a1bSJed Brown   PetscReal      hx;
174c4762a1bSJed Brown   Field          *f;
175c4762a1bSJed Brown   const Field    *x;
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   PetscFunctionBeginUser;
1789566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
1799566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
180c4762a1bSJed Brown   hx   = 1.0/(PetscReal)info.mx;
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /*
183c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
184c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
185c4762a1bSJed Brown      By placing code between these two statements, computations can be
186c4762a1bSJed Brown      done while messages are in transition.
187c4762a1bSJed Brown   */
1889566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&Xloc));
1899566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc));
1909566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* Get pointers to vector data */
1939566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Xloc,(void*)&x));
1949566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
197c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
198c4762a1bSJed Brown     const PetscReal *a  = user->a;
199c4762a1bSJed Brown     PetscReal       u0t[2];
200c4762a1bSJed Brown     u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12*t),4);
201c4762a1bSJed Brown     u0t[1] = 0.0;
202c4762a1bSJed Brown     for (j=0; j<2; j++) {
203c4762a1bSJed Brown       if (i == 0)              f[i][j] = a[j]/hx*(1./3*u0t[j] + 0.5*x[i][j] - x[i+1][j] + 1./6*x[i+2][j]);
204c4762a1bSJed Brown       else if (i == 1)         f[i][j] = a[j]/hx*(-1./12*u0t[j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
205c4762a1bSJed Brown       else if (i == info.mx-2) f[i][j] = a[j]/hx*(-1./6*x[i-2][j] + x[i-1][j] - 0.5*x[i][j] - 1./3*x[i+1][j]);
206c4762a1bSJed Brown       else if (i == info.mx-1) f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]);
207c4762a1bSJed Brown       else                     f[i][j] = a[j]/hx*(-1./12*x[i-2][j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
208c4762a1bSJed Brown     }
209c4762a1bSJed Brown   }
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* Restore vectors */
2129566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Xloc,(void*)&x));
2139566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
2149566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&Xloc));
215c4762a1bSJed Brown   PetscFunctionReturn(0);
216c4762a1bSJed Brown }
217c4762a1bSJed Brown 
218c4762a1bSJed Brown /* --------------------------------------------------------------------- */
219c4762a1bSJed Brown /*
220c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
221c4762a1bSJed Brown */
222c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
223c4762a1bSJed Brown {
224c4762a1bSJed Brown   User           user = (User)ptr;
225c4762a1bSJed Brown   DMDALocalInfo  info;
226c4762a1bSJed Brown   PetscInt       i;
227c4762a1bSJed Brown   DM             da;
228c4762a1bSJed Brown   const Field    *x,*xdot;
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   PetscFunctionBeginUser;
2319566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
2329566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   /* Get pointers to vector data */
2359566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,X,(void*)&x));
2369566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Xdot,(void*)&xdot));
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
239c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
240c4762a1bSJed Brown     const PetscReal *k = user->k;
241c4762a1bSJed Brown     PetscScalar     v[2][2];
242c4762a1bSJed Brown 
243c4762a1bSJed Brown     v[0][0] = a + k[0]; v[0][1] =  -k[1];
244c4762a1bSJed Brown     v[1][0] =    -k[0]; v[1][1] = a+k[1];
2459566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES));
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /* Restore vectors */
2499566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,X,(void*)&x));
2509566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot));
251c4762a1bSJed Brown 
2529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
254c4762a1bSJed Brown   if (J != Jpre) {
2559566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
2569566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
257c4762a1bSJed Brown   }
258c4762a1bSJed Brown   PetscFunctionReturn(0);
259c4762a1bSJed Brown }
260c4762a1bSJed Brown 
261c4762a1bSJed Brown PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
262c4762a1bSJed Brown {
263c4762a1bSJed Brown   User           user = (User)ctx;
264c4762a1bSJed Brown   DM             da;
265c4762a1bSJed Brown   PetscInt       i;
266c4762a1bSJed Brown   DMDALocalInfo  info;
267c4762a1bSJed Brown   Field          *x;
268c4762a1bSJed Brown   PetscReal      hx;
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   PetscFunctionBeginUser;
2719566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
2729566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
273c4762a1bSJed Brown   hx   = 1.0/(PetscReal)info.mx;
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /* Get pointers to vector data */
2769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
279c4762a1bSJed Brown   for (i=info.xs; i<info.xs+info.xm; i++) {
280c4762a1bSJed Brown     PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0;
281c4762a1bSJed Brown     x[i][0] = 1 + user->s[1]*r;
282c4762a1bSJed Brown     x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik;
283c4762a1bSJed Brown   }
2849566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,X,&x));
285c4762a1bSJed Brown   PetscFunctionReturn(0);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown /*TEST
289c4762a1bSJed Brown 
290c4762a1bSJed Brown     test:
291c4762a1bSJed Brown       args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_arkimex_type 4 -ts_adapt_type none -ts_dt .005 -ts_max_time .1
292c4762a1bSJed Brown       requires: !single
293c4762a1bSJed Brown 
294c4762a1bSJed Brown     test:
295c4762a1bSJed Brown       suffix: 2
296c4762a1bSJed Brown       args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_dt 1e-3 -ts_adapt_type none -ts_dt .005 -ts_max_time .1
297c4762a1bSJed Brown       nsize: 2
298c4762a1bSJed Brown 
299c4762a1bSJed Brown     test:
300c4762a1bSJed Brown       suffix: 3
301c4762a1bSJed Brown       args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type ra34pw2 -ts_dt 5e-3 -ts_max_time .1  -ts_adapt_type none
302c4762a1bSJed Brown       nsize: 2
303c4762a1bSJed Brown 
304c4762a1bSJed Brown     test:
305c4762a1bSJed Brown       suffix: 4
306c4762a1bSJed Brown       args: -ts_type eimex -da_grid_x 200 -ts_eimex_order_adapt true -ts_dt 0.001 -ts_monitor -ts_max_steps 100
307c4762a1bSJed Brown       filter: sed "s/ITS/TIME/g"
308c4762a1bSJed Brown       nsize: 2
309c4762a1bSJed Brown 
310c4762a1bSJed Brown TEST*/
311