xref: /petsc/src/ts/tutorials/ex20adj.c (revision 11486bccf1aeea1ca5536228f99d437b39bdaca6)
1 static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n";
2 
3 /* ------------------------------------------------------------------------
4 
5    This program solves the van der Pol DAE ODE equivalent
6       [ u_1' ] = [          u_2                ]  (2)
7       [ u_2' ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]
8    on the domain 0 <= x <= 1, with the boundary conditions
9        u_1(0) = 2, u_2(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
10    and
11        \mu = 10^6 ( y'(0) ~ -0.6666665432100101).,
12    and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with the implicit theta method and its discrete adjoint.
13 
14    Notes:
15    This code demonstrates the TSAdjoint interface to a DAE system.
16 
17    The user provides the implicit right-hand-side function
18    [ F(u',u,t) ] = [u' - f(u,t)] = [ u_1'] - [        u_2             ]
19                                    [ u_2']   [ \mu ((1-u_1^2)u_2-u_1) ]
20 
21    and the Jacobian of F (from the PETSc user manual)
22 
23               dF   dF
24    J(F) = a * -- + --
25               du'  du
26 
27    and the JacobianP of the explicit right-hand side of (2) f(u,t) ( which is equivalent to -F(0,u,t)).
28    df   [       0               ]
29    -- = [                       ]
30    dp   [ (1 - u_1^2) u_2 - u_1 ].
31 
32    See ex20.c for more details on the Jacobian.
33 
34   ------------------------------------------------------------------------- */
35 #include <petscts.h>
36 #include <petsctao.h>
37 
38 typedef struct _n_User *User;
39 struct _n_User {
40   PetscReal mu;
41   PetscReal next_output;
42 
43   /* Sensitivity analysis support */
44   PetscInt  steps;
45   PetscReal ftime;
46   Mat       A;                   /* Jacobian matrix */
47   Mat       Jacp;                /* JacobianP matrix */
48   Vec       U,lambda[2],mup[2];  /* adjoint variables */
49 };
50 
51 /* ----------------------- Explicit form of the ODE  -------------------- */
52 
53 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
54 {
55   User              user = (User)ctx;
56   PetscScalar       *f;
57   const PetscScalar *u;
58 
59   PetscFunctionBeginUser;
60   PetscCall(VecGetArrayRead(U,&u));
61   PetscCall(VecGetArray(F,&f));
62   f[0] = u[1];
63   f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]);
64   PetscCall(VecRestoreArrayRead(U,&u));
65   PetscCall(VecRestoreArray(F,&f));
66   PetscFunctionReturn(0);
67 }
68 
69 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
70 {
71   User              user = (User)ctx;
72   PetscReal         mu   = user->mu;
73   PetscInt          rowcol[] = {0,1};
74   PetscScalar       J[2][2];
75   const PetscScalar *u;
76 
77   PetscFunctionBeginUser;
78   PetscCall(VecGetArrayRead(U,&u));
79   J[0][0] = 0;
80   J[1][0] = -mu*(2.0*u[1]*u[0]+1.);
81   J[0][1] = 1.0;
82   J[1][1] = mu*(1.0-u[0]*u[0]);
83   PetscCall(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
84   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
85   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
86   if (A != B) {
87     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
88     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
89   }
90   PetscCall(VecRestoreArrayRead(U,&u));
91   PetscFunctionReturn(0);
92 }
93 
94 /* ----------------------- Implicit form of the ODE  -------------------- */
95 
96 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
97 {
98   User              user = (User)ctx;
99   const PetscScalar *u,*udot;
100   PetscScalar       *f;
101 
102   PetscFunctionBeginUser;
103   PetscCall(VecGetArrayRead(U,&u));
104   PetscCall(VecGetArrayRead(Udot,&udot));
105   PetscCall(VecGetArray(F,&f));
106   f[0] = udot[0] - u[1];
107   f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]);
108   PetscCall(VecRestoreArrayRead(U,&u));
109   PetscCall(VecRestoreArrayRead(Udot,&udot));
110   PetscCall(VecRestoreArray(F,&f));
111   PetscFunctionReturn(0);
112 }
113 
114 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
115 {
116   User              user     = (User)ctx;
117   PetscInt          rowcol[] = {0,1};
118   PetscScalar       J[2][2];
119   const PetscScalar *u;
120 
121   PetscFunctionBeginUser;
122   PetscCall(VecGetArrayRead(U,&u));
123 
124   J[0][0] = a;     J[0][1] =  -1.0;
125   J[1][0] = user->mu*(2.0*u[0]*u[1] + 1.0);   J[1][1] = a - user->mu*(1.0-u[0]*u[0]);
126 
127   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
128   PetscCall(VecRestoreArrayRead(U,&u));
129 
130   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
131   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
132   if (B && A != B) {
133     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
134     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
135   }
136   PetscFunctionReturn(0);
137 }
138 
139 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
140 {
141   PetscInt          row[] = {0,1},col[]={0};
142   PetscScalar       J[2][1];
143   const PetscScalar *u;
144 
145   PetscFunctionBeginUser;
146   PetscCall(VecGetArrayRead(U,&u));
147   J[0][0] = 0;
148   J[1][0] = (1.-u[0]*u[0])*u[1]-u[0];
149   PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
150   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
151   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
152   PetscCall(VecRestoreArrayRead(U,&u));
153   PetscFunctionReturn(0);
154 }
155 
156 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
157 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx)
158 {
159   const PetscScalar *u;
160   PetscReal         tfinal, dt;
161   User              user = (User)ctx;
162   Vec               interpolatedU;
163 
164   PetscFunctionBeginUser;
165   PetscCall(TSGetTimeStep(ts,&dt));
166   PetscCall(TSGetMaxTime(ts,&tfinal));
167 
168   while (user->next_output <= t && user->next_output <= tfinal) {
169     PetscCall(VecDuplicate(U,&interpolatedU));
170     PetscCall(TSInterpolate(ts,user->next_output,interpolatedU));
171     PetscCall(VecGetArrayRead(interpolatedU,&u));
172     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%g] %" PetscInt_FMT " TS %g (dt = %g) X %g %g\n",
173                           (double)user->next_output,step,(double)t,(double)dt,
174                           (double)PetscRealPart(u[0]),(double)PetscRealPart(u[1])));
175     PetscCall(VecRestoreArrayRead(interpolatedU,&u));
176     PetscCall(VecDestroy(&interpolatedU));
177     user->next_output += 0.1;
178   }
179   PetscFunctionReturn(0);
180 }
181 
182 int main(int argc,char **argv)
183 {
184   TS             ts;
185   PetscBool      monitor = PETSC_FALSE,implicitform = PETSC_TRUE;
186   PetscScalar    *x_ptr,*y_ptr,derp;
187   PetscMPIInt    size;
188   struct _n_User user;
189 
190   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191      Initialize program
192      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193   PetscFunctionBeginUser;
194   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
195   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
196   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
197 
198   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199     Set runtime options
200     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201   user.next_output = 0.0;
202   user.mu          = 1.0e3;
203   user.steps       = 0;
204   user.ftime       = 0.5;
205   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
206   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
207   PetscCall(PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL));
208 
209   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210     Create necessary matrix and vectors, solve same ODE on every process
211     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.A));
213   PetscCall(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
214   PetscCall(MatSetFromOptions(user.A));
215   PetscCall(MatSetUp(user.A));
216   PetscCall(MatCreateVecs(user.A,&user.U,NULL));
217 
218   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
219   PetscCall(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
220   PetscCall(MatSetFromOptions(user.Jacp));
221   PetscCall(MatSetUp(user.Jacp));
222 
223   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224      Create timestepping solver context
225    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
227   PetscCall(TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
228   if (implicitform) {
229     PetscCall(TSSetIFunction(ts,NULL,IFunction,&user));
230     PetscCall(TSSetIJacobian(ts,user.A,user.A,IJacobian,&user));
231     PetscCall(TSSetType(ts,TSCN));
232   } else {
233     PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user));
234     PetscCall(TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user));
235     PetscCall(TSSetType(ts,TSRK));
236   }
237   PetscCall(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user));
238   PetscCall(TSSetMaxTime(ts,user.ftime));
239   PetscCall(TSSetTimeStep(ts,0.001));
240   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
241   if (monitor) PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
242 
243   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244      Set initial conditions
245    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246   PetscCall(VecGetArray(user.U,&x_ptr));
247   x_ptr[0] = 2.0;
248   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
249   PetscCall(VecRestoreArray(user.U,&x_ptr));
250   PetscCall(TSSetTimeStep(ts,0.001));
251 
252   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253     Save trajectory of solution so that TSAdjointSolve() may be used
254    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255   PetscCall(TSSetSaveTrajectory(ts));
256 
257   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258      Set runtime options
259    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260   PetscCall(TSSetFromOptions(ts));
261 
262   PetscCall(TSSolve(ts,user.U));
263   PetscCall(TSGetSolveTime(ts,&user.ftime));
264   PetscCall(TSGetStepNumber(ts,&user.steps));
265 
266   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267      Adjoint model starts here
268      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
269   PetscCall(MatCreateVecs(user.A,&user.lambda[0],NULL));
270   /* Set initial conditions for the adjoint integration */
271   PetscCall(VecGetArray(user.lambda[0],&y_ptr));
272   y_ptr[0] = 1.0; y_ptr[1] = 0.0;
273   PetscCall(VecRestoreArray(user.lambda[0],&y_ptr));
274   PetscCall(MatCreateVecs(user.A,&user.lambda[1],NULL));
275   PetscCall(VecGetArray(user.lambda[1],&y_ptr));
276   y_ptr[0] = 0.0; y_ptr[1] = 1.0;
277   PetscCall(VecRestoreArray(user.lambda[1],&y_ptr));
278 
279   PetscCall(MatCreateVecs(user.Jacp,&user.mup[0],NULL));
280   PetscCall(VecGetArray(user.mup[0],&x_ptr));
281   x_ptr[0] = 0.0;
282   PetscCall(VecRestoreArray(user.mup[0],&x_ptr));
283   PetscCall(MatCreateVecs(user.Jacp,&user.mup[1],NULL));
284   PetscCall(VecGetArray(user.mup[1],&x_ptr));
285   x_ptr[0] = 0.0;
286   PetscCall(VecRestoreArray(user.mup[1],&x_ptr));
287 
288   PetscCall(TSSetCostGradients(ts,2,user.lambda,user.mup));
289 
290   PetscCall(TSAdjointSolve(ts));
291 
292   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[y(tf)]/d[y0]  d[y(tf)]/d[z0]\n"));
293   PetscCall(VecView(user.lambda[0],PETSC_VIEWER_STDOUT_WORLD));
294   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[z(tf)]/d[y0]  d[z(tf)]/d[z0]\n"));
295   PetscCall(VecView(user.lambda[1],PETSC_VIEWER_STDOUT_WORLD));
296 
297   PetscCall(VecGetArray(user.mup[0],&x_ptr));
298   PetscCall(VecGetArray(user.lambda[0],&y_ptr));
299   derp = y_ptr[1]*(-10.0/(81.0*user.mu*user.mu)+2.0*292.0/(2187.0*user.mu*user.mu*user.mu))+x_ptr[0];
300   PetscCall(VecRestoreArray(user.mup[0],&x_ptr));
301   PetscCall(VecRestoreArray(user.lambda[0],&y_ptr));
302   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp)));
303 
304   PetscCall(VecGetArray(user.mup[1],&x_ptr));
305   PetscCall(VecGetArray(user.lambda[1],&y_ptr));
306   derp = y_ptr[1]*(-10.0/(81.0*user.mu*user.mu)+2.0*292.0/(2187.0*user.mu*user.mu*user.mu))+x_ptr[0];
307   PetscCall(VecRestoreArray(user.mup[1],&x_ptr));
308   PetscCall(VecRestoreArray(user.lambda[1],&y_ptr));
309   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp)));
310 
311   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
312      Free work space.  All PETSc objects should be destroyed when they
313      are no longer needed.
314    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
315   PetscCall(MatDestroy(&user.A));
316   PetscCall(MatDestroy(&user.Jacp));
317   PetscCall(VecDestroy(&user.U));
318   PetscCall(VecDestroy(&user.lambda[0]));
319   PetscCall(VecDestroy(&user.lambda[1]));
320   PetscCall(VecDestroy(&user.mup[0]));
321   PetscCall(VecDestroy(&user.mup[1]));
322   PetscCall(TSDestroy(&ts));
323 
324   PetscCall(PetscFinalize());
325   return 0;
326 }
327 
328 /*TEST
329 
330     test:
331       requires: revolve
332       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -viewer_binary_skip_info -ts_dt 0.001 -mu 100000
333 
334     test:
335       suffix: 2
336       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only
337 
338     test:
339       suffix: 3
340       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 0
341       output_file: output/ex20adj_2.out
342 
343     test:
344       suffix: 4
345       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
346       output_file: output/ex20adj_2.out
347 
348     test:
349       suffix: 5
350       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
351       output_file: output/ex20adj_2.out
352 
353     test:
354       suffix: 6
355       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0
356       output_file: output/ex20adj_2.out
357 
358     test:
359       suffix: 7
360       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
361       output_file: output/ex20adj_2.out
362 
363     test:
364       suffix: 8
365       requires: revolve !cams
366       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only -ts_trajectory_monitor
367       output_file: output/ex20adj_3.out
368 
369     test:
370       suffix: 9
371       requires: revolve !cams
372       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only 0 -ts_trajectory_monitor
373       output_file: output/ex20adj_4.out
374 
375     test:
376       requires: revolve
377       suffix: 10
378       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only
379       output_file: output/ex20adj_2.out
380 
381     test:
382       requires: revolve
383       suffix: 11
384       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only 0
385       output_file: output/ex20adj_2.out
386 
387     test:
388       suffix: 12
389       requires: revolve
390       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only
391       output_file: output/ex20adj_2.out
392 
393     test:
394       suffix: 13
395       requires: revolve
396       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0
397       output_file: output/ex20adj_2.out
398 
399     test:
400       suffix: 14
401       requires: revolve
402       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
403       output_file: output/ex20adj_2.out
404 
405     test:
406       suffix: 15
407       requires: revolve
408       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0
409       output_file: output/ex20adj_2.out
410 
411     test:
412       suffix: 16
413       requires: revolve
414       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
415       output_file: output/ex20adj_2.out
416 
417     test:
418       suffix: 17
419       requires: revolve
420       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
421       output_file: output/ex20adj_2.out
422 
423     test:
424       suffix: 18
425       requires: revolve
426       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
427       output_file: output/ex20adj_2.out
428 
429     test:
430       suffix: 19
431       requires: revolve
432       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
433       output_file: output/ex20adj_2.out
434 
435     test:
436       suffix: 20
437       requires: revolve
438       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0
439       output_file: output/ex20adj_2.out
440 
441     test:
442       suffix: 21
443       requires: revolve
444       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
445       output_file: output/ex20adj_2.out
446 
447     test:
448       suffix: 22
449       args: -ts_type beuler -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only
450       output_file: output/ex20adj_2.out
451 
452     test:
453       suffix: 23
454       requires: cams
455       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_units_ram 5 -ts_trajectory_solution_only -ts_trajectory_monitor -ts_trajectory_memory_type cams
456       output_file: output/ex20adj_5.out
457 
458     test:
459       suffix: 24
460       requires: cams
461       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_units_ram 5 -ts_trajectory_solution_only 0 -ts_trajectory_monitor -ts_trajectory_memory_type cams
462       output_file: output/ex20adj_6.out
463 
464 TEST*/
465