xref: /petsc/src/ts/tutorials/ex41.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
17*9371c9d4SSatish Balay PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) {
18c4762a1bSJed Brown   const PetscScalar *u;
19c4762a1bSJed Brown 
207510d9b0SBarry Smith   PetscFunctionBeginUser;
21c4762a1bSJed Brown   /* Event for ball height */
229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
23c4762a1bSJed Brown   fvalue[0] = u[0];
249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
25c4762a1bSJed Brown   PetscFunctionReturn(0);
26c4762a1bSJed Brown }
27c4762a1bSJed Brown 
28*9371c9d4SSatish Balay PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) {
29c4762a1bSJed Brown   PetscScalar *u;
30c4762a1bSJed Brown   PetscMPIInt  rank;
31c4762a1bSJed Brown 
327510d9b0SBarry Smith   PetscFunctionBeginUser;
339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
34c4762a1bSJed Brown   if (nevents) {
359566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball hit the ground at t = %5.2f seconds -> Processor[%d]\n", (double)t, rank));
36c4762a1bSJed Brown     /* Set new initial conditions with .9 attenuation */
379566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
38c4762a1bSJed Brown     u[0] = 1.0 * rank;
39c4762a1bSJed Brown     u[1] = -0.9 * u[1];
409566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown   PetscFunctionReturn(0);
43c4762a1bSJed Brown }
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /*
46c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in explicit form: U_t = F(U)
47c4762a1bSJed Brown */
48*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) {
49c4762a1bSJed Brown   PetscScalar       *f;
50c4762a1bSJed Brown   const PetscScalar *u;
51c4762a1bSJed Brown 
527510d9b0SBarry Smith   PetscFunctionBeginUser;
53c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   f[0] = u[1];
58c4762a1bSJed Brown   f[1] = -9.8;
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
62c4762a1bSJed Brown   PetscFunctionReturn(0);
63c4762a1bSJed Brown }
64c4762a1bSJed Brown 
65c4762a1bSJed Brown /*
66c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning the Jacobian.
67c4762a1bSJed Brown */
68*9371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) {
69c4762a1bSJed Brown   PetscInt           rowcol[2], rstart;
70c4762a1bSJed Brown   PetscScalar        J[2][2];
71c4762a1bSJed Brown   const PetscScalar *u;
72c4762a1bSJed Brown 
737510d9b0SBarry Smith   PetscFunctionBeginUser;
749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
77*9371c9d4SSatish Balay   rowcol[0] = rstart;
78*9371c9d4SSatish Balay   rowcol[1] = rstart + 1;
79c4762a1bSJed Brown 
80*9371c9d4SSatish Balay   J[0][0] = 0.0;
81*9371c9d4SSatish Balay   J[0][1] = 1.0;
82*9371c9d4SSatish Balay   J[1][0] = 0.0;
83*9371c9d4SSatish Balay   J[1][1] = 0.0;
849566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
89c4762a1bSJed Brown   if (A != B) {
909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
92c4762a1bSJed Brown   }
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown /*
97c4762a1bSJed Brown      Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0
98c4762a1bSJed Brown */
99*9371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
100c4762a1bSJed Brown   PetscScalar       *f;
101c4762a1bSJed Brown   const PetscScalar *u, *udot;
102c4762a1bSJed Brown 
1037510d9b0SBarry Smith   PetscFunctionBeginUser;
104c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   f[0] = udot[0] - u[1];
110c4762a1bSJed Brown   f[1] = udot[1] + 9.8;
111c4762a1bSJed Brown 
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
115c4762a1bSJed Brown   PetscFunctionReturn(0);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown /*
119c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
120c4762a1bSJed Brown */
121*9371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) {
122c4762a1bSJed Brown   PetscInt           rowcol[2], rstart;
123c4762a1bSJed Brown   PetscScalar        J[2][2];
124c4762a1bSJed Brown   const PetscScalar *u, *udot;
125c4762a1bSJed Brown 
1267510d9b0SBarry Smith   PetscFunctionBeginUser;
1279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
129c4762a1bSJed Brown 
1309566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
131*9371c9d4SSatish Balay   rowcol[0] = rstart;
132*9371c9d4SSatish Balay   rowcol[1] = rstart + 1;
133c4762a1bSJed Brown 
134*9371c9d4SSatish Balay   J[0][0] = a;
135*9371c9d4SSatish Balay   J[0][1] = -1.0;
136*9371c9d4SSatish Balay   J[1][0] = 0.0;
137*9371c9d4SSatish Balay   J[1][1] = a;
1389566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
145c4762a1bSJed Brown   if (A != B) {
1469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
148c4762a1bSJed Brown   }
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152*9371c9d4SSatish Balay int main(int argc, char **argv) {
153c4762a1bSJed Brown   TS           ts; /* ODE integrator */
154c4762a1bSJed Brown   Vec          U;  /* solution will be stored here */
155c4762a1bSJed Brown   PetscMPIInt  rank;
156c4762a1bSJed Brown   PetscInt     n = 2;
157c4762a1bSJed Brown   PetscScalar *u;
158c4762a1bSJed Brown   PetscInt     direction = -1;
159c4762a1bSJed Brown   PetscBool    terminate = PETSC_FALSE;
160c4762a1bSJed Brown   PetscBool    rhs_form = PETSC_FALSE, hist = PETSC_TRUE;
161c4762a1bSJed Brown   TSAdapt      adapt;
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164c4762a1bSJed Brown      Initialize program
165c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166327415f7SBarry Smith   PetscFunctionBeginUser;
1679566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown      Create timestepping solver context
172c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1739566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1749566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177c4762a1bSJed Brown      Set ODE routines
178c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1799566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
180c4762a1bSJed Brown   /* Users are advised against the following branching and code duplication.
181c4762a1bSJed Brown      For problems without a mass matrix like the one at hand, the RHSFunction
182c4762a1bSJed Brown      (and companion RHSJacobian) interface is enough to support both explicit
183c4762a1bSJed Brown      and implicit timesteppers. This tutorial example also deals with the
184c4762a1bSJed Brown      IFunction/IJacobian interface for demonstration and testing purposes. */
1859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-rhs-form", &rhs_form, NULL));
186c4762a1bSJed Brown   if (rhs_form) {
1879566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL));
1889566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, NULL));
189c4762a1bSJed Brown   } else {
1909566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, NULL));
1919566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, NULL));
192c4762a1bSJed Brown   }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195c4762a1bSJed Brown      Set initial conditions
196c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1979566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &U));
1989566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(U, n, PETSC_DETERMINE));
1999566063dSJacob Faibussowitsch   PetscCall(VecSetUp(U));
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
201c4762a1bSJed Brown   u[0] = 1.0 * rank;
202c4762a1bSJed Brown   u[1] = 20.0;
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
2049566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207c4762a1bSJed Brown      Set solver options
208c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2099566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
2109566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 30.0));
2119566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2129566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.1));
213a5b23f4aSJose E. Roman   /* The adaptive time step controller could take very large timesteps resulting in
214a5b23f4aSJose E. Roman      the same event occurring multiple times in the same interval. A maximum step size
215c4762a1bSJed Brown      limit is enforced here to avoid this issue. */
2169566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
2179566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetType(adapt, TSADAPTBASIC));
2189566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* Set direction and terminate flag for the event */
2219566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 1, &direction, &terminate, EventFunction, PostEventFunction, NULL));
222c4762a1bSJed Brown 
2239566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226c4762a1bSJed Brown      Run timestepping solver
227c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2289566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   if (hist) { /* replay following history */
231c4762a1bSJed Brown     TSTrajectory tj;
232c4762a1bSJed Brown     PetscReal    tf, t0, dt;
233c4762a1bSJed Brown 
2349566063dSJacob Faibussowitsch     PetscCall(TSGetTime(ts, &tf));
2359566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts, tf));
2369566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, 0));
2379566063dSJacob Faibussowitsch     PetscCall(TSRestartStep(ts));
2389566063dSJacob Faibussowitsch     PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2399566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
2409566063dSJacob Faibussowitsch     PetscCall(TSSetEventHandler(ts, 1, &direction, &terminate, EventFunction, PostEventFunction, NULL));
2419566063dSJacob Faibussowitsch     PetscCall(TSGetAdapt(ts, &adapt));
2429566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetType(adapt, TSADAPTHISTORY));
2439566063dSJacob Faibussowitsch     PetscCall(TSGetTrajectory(ts, &tj));
2449566063dSJacob Faibussowitsch     PetscCall(TSAdaptHistorySetTrajectory(adapt, tj, PETSC_FALSE));
2459566063dSJacob Faibussowitsch     PetscCall(TSAdaptHistoryGetStep(adapt, 0, &t0, &dt));
246c4762a1bSJed Brown     /* this example fails with single (or smaller) precision */
247c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL__FP16)
2489566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetType(adapt, TSADAPTBASIC));
2499566063dSJacob Faibussowitsch     PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5));
2509566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
251c4762a1bSJed Brown #endif
2529566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, t0));
2539566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
2549566063dSJacob Faibussowitsch     PetscCall(TSResetTrajectory(ts));
2559566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
256c4762a1bSJed Brown     u[0] = 1.0 * rank;
257c4762a1bSJed Brown     u[1] = 20.0;
2589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
2599566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, U));
260c4762a1bSJed Brown   }
261c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
263c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2659566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
266c4762a1bSJed Brown 
2679566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
268b122ec5aSJacob Faibussowitsch   return 0;
269c4762a1bSJed Brown }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown /*TEST
272c4762a1bSJed Brown 
273c4762a1bSJed Brown    test:
274c4762a1bSJed Brown       suffix: a
275c4762a1bSJed Brown       nsize: 2
276c4762a1bSJed Brown       args: -ts_trajectory_type memory -snes_stol 1e-4
277c4762a1bSJed Brown       filter: sort -b
278c4762a1bSJed Brown 
279c4762a1bSJed Brown    test:
280c4762a1bSJed Brown       suffix: b
281c4762a1bSJed Brown       nsize: 2
282c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type arkimex -snes_stol 1e-4
283c4762a1bSJed Brown       filter: sort -b
284c4762a1bSJed Brown 
285c4762a1bSJed Brown    test:
286c4762a1bSJed Brown       suffix: c
287c4762a1bSJed Brown       nsize: 2
288c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4
289c4762a1bSJed Brown       filter: sort -b
290c4762a1bSJed Brown 
291c4762a1bSJed Brown    test:
292c4762a1bSJed Brown       suffix: d
293c4762a1bSJed Brown       nsize: 2
294c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4
295c4762a1bSJed Brown       filter: sort -b
296c4762a1bSJed Brown 
297c4762a1bSJed Brown    test:
298c4762a1bSJed Brown       suffix: e
299c4762a1bSJed Brown       nsize: 2
300c4762a1bSJed Brown       args: -ts_trajectory_type memory -ts_type bdf -ts_adapt_dt_max 0.015 -ts_max_steps 3000
301c4762a1bSJed Brown       filter: sort -b
302c4762a1bSJed Brown 
303c4762a1bSJed Brown    test:
304c4762a1bSJed Brown       suffix: f
305c4762a1bSJed Brown       nsize: 2
306c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 3bs
307c4762a1bSJed Brown       filter: sort -b
308c4762a1bSJed Brown 
309c4762a1bSJed Brown    test:
310c4762a1bSJed Brown       suffix: g
311c4762a1bSJed Brown       nsize: 2
312c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 5bs
313c4762a1bSJed Brown       filter: sort -b
314c4762a1bSJed Brown 
315c4762a1bSJed Brown    test:
316c4762a1bSJed Brown       suffix: h
317c4762a1bSJed Brown       nsize: 2
318c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 6vr
319c4762a1bSJed Brown       filter: sort -b
320c4762a1bSJed Brown       output_file: output/ex41_g.out
321c4762a1bSJed Brown 
322c4762a1bSJed Brown    test:
323c4762a1bSJed Brown       suffix: i
324c4762a1bSJed Brown       nsize: 2
325c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 7vr
326c4762a1bSJed Brown       filter: sort -b
327c4762a1bSJed Brown       output_file: output/ex41_g.out
328c4762a1bSJed Brown 
329c4762a1bSJed Brown    test:
330c4762a1bSJed Brown       suffix: j
331c4762a1bSJed Brown       nsize: 2
332c4762a1bSJed Brown       args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 8vr
333c4762a1bSJed Brown       filter: sort -b
334c4762a1bSJed Brown       output_file: output/ex41_g.out
335c4762a1bSJed Brown 
336c4762a1bSJed Brown TEST*/
337