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