xref: /petsc/src/ts/tutorials/ex41.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Parallel bouncing ball example to test TS event feature.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   The dynamics of the bouncing ball is described by the ODE
5c4762a1bSJed Brown                   u1_t = u2
6c4762a1bSJed Brown                   u2_t = -9.8
7c4762a1bSJed Brown 
8c4762a1bSJed Brown   Each processor is assigned one ball.
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   The event function routine checks for the ball hitting the
11c4762a1bSJed Brown   ground (u1 = 0). Every time the ball hits the ground, its velocity u2 is attenuated by
12c4762a1bSJed Brown   a factor of 0.9 and its height set to 1.0*rank.
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown #include <petscts.h>
16c4762a1bSJed Brown 
EventFunction(TS ts,PetscReal t,Vec U,PetscReal * fvalue,PetscCtx ctx)17*2a8381b2SBarry Smith PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, PetscCtx ctx)
18d71ae5a4SJacob Faibussowitsch {
19c4762a1bSJed Brown   const PetscScalar *u;
20c4762a1bSJed Brown 
217510d9b0SBarry Smith   PetscFunctionBeginUser;
22c4762a1bSJed Brown   /* Event for ball height */
239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
24fa9584fbSIlya Fursov   fvalue[0] = PetscRealPart(u[0]);
259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown 
PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,PetscCtx ctx)29*2a8381b2SBarry Smith PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, PetscCtx ctx)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown   PetscScalar *u;
32c4762a1bSJed Brown   PetscMPIInt  rank;
33c4762a1bSJed Brown 
347510d9b0SBarry Smith   PetscFunctionBeginUser;
359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
36c4762a1bSJed Brown   if (nevents) {
379566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball hit the ground at t = %5.2f seconds -> Processor[%d]\n", (double)t, rank));
38c4762a1bSJed Brown     /* Set new initial conditions with .9 attenuation */
399566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
40c4762a1bSJed Brown     u[0] = 1.0 * rank;
41c4762a1bSJed Brown     u[1] = -0.9 * u[1];
429566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
43c4762a1bSJed Brown   }
443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /*
48c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in explicit form: U_t = F(U)
49c4762a1bSJed Brown */
RHSFunction(TS ts,PetscReal t,Vec U,Vec F,PetscCtx ctx)50*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, PetscCtx ctx)
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown   PetscScalar       *f;
53c4762a1bSJed Brown   const PetscScalar *u;
54c4762a1bSJed Brown 
557510d9b0SBarry Smith   PetscFunctionBeginUser;
56c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   f[0] = u[1];
61c4762a1bSJed Brown   f[1] = -9.8;
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown /*
69c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning the Jacobian.
70c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,PetscCtx ctx)71*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx)
72d71ae5a4SJacob Faibussowitsch {
73c4762a1bSJed Brown   PetscInt           rowcol[2], rstart;
74c4762a1bSJed Brown   PetscScalar        J[2][2];
75c4762a1bSJed Brown   const PetscScalar *u;
76c4762a1bSJed Brown 
777510d9b0SBarry Smith   PetscFunctionBeginUser;
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
819371c9d4SSatish Balay   rowcol[0] = rstart;
829371c9d4SSatish Balay   rowcol[1] = rstart + 1;
83c4762a1bSJed Brown 
849371c9d4SSatish Balay   J[0][0] = 0.0;
859371c9d4SSatish Balay   J[0][1] = 1.0;
869371c9d4SSatish Balay   J[1][0] = 0.0;
879371c9d4SSatish Balay   J[1][1] = 0.0;
889566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown   if (A != B) {
949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
959566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
96c4762a1bSJed Brown   }
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown /*
101c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0
102c4762a1bSJed Brown */
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)103*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
104d71ae5a4SJacob Faibussowitsch {
105c4762a1bSJed Brown   PetscScalar       *f;
106c4762a1bSJed Brown   const PetscScalar *u, *udot;
107c4762a1bSJed Brown 
1087510d9b0SBarry Smith   PetscFunctionBeginUser;
109c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   f[0] = udot[0] - u[1];
115c4762a1bSJed Brown   f[1] = udot[1] + 9.8;
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown /*
124c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
125c4762a1bSJed Brown */
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)126*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
127d71ae5a4SJacob Faibussowitsch {
128c4762a1bSJed Brown   PetscInt           rowcol[2], rstart;
129c4762a1bSJed Brown   PetscScalar        J[2][2];
130c4762a1bSJed Brown   const PetscScalar *u, *udot;
131c4762a1bSJed Brown 
1327510d9b0SBarry Smith   PetscFunctionBeginUser;
1339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
135c4762a1bSJed Brown 
1369566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
1379371c9d4SSatish Balay   rowcol[0] = rstart;
1389371c9d4SSatish Balay   rowcol[1] = rstart + 1;
139c4762a1bSJed Brown 
1409371c9d4SSatish Balay   J[0][0] = a;
1419371c9d4SSatish Balay   J[0][1] = -1.0;
1429371c9d4SSatish Balay   J[1][0] = 0.0;
1439371c9d4SSatish Balay   J[1][1] = a;
1449566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
148c4762a1bSJed Brown 
1499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
151c4762a1bSJed Brown   if (A != B) {
1529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
154c4762a1bSJed Brown   }
1553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
156c4762a1bSJed Brown }
157c4762a1bSJed Brown 
main(int argc,char ** argv)158d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
159d71ae5a4SJacob Faibussowitsch {
160c4762a1bSJed Brown   TS           ts; /* ODE integrator */
161c4762a1bSJed Brown   Vec          U;  /* solution will be stored here */
162c4762a1bSJed Brown   PetscMPIInt  rank;
163c4762a1bSJed Brown   PetscInt     n = 2;
164c4762a1bSJed Brown   PetscScalar *u;
165c4762a1bSJed Brown   PetscInt     direction = -1;
166c4762a1bSJed Brown   PetscBool    terminate = PETSC_FALSE;
167c4762a1bSJed Brown   PetscBool    rhs_form = PETSC_FALSE, hist = PETSC_TRUE;
168c4762a1bSJed Brown   TSAdapt      adapt;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown      Initialize program
172c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173327415f7SBarry Smith   PetscFunctionBeginUser;
174c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown      Create timestepping solver context
179c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1809566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1819566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184c4762a1bSJed Brown      Set ODE routines
185c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1869566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
187c4762a1bSJed Brown   /* Users are advised against the following branching and code duplication.
188c4762a1bSJed Brown      For problems without a mass matrix like the one at hand, the RHSFunction
189c4762a1bSJed Brown      (and companion RHSJacobian) interface is enough to support both explicit
190c4762a1bSJed Brown      and implicit timesteppers. This tutorial example also deals with the
191c4762a1bSJed Brown      IFunction/IJacobian interface for demonstration and testing purposes. */
1929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-rhs-form", &rhs_form, NULL));
193c4762a1bSJed Brown   if (rhs_form) {
1949566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
1959566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, NULL));
196c4762a1bSJed Brown   } else {
1979566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, NULL));
1989566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, NULL));
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown      Set initial conditions
203c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2049566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &U));
2059566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(U, n, PETSC_DETERMINE));
2069566063dSJacob Faibussowitsch   PetscCall(VecSetUp(U));
2079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
208c4762a1bSJed Brown   u[0] = 1.0 * rank;
209c4762a1bSJed Brown   u[1] = 20.0;
2109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
2119566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214c4762a1bSJed Brown      Set solver options
215c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2169566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
2179566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 30.0));
2189566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2199566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.1));
220a5b23f4aSJose E. Roman   /* The adaptive time step controller could take very large timesteps resulting in
221a5b23f4aSJose E. Roman      the same event occurring multiple times in the same interval. A maximum step size
222c4762a1bSJed Brown      limit is enforced here to avoid this issue. */
2239566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
2249566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetType(adapt, TSADAPTBASIC));
2259566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* Set direction and terminate flag for the event */
2289566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 1, &direction, &terminate, EventFunction, PostEventFunction, NULL));
229c4762a1bSJed Brown 
2309566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233c4762a1bSJed Brown      Run timestepping solver
234c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2359566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   if (hist) { /* replay following history */
238c4762a1bSJed Brown     TSTrajectory tj;
239c4762a1bSJed Brown     PetscReal    tf, t0, dt;
240c4762a1bSJed Brown 
2419566063dSJacob Faibussowitsch     PetscCall(TSGetTime(ts, &tf));
2429566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, tf));
2439566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, 0));
2449566063dSJacob Faibussowitsch     PetscCall(TSRestartStep(ts));
2459566063dSJacob Faibussowitsch     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2469566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
2479566063dSJacob Faibussowitsch     PetscCall(TSSetEventHandler(ts, 1, &direction, &terminate, EventFunction, PostEventFunction, NULL));
2489566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
2499566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetType(adapt, TSADAPTHISTORY));
2509566063dSJacob Faibussowitsch     PetscCall(TSGetTrajectory(ts, &tj));
2519566063dSJacob Faibussowitsch     PetscCall(TSAdaptHistorySetTrajectory(adapt, tj, PETSC_FALSE));
2529566063dSJacob Faibussowitsch     PetscCall(TSAdaptHistoryGetStep(adapt, 0, &t0, &dt));
253c4762a1bSJed Brown     /* this example fails with single (or smaller) precision */
254cc404788SPierre Jolivet #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL___FP16)
2559566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetType(adapt, TSADAPTBASIC));
2569566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
2579566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
258c4762a1bSJed Brown #endif
2599566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, t0));
2609566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
2619566063dSJacob Faibussowitsch     PetscCall(TSResetTrajectory(ts));
2629566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
263c4762a1bSJed Brown     u[0] = 1.0 * rank;
264c4762a1bSJed Brown     u[1] = 20.0;
2659566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
2669566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, U));
267c4762a1bSJed Brown   }
268c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
270c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2729566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
273c4762a1bSJed Brown 
2749566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
275b122ec5aSJacob Faibussowitsch   return 0;
276c4762a1bSJed Brown }
277c4762a1bSJed Brown 
278c4762a1bSJed Brown /*TEST
279c4762a1bSJed Brown 
280c4762a1bSJed Brown    test:
281c4762a1bSJed Brown       suffix: a
282c4762a1bSJed Brown       nsize: 2
283c4762a1bSJed Brown       args: -ts_trajectory_type memory -snes_stol 1e-4
284c4762a1bSJed Brown       filter: sort -b
285c4762a1bSJed Brown 
286c4762a1bSJed Brown    test:
287c4762a1bSJed Brown       suffix: b
288c4762a1bSJed Brown       nsize: 2
289c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type arkimex -snes_stol 1e-4
290c4762a1bSJed Brown       filter: sort -b
291c4762a1bSJed Brown 
292c4762a1bSJed Brown    test:
293c4762a1bSJed Brown       suffix: c
294c4762a1bSJed Brown       nsize: 2
295c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4
296c4762a1bSJed Brown       filter: sort -b
297c4762a1bSJed Brown 
298c4762a1bSJed Brown    test:
299c4762a1bSJed Brown       suffix: d
300c4762a1bSJed Brown       nsize: 2
301c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4
302c4762a1bSJed Brown       filter: sort -b
303c4762a1bSJed Brown 
304c4762a1bSJed Brown    test:
305c4762a1bSJed Brown       suffix: e
306c4762a1bSJed Brown       nsize: 2
307c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type bdf -ts_adapt_dt_max 0.015 -ts_max_steps 3000
308c4762a1bSJed Brown       filter: sort -b
309c4762a1bSJed Brown 
310c4762a1bSJed Brown    test:
311c4762a1bSJed Brown       suffix: f
312c4762a1bSJed Brown       nsize: 2
313c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 3bs
314c4762a1bSJed Brown       filter: sort -b
315c4762a1bSJed Brown 
316c4762a1bSJed Brown    test:
317c4762a1bSJed Brown       suffix: g
318c4762a1bSJed Brown       nsize: 2
319c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 5bs
320c4762a1bSJed Brown       filter: sort -b
321c4762a1bSJed Brown 
322c4762a1bSJed Brown    test:
323c4762a1bSJed Brown       suffix: h
324c4762a1bSJed Brown       nsize: 2
325c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 6vr
326c4762a1bSJed Brown       filter: sort -b
327c4762a1bSJed Brown       output_file: output/ex41_g.out
328c4762a1bSJed Brown 
329c4762a1bSJed Brown    test:
330c4762a1bSJed Brown       suffix: i
331c4762a1bSJed Brown       nsize: 2
332c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 7vr
333c4762a1bSJed Brown       filter: sort -b
334c4762a1bSJed Brown       output_file: output/ex41_g.out
335c4762a1bSJed Brown 
336c4762a1bSJed Brown    test:
337c4762a1bSJed Brown       suffix: j
338c4762a1bSJed Brown       nsize: 2
339c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 8vr
340c4762a1bSJed Brown       filter: sort -b
341c4762a1bSJed Brown       output_file: output/ex41_g.out
342c4762a1bSJed Brown 
343c4762a1bSJed Brown TEST*/
344