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