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