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
main(int argc,char ** argv)39d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
40d71ae5a4SJacob Faibussowitsch {
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
52327415f7SBarry Smith PetscFunctionBeginUser;
53c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
54c4762a1bSJed Brown
55c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
57c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
589566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 2, 2, NULL, &da));
599566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
609566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
61c4762a1bSJed Brown
62c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63c4762a1bSJed Brown Extract global vectors from DMDA;
64c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
659566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &X));
66c4762a1bSJed Brown
67c4762a1bSJed Brown /* Initialize user application context */
68d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Advection-reaction options", "");
69c4762a1bSJed Brown {
709371c9d4SSatish Balay user.a[0] = 1;
719371c9d4SSatish Balay PetscCall(PetscOptionsReal("-a0", "Advection rate 0", "", user.a[0], &user.a[0], NULL));
729371c9d4SSatish Balay user.a[1] = 0;
739371c9d4SSatish Balay PetscCall(PetscOptionsReal("-a1", "Advection rate 1", "", user.a[1], &user.a[1], NULL));
749371c9d4SSatish Balay user.k[0] = 1e6;
759371c9d4SSatish Balay PetscCall(PetscOptionsReal("-k0", "Reaction rate 0", "", user.k[0], &user.k[0], NULL));
769371c9d4SSatish Balay user.k[1] = 2 * user.k[0];
779371c9d4SSatish Balay PetscCall(PetscOptionsReal("-k1", "Reaction rate 1", "", user.k[1], &user.k[1], NULL));
789371c9d4SSatish Balay user.s[0] = 0;
799371c9d4SSatish Balay PetscCall(PetscOptionsReal("-s0", "Source 0", "", user.s[0], &user.s[0], NULL));
809371c9d4SSatish Balay user.s[1] = 1;
819371c9d4SSatish Balay PetscCall(PetscOptionsReal("-s1", "Source 1", "", user.s[1], &user.s[1], NULL));
82c4762a1bSJed Brown }
83d0609cedSBarry Smith PetscOptionsEnd();
84c4762a1bSJed Brown
85c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86c4762a1bSJed Brown Create timestepping solver context
87c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
889566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
899566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
909566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX));
919566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user));
929566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &user));
939566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATAIJ));
949566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &J));
959566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user));
96c4762a1bSJed Brown
97c4762a1bSJed Brown /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since
98c4762a1bSJed Brown * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use
99c4762a1bSJed Brown * SNESSetType(snes,SNESKSPONLY). */
1009566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
1019566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch));
1029566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
103c4762a1bSJed Brown
104c4762a1bSJed Brown ftime = .1;
1059566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime));
1069566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
107c4762a1bSJed Brown
108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109c4762a1bSJed Brown Set initial conditions
110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1119566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts, X, &user));
1129566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X));
1139566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &mx));
114c4762a1bSJed Brown dt = .1 * PetscMax(user.a[0], user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
1159566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt));
116c4762a1bSJed Brown
117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118c4762a1bSJed Brown Set runtime options
119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1209566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
121c4762a1bSJed Brown
122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123c4762a1bSJed Brown Solve nonlinear system
124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1259566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
1269566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
1279566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
1289566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason));
12963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps));
130c4762a1bSJed Brown
131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132c4762a1bSJed Brown Free work space.
133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
1359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
1369566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
1389566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
139b122ec5aSJacob Faibussowitsch return 0;
140c4762a1bSJed Brown }
141c4762a1bSJed Brown
FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void * ptr)142d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ptr)
143d71ae5a4SJacob Faibussowitsch {
144c4762a1bSJed Brown User user = (User)ptr;
145c4762a1bSJed Brown DM da;
146c4762a1bSJed Brown DMDALocalInfo info;
147c4762a1bSJed Brown PetscInt i;
148c4762a1bSJed Brown Field *f;
149c4762a1bSJed Brown const Field *x, *xdot;
150c4762a1bSJed Brown
151c4762a1bSJed Brown PetscFunctionBeginUser;
1529566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
1539566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
154c4762a1bSJed Brown
155c4762a1bSJed Brown /* Get pointers to vector data */
1569566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&x));
1579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xdot, (void *)&xdot));
1589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
159c4762a1bSJed Brown
160c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
161c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) {
162c4762a1bSJed Brown f[i][0] = xdot[i][0] + user->k[0] * x[i][0] - user->k[1] * x[i][1] - user->s[0];
163c4762a1bSJed Brown f[i][1] = xdot[i][1] - user->k[0] * x[i][0] + user->k[1] * x[i][1] - user->s[1];
164c4762a1bSJed Brown }
165c4762a1bSJed Brown
166c4762a1bSJed Brown /* Restore vectors */
1679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&x));
1689566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xdot, (void *)&xdot));
1699566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
171c4762a1bSJed Brown }
172c4762a1bSJed Brown
FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void * ptr)173d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr)
174d71ae5a4SJacob Faibussowitsch {
175c4762a1bSJed Brown User user = (User)ptr;
176c4762a1bSJed Brown DM da;
177c4762a1bSJed Brown Vec Xloc;
178c4762a1bSJed Brown DMDALocalInfo info;
179c4762a1bSJed Brown PetscInt i, j;
180c4762a1bSJed Brown PetscReal hx;
181c4762a1bSJed Brown Field *f;
182c4762a1bSJed Brown const Field *x;
183c4762a1bSJed Brown
184c4762a1bSJed Brown PetscFunctionBeginUser;
1859566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
1869566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
187c4762a1bSJed Brown hx = 1.0 / (PetscReal)info.mx;
188c4762a1bSJed Brown
189c4762a1bSJed Brown /*
190c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
191c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
192c4762a1bSJed Brown By placing code between these two statements, computations can be
193c4762a1bSJed Brown done while messages are in transition.
194c4762a1bSJed Brown */
1959566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc));
1969566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
1979566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
198c4762a1bSJed Brown
199c4762a1bSJed Brown /* Get pointers to vector data */
2009566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xloc, (void *)&x));
2019566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
202c4762a1bSJed Brown
203c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
204c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) {
205c4762a1bSJed Brown const PetscReal *a = user->a;
206c4762a1bSJed Brown PetscReal u0t[2];
207c4762a1bSJed Brown u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12 * t), 4);
208c4762a1bSJed Brown u0t[1] = 0.0;
209c4762a1bSJed Brown for (j = 0; j < 2; j++) {
210c4762a1bSJed 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]);
211c4762a1bSJed 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]);
212c4762a1bSJed 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]);
213c4762a1bSJed Brown else if (i == info.mx - 1) f[i][j] = a[j] / hx * (-x[i][j] + x[i - 1][j]);
214c4762a1bSJed 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]);
215c4762a1bSJed Brown }
216c4762a1bSJed Brown }
217c4762a1bSJed Brown
218c4762a1bSJed Brown /* Restore vectors */
2199566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xloc, (void *)&x));
2209566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
2219566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc));
2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
223c4762a1bSJed Brown }
224c4762a1bSJed Brown
225c4762a1bSJed Brown /* --------------------------------------------------------------------- */
226c4762a1bSJed Brown /*
227c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
228c4762a1bSJed Brown */
FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void * ptr)229d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ptr)
230d71ae5a4SJacob Faibussowitsch {
231c4762a1bSJed Brown User user = (User)ptr;
232c4762a1bSJed Brown DMDALocalInfo info;
233c4762a1bSJed Brown PetscInt i;
234c4762a1bSJed Brown DM da;
235c4762a1bSJed Brown const Field *x, *xdot;
236c4762a1bSJed Brown
237c4762a1bSJed Brown PetscFunctionBeginUser;
2389566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
2399566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
240c4762a1bSJed Brown
241c4762a1bSJed Brown /* Get pointers to vector data */
2429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&x));
2439566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xdot, (void *)&xdot));
244c4762a1bSJed Brown
245c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
246c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) {
247c4762a1bSJed Brown const PetscReal *k = user->k;
248c4762a1bSJed Brown PetscScalar v[2][2];
249c4762a1bSJed Brown
2509371c9d4SSatish Balay v[0][0] = a + k[0];
2519371c9d4SSatish Balay v[0][1] = -k[1];
2529371c9d4SSatish Balay v[1][0] = -k[0];
2539371c9d4SSatish Balay v[1][1] = a + k[1];
2549566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre, 1, &i, 1, &i, &v[0][0], INSERT_VALUES));
255c4762a1bSJed Brown }
256c4762a1bSJed Brown
257c4762a1bSJed Brown /* Restore vectors */
2589566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&x));
2599566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xdot, (void *)&xdot));
260c4762a1bSJed Brown
2619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
2629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
263c4762a1bSJed Brown if (J != Jpre) {
2649566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2659566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
266c4762a1bSJed Brown }
2673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
268c4762a1bSJed Brown }
269c4762a1bSJed Brown
FormInitialSolution(TS ts,Vec X,PetscCtx ctx)270*2a8381b2SBarry Smith PetscErrorCode FormInitialSolution(TS ts, Vec X, PetscCtx ctx)
271d71ae5a4SJacob Faibussowitsch {
272c4762a1bSJed Brown User user = (User)ctx;
273c4762a1bSJed Brown DM da;
274c4762a1bSJed Brown PetscInt i;
275c4762a1bSJed Brown DMDALocalInfo info;
276c4762a1bSJed Brown Field *x;
277c4762a1bSJed Brown PetscReal hx;
278c4762a1bSJed Brown
279c4762a1bSJed Brown PetscFunctionBeginUser;
2809566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
2819566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
282c4762a1bSJed Brown hx = 1.0 / (PetscReal)info.mx;
283c4762a1bSJed Brown
284c4762a1bSJed Brown /* Get pointers to vector data */
2859566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, X, &x));
286c4762a1bSJed Brown
287c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
288c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) {
289c4762a1bSJed Brown PetscReal r = (i + 1) * hx, ik = user->k[1] != 0.0 ? 1.0 / user->k[1] : 1.0;
290c4762a1bSJed Brown x[i][0] = 1 + user->s[1] * r;
291c4762a1bSJed Brown x[i][1] = user->k[0] * ik * x[i][0] + user->s[1] * ik;
292c4762a1bSJed Brown }
2939566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, X, &x));
2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
295c4762a1bSJed Brown }
296c4762a1bSJed Brown
297c4762a1bSJed Brown /*TEST
298c4762a1bSJed Brown
299c4762a1bSJed Brown test:
300188af4bfSBarry Smith args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_arkimex_type 4 -ts_adapt_type none -ts_time_step .005 -ts_max_time .1
301c4762a1bSJed Brown requires: !single
302c4762a1bSJed Brown
303c4762a1bSJed Brown test:
304c4762a1bSJed Brown suffix: 2
305188af4bfSBarry Smith args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_time_step 1e-3 -ts_adapt_type none -ts_time_step .005 -ts_max_time .1
306c4762a1bSJed Brown nsize: 2
307c4762a1bSJed Brown
308c4762a1bSJed Brown test:
309c4762a1bSJed Brown suffix: 3
310188af4bfSBarry Smith args: -nox -da_grid_x 200 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type ra34pw2 -ts_time_step 5e-3 -ts_max_time .1 -ts_adapt_type none
311c4762a1bSJed Brown nsize: 2
312c4762a1bSJed Brown
313c4762a1bSJed Brown test:
314c4762a1bSJed Brown suffix: 4
315188af4bfSBarry Smith args: -ts_type eimex -da_grid_x 200 -ts_eimex_order_adapt true -ts_time_step 0.001 -ts_monitor -ts_max_steps 100
316c4762a1bSJed Brown filter: sed "s/ITS/TIME/g"
317c4762a1bSJed Brown nsize: 2
318c4762a1bSJed Brown
319c4762a1bSJed Brown TEST*/
320