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 39*9371c9d4SSatish Balay int main(int argc, char **argv) { 40c4762a1bSJed Brown TS ts; /* time integrator */ 41c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 42c4762a1bSJed Brown SNESLineSearch linesearch; /* line search */ 43c4762a1bSJed Brown Vec X; /* solution, residual vectors */ 44c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 45c4762a1bSJed Brown PetscInt steps, mx; 46c4762a1bSJed Brown DM da; 47c4762a1bSJed Brown PetscReal ftime, dt; 48c4762a1bSJed Brown struct _User user; /* user-defined work context */ 49c4762a1bSJed Brown TSConvergedReason reason; 50c4762a1bSJed Brown 51327415f7SBarry Smith PetscFunctionBeginUser; 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 */ 67d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Advection-reaction options", ""); 68c4762a1bSJed Brown { 69*9371c9d4SSatish Balay user.a[0] = 1; 70*9371c9d4SSatish Balay PetscCall(PetscOptionsReal("-a0", "Advection rate 0", "", user.a[0], &user.a[0], NULL)); 71*9371c9d4SSatish Balay user.a[1] = 0; 72*9371c9d4SSatish Balay PetscCall(PetscOptionsReal("-a1", "Advection rate 1", "", user.a[1], &user.a[1], NULL)); 73*9371c9d4SSatish Balay user.k[0] = 1e6; 74*9371c9d4SSatish Balay PetscCall(PetscOptionsReal("-k0", "Reaction rate 0", "", user.k[0], &user.k[0], NULL)); 75*9371c9d4SSatish Balay user.k[1] = 2 * user.k[0]; 76*9371c9d4SSatish Balay PetscCall(PetscOptionsReal("-k1", "Reaction rate 1", "", user.k[1], &user.k[1], NULL)); 77*9371c9d4SSatish Balay user.s[0] = 0; 78*9371c9d4SSatish Balay PetscCall(PetscOptionsReal("-s0", "Source 0", "", user.s[0], &user.s[0], NULL)); 79*9371c9d4SSatish Balay user.s[1] = 1; 80*9371c9d4SSatish Balay PetscCall(PetscOptionsReal("-s1", "Source 1", "", user.s[1], &user.s[1], NULL)); 81c4762a1bSJed Brown } 82d0609cedSBarry Smith PetscOptionsEnd(); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 85c4762a1bSJed Brown Create timestepping solver context 86c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 879566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 889566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 899566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX)); 909566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user)); 919566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &user)); 929566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATAIJ)); 939566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &J)); 949566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since 97c4762a1bSJed Brown * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use 98c4762a1bSJed Brown * SNESSetType(snes,SNESKSPONLY). */ 999566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 1009566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 1019566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown ftime = .1; 1049566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1059566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108c4762a1bSJed Brown Set initial conditions 109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1109566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts, X, &user)); 1119566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &mx)); 113c4762a1bSJed Brown dt = .1 * PetscMax(user.a[0], user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */ 1149566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117c4762a1bSJed Brown Set runtime options 118c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1199566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122c4762a1bSJed Brown Solve nonlinear system 123c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1249566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 1259566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 1269566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 1279566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 12863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131c4762a1bSJed Brown Free work space. 132c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1359566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1379566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 138b122ec5aSJacob Faibussowitsch return 0; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141*9371c9d4SSatish Balay static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ptr) { 142c4762a1bSJed Brown User user = (User)ptr; 143c4762a1bSJed Brown DM da; 144c4762a1bSJed Brown DMDALocalInfo info; 145c4762a1bSJed Brown PetscInt i; 146c4762a1bSJed Brown Field *f; 147c4762a1bSJed Brown const Field *x, *xdot; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBeginUser; 1509566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1519566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* Get pointers to vector data */ 1549566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&x)); 1559566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xdot, (void *)&xdot)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 159c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 160c4762a1bSJed Brown f[i][0] = xdot[i][0] + user->k[0] * x[i][0] - user->k[1] * x[i][1] - user->s[0]; 161c4762a1bSJed Brown f[i][1] = xdot[i][1] - user->k[0] * x[i][0] + user->k[1] * x[i][1] - user->s[1]; 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* Restore vectors */ 1659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&x)); 1669566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xdot, (void *)&xdot)); 1679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 168c4762a1bSJed Brown PetscFunctionReturn(0); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171*9371c9d4SSatish Balay static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr) { 172c4762a1bSJed Brown User user = (User)ptr; 173c4762a1bSJed Brown DM da; 174c4762a1bSJed Brown Vec Xloc; 175c4762a1bSJed Brown DMDALocalInfo info; 176c4762a1bSJed Brown PetscInt i, j; 177c4762a1bSJed Brown PetscReal hx; 178c4762a1bSJed Brown Field *f; 179c4762a1bSJed Brown const Field *x; 180c4762a1bSJed Brown 181c4762a1bSJed Brown PetscFunctionBeginUser; 1829566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1839566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 184c4762a1bSJed Brown hx = 1.0 / (PetscReal)info.mx; 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 188c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 189c4762a1bSJed Brown By placing code between these two statements, computations can be 190c4762a1bSJed Brown done while messages are in transition. 191c4762a1bSJed Brown */ 1929566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 1939566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 1949566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* Get pointers to vector data */ 1979566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xloc, (void *)&x)); 1989566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 201c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 202c4762a1bSJed Brown const PetscReal *a = user->a; 203c4762a1bSJed Brown PetscReal u0t[2]; 204c4762a1bSJed Brown u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12 * t), 4); 205c4762a1bSJed Brown u0t[1] = 0.0; 206c4762a1bSJed Brown for (j = 0; j < 2; j++) { 207c4762a1bSJed 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]); 208c4762a1bSJed 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]); 209c4762a1bSJed 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]); 210c4762a1bSJed Brown else if (i == info.mx - 1) f[i][j] = a[j] / hx * (-x[i][j] + x[i - 1][j]); 211c4762a1bSJed 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]); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* Restore vectors */ 2169566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xloc, (void *)&x)); 2179566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 2189566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 219c4762a1bSJed Brown PetscFunctionReturn(0); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* --------------------------------------------------------------------- */ 223c4762a1bSJed Brown /* 224c4762a1bSJed Brown IJacobian - Compute IJacobian = dF/dU + a dF/dUdot 225c4762a1bSJed Brown */ 226*9371c9d4SSatish Balay PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ptr) { 227c4762a1bSJed Brown User user = (User)ptr; 228c4762a1bSJed Brown DMDALocalInfo info; 229c4762a1bSJed Brown PetscInt i; 230c4762a1bSJed Brown DM da; 231c4762a1bSJed Brown const Field *x, *xdot; 232c4762a1bSJed Brown 233c4762a1bSJed Brown PetscFunctionBeginUser; 2349566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2359566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* Get pointers to vector data */ 2389566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, X, (void *)&x)); 2399566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xdot, (void *)&xdot)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 242c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 243c4762a1bSJed Brown const PetscReal *k = user->k; 244c4762a1bSJed Brown PetscScalar v[2][2]; 245c4762a1bSJed Brown 246*9371c9d4SSatish Balay v[0][0] = a + k[0]; 247*9371c9d4SSatish Balay v[0][1] = -k[1]; 248*9371c9d4SSatish Balay v[1][0] = -k[0]; 249*9371c9d4SSatish Balay v[1][1] = a + k[1]; 2509566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(Jpre, 1, &i, 1, &i, &v[0][0], INSERT_VALUES)); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* Restore vectors */ 2549566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, X, (void *)&x)); 2559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xdot, (void *)&xdot)); 256c4762a1bSJed Brown 2579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 2589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 259c4762a1bSJed Brown if (J != Jpre) { 2609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 2619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 262c4762a1bSJed Brown } 263c4762a1bSJed Brown PetscFunctionReturn(0); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266*9371c9d4SSatish Balay PetscErrorCode FormInitialSolution(TS ts, Vec X, void *ctx) { 267c4762a1bSJed Brown User user = (User)ctx; 268c4762a1bSJed Brown DM da; 269c4762a1bSJed Brown PetscInt i; 270c4762a1bSJed Brown DMDALocalInfo info; 271c4762a1bSJed Brown Field *x; 272c4762a1bSJed Brown PetscReal hx; 273c4762a1bSJed Brown 274c4762a1bSJed Brown PetscFunctionBeginUser; 2759566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2769566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 277c4762a1bSJed Brown hx = 1.0 / (PetscReal)info.mx; 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* Get pointers to vector data */ 2809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, X, &x)); 281c4762a1bSJed Brown 282c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 283c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 284c4762a1bSJed Brown PetscReal r = (i + 1) * hx, ik = user->k[1] != 0.0 ? 1.0 / user->k[1] : 1.0; 285c4762a1bSJed Brown x[i][0] = 1 + user->s[1] * r; 286c4762a1bSJed Brown x[i][1] = user->k[0] * ik * x[i][0] + user->s[1] * ik; 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, X, &x)); 289c4762a1bSJed Brown PetscFunctionReturn(0); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown /*TEST 293c4762a1bSJed Brown 294c4762a1bSJed Brown test: 295c4762a1bSJed 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 296c4762a1bSJed Brown requires: !single 297c4762a1bSJed Brown 298c4762a1bSJed Brown test: 299c4762a1bSJed Brown suffix: 2 300c4762a1bSJed 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 301c4762a1bSJed Brown nsize: 2 302c4762a1bSJed Brown 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 3 305c4762a1bSJed 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 306c4762a1bSJed Brown nsize: 2 307c4762a1bSJed Brown 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: 4 310c4762a1bSJed Brown args: -ts_type eimex -da_grid_x 200 -ts_eimex_order_adapt true -ts_dt 0.001 -ts_monitor -ts_max_steps 100 311c4762a1bSJed Brown filter: sed "s/ITS/TIME/g" 312c4762a1bSJed Brown nsize: 2 313c4762a1bSJed Brown 314c4762a1bSJed Brown TEST*/ 315