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; 125 J[0][1] = -1.0; 126 J[1][0] = user->mu * (2.0 * u[0] * u[1] + 1.0); 127 J[1][1] = a - user->mu * (1.0 - u[0] * u[0]); 128 129 PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 130 PetscCall(VecRestoreArrayRead(U, &u)); 131 132 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 133 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 134 if (B && A != B) { 135 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 136 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 137 } 138 PetscFunctionReturn(0); 139 } 140 141 static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec U, Mat A, void *ctx) 142 { 143 PetscInt row[] = {0, 1}, col[] = {0}; 144 PetscScalar J[2][1]; 145 const PetscScalar *u; 146 147 PetscFunctionBeginUser; 148 PetscCall(VecGetArrayRead(U, &u)); 149 J[0][0] = 0; 150 J[1][0] = (1. - u[0] * u[0]) * u[1] - u[0]; 151 PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 152 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 153 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 154 PetscCall(VecRestoreArrayRead(U, &u)); 155 PetscFunctionReturn(0); 156 } 157 158 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 159 static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) 160 { 161 const PetscScalar *u; 162 PetscReal tfinal, dt; 163 User user = (User)ctx; 164 Vec interpolatedU; 165 166 PetscFunctionBeginUser; 167 PetscCall(TSGetTimeStep(ts, &dt)); 168 PetscCall(TSGetMaxTime(ts, &tfinal)); 169 170 while (user->next_output <= t && user->next_output <= tfinal) { 171 PetscCall(VecDuplicate(U, &interpolatedU)); 172 PetscCall(TSInterpolate(ts, user->next_output, interpolatedU)); 173 PetscCall(VecGetArrayRead(interpolatedU, &u)); 174 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%g] %" PetscInt_FMT " TS %g (dt = %g) X %g %g\n", (double)user->next_output, step, (double)t, (double)dt, (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; 273 y_ptr[1] = 0.0; 274 PetscCall(VecRestoreArray(user.lambda[0], &y_ptr)); 275 PetscCall(MatCreateVecs(user.A, &user.lambda[1], NULL)); 276 PetscCall(VecGetArray(user.lambda[1], &y_ptr)); 277 y_ptr[0] = 0.0; 278 y_ptr[1] = 1.0; 279 PetscCall(VecRestoreArray(user.lambda[1], &y_ptr)); 280 281 PetscCall(MatCreateVecs(user.Jacp, &user.mup[0], NULL)); 282 PetscCall(VecGetArray(user.mup[0], &x_ptr)); 283 x_ptr[0] = 0.0; 284 PetscCall(VecRestoreArray(user.mup[0], &x_ptr)); 285 PetscCall(MatCreateVecs(user.Jacp, &user.mup[1], NULL)); 286 PetscCall(VecGetArray(user.mup[1], &x_ptr)); 287 x_ptr[0] = 0.0; 288 PetscCall(VecRestoreArray(user.mup[1], &x_ptr)); 289 290 PetscCall(TSSetCostGradients(ts, 2, user.lambda, user.mup)); 291 292 PetscCall(TSAdjointSolve(ts)); 293 294 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: d[y(tf)]/d[y0] d[y(tf)]/d[z0]\n")); 295 PetscCall(VecView(user.lambda[0], PETSC_VIEWER_STDOUT_WORLD)); 296 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: d[z(tf)]/d[y0] d[z(tf)]/d[z0]\n")); 297 PetscCall(VecView(user.lambda[1], PETSC_VIEWER_STDOUT_WORLD)); 298 299 PetscCall(VecGetArray(user.mup[0], &x_ptr)); 300 PetscCall(VecGetArray(user.lambda[0], &y_ptr)); 301 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]; 302 PetscCall(VecRestoreArray(user.mup[0], &x_ptr)); 303 PetscCall(VecRestoreArray(user.lambda[0], &y_ptr)); 304 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n%g\n", (double)PetscRealPart(derp))); 305 306 PetscCall(VecGetArray(user.mup[1], &x_ptr)); 307 PetscCall(VecGetArray(user.lambda[1], &y_ptr)); 308 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]; 309 PetscCall(VecRestoreArray(user.mup[1], &x_ptr)); 310 PetscCall(VecRestoreArray(user.lambda[1], &y_ptr)); 311 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n%g\n", (double)PetscRealPart(derp))); 312 313 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 314 Free work space. All PETSc objects should be destroyed when they 315 are no longer needed. 316 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 317 PetscCall(MatDestroy(&user.A)); 318 PetscCall(MatDestroy(&user.Jacp)); 319 PetscCall(VecDestroy(&user.U)); 320 PetscCall(VecDestroy(&user.lambda[0])); 321 PetscCall(VecDestroy(&user.lambda[1])); 322 PetscCall(VecDestroy(&user.mup[0])); 323 PetscCall(VecDestroy(&user.mup[1])); 324 PetscCall(TSDestroy(&ts)); 325 326 PetscCall(PetscFinalize()); 327 return 0; 328 } 329 330 /*TEST 331 332 test: 333 requires: revolve 334 args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -viewer_binary_skip_info -ts_dt 0.001 -mu 100000 335 336 test: 337 suffix: 2 338 args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 339 340 test: 341 suffix: 3 342 args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 0 343 output_file: output/ex20adj_2.out 344 345 test: 346 suffix: 4 347 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 348 output_file: output/ex20adj_2.out 349 350 test: 351 suffix: 5 352 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 353 output_file: output/ex20adj_2.out 354 355 test: 356 suffix: 6 357 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 358 output_file: output/ex20adj_2.out 359 360 test: 361 suffix: 7 362 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 363 output_file: output/ex20adj_2.out 364 365 test: 366 suffix: 8 367 requires: revolve !cams 368 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 369 output_file: output/ex20adj_3.out 370 371 test: 372 suffix: 9 373 requires: revolve !cams 374 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 375 output_file: output/ex20adj_4.out 376 377 test: 378 requires: revolve 379 suffix: 10 380 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 381 output_file: output/ex20adj_2.out 382 383 test: 384 requires: revolve 385 suffix: 11 386 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 387 output_file: output/ex20adj_2.out 388 389 test: 390 suffix: 12 391 requires: revolve 392 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 393 output_file: output/ex20adj_2.out 394 395 test: 396 suffix: 13 397 requires: revolve 398 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 399 output_file: output/ex20adj_2.out 400 401 test: 402 suffix: 14 403 requires: revolve 404 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 405 output_file: output/ex20adj_2.out 406 407 test: 408 suffix: 15 409 requires: revolve 410 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 411 output_file: output/ex20adj_2.out 412 413 test: 414 suffix: 16 415 requires: revolve 416 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 417 output_file: output/ex20adj_2.out 418 419 test: 420 suffix: 17 421 requires: revolve 422 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 423 output_file: output/ex20adj_2.out 424 425 test: 426 suffix: 18 427 requires: revolve 428 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 429 output_file: output/ex20adj_2.out 430 431 test: 432 suffix: 19 433 requires: revolve 434 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 435 output_file: output/ex20adj_2.out 436 437 test: 438 suffix: 20 439 requires: revolve 440 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 441 output_file: output/ex20adj_2.out 442 443 test: 444 suffix: 21 445 requires: revolve 446 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 447 output_file: output/ex20adj_2.out 448 449 test: 450 suffix: 22 451 args: -ts_type beuler -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 452 output_file: output/ex20adj_2.out 453 454 test: 455 suffix: 23 456 requires: cams 457 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 458 output_file: output/ex20adj_5.out 459 460 test: 461 suffix: 24 462 requires: cams 463 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 464 output_file: output/ex20adj_6.out 465 466 TEST*/ 467