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 PetscErrorCode ierr; 166 const PetscScalar *u; 167 PetscReal tfinal, dt; 168 User user = (User)ctx; 169 Vec interpolatedU; 170 171 PetscFunctionBeginUser; 172 CHKERRQ(TSGetTimeStep(ts,&dt)); 173 CHKERRQ(TSGetMaxTime(ts,&tfinal)); 174 175 while (user->next_output <= t && user->next_output <= tfinal) { 176 CHKERRQ(VecDuplicate(U,&interpolatedU)); 177 CHKERRQ(TSInterpolate(ts,user->next_output,interpolatedU)); 178 CHKERRQ(VecGetArrayRead(interpolatedU,&u)); 179 ierr = PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n", 180 (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(u[0]), 181 (double)PetscRealPart(u[1]));CHKERRQ(ierr); 182 CHKERRQ(VecRestoreArrayRead(interpolatedU,&u)); 183 CHKERRQ(VecDestroy(&interpolatedU)); 184 user->next_output += 0.1; 185 } 186 PetscFunctionReturn(0); 187 } 188 189 int main(int argc,char **argv) 190 { 191 TS ts; 192 PetscBool monitor = PETSC_FALSE,implicitform = PETSC_TRUE; 193 PetscScalar *x_ptr,*y_ptr,derp; 194 PetscMPIInt size; 195 struct _n_User user; 196 PetscErrorCode ierr; 197 198 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199 Initialize program 200 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 201 ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 202 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 203 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 204 205 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206 Set runtime options 207 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 208 user.next_output = 0.0; 209 user.mu = 1.0e3; 210 user.steps = 0; 211 user.ftime = 0.5; 212 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 213 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 214 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL)); 215 216 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 217 Create necessary matrix and vectors, solve same ODE on every process 218 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 219 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.A)); 220 CHKERRQ(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 221 CHKERRQ(MatSetFromOptions(user.A)); 222 CHKERRQ(MatSetUp(user.A)); 223 CHKERRQ(MatCreateVecs(user.A,&user.U,NULL)); 224 225 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jacp)); 226 CHKERRQ(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1)); 227 CHKERRQ(MatSetFromOptions(user.Jacp)); 228 CHKERRQ(MatSetUp(user.Jacp)); 229 230 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 231 Create timestepping solver context 232 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 233 CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 234 CHKERRQ(TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 235 if (implicitform) { 236 CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&user)); 237 CHKERRQ(TSSetIJacobian(ts,user.A,user.A,IJacobian,&user)); 238 CHKERRQ(TSSetType(ts,TSCN)); 239 } else { 240 CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); 241 CHKERRQ(TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user)); 242 CHKERRQ(TSSetType(ts,TSRK)); 243 } 244 CHKERRQ(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user)); 245 CHKERRQ(TSSetMaxTime(ts,user.ftime)); 246 CHKERRQ(TSSetTimeStep(ts,0.001)); 247 CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 248 if (monitor) { 249 CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL)); 250 } 251 252 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 253 Set initial conditions 254 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 255 CHKERRQ(VecGetArray(user.U,&x_ptr)); 256 x_ptr[0] = 2.0; 257 x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu); 258 CHKERRQ(VecRestoreArray(user.U,&x_ptr)); 259 CHKERRQ(TSSetTimeStep(ts,0.001)); 260 261 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 262 Save trajectory of solution so that TSAdjointSolve() may be used 263 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 264 CHKERRQ(TSSetSaveTrajectory(ts)); 265 266 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 267 Set runtime options 268 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 269 CHKERRQ(TSSetFromOptions(ts)); 270 271 CHKERRQ(TSSolve(ts,user.U)); 272 CHKERRQ(TSGetSolveTime(ts,&user.ftime)); 273 CHKERRQ(TSGetStepNumber(ts,&user.steps)); 274 275 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 276 Adjoint model starts here 277 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 278 CHKERRQ(MatCreateVecs(user.A,&user.lambda[0],NULL)); 279 /* Set initial conditions for the adjoint integration */ 280 CHKERRQ(VecGetArray(user.lambda[0],&y_ptr)); 281 y_ptr[0] = 1.0; y_ptr[1] = 0.0; 282 CHKERRQ(VecRestoreArray(user.lambda[0],&y_ptr)); 283 CHKERRQ(MatCreateVecs(user.A,&user.lambda[1],NULL)); 284 CHKERRQ(VecGetArray(user.lambda[1],&y_ptr)); 285 y_ptr[0] = 0.0; y_ptr[1] = 1.0; 286 CHKERRQ(VecRestoreArray(user.lambda[1],&y_ptr)); 287 288 CHKERRQ(MatCreateVecs(user.Jacp,&user.mup[0],NULL)); 289 CHKERRQ(VecGetArray(user.mup[0],&x_ptr)); 290 x_ptr[0] = 0.0; 291 CHKERRQ(VecRestoreArray(user.mup[0],&x_ptr)); 292 CHKERRQ(MatCreateVecs(user.Jacp,&user.mup[1],NULL)); 293 CHKERRQ(VecGetArray(user.mup[1],&x_ptr)); 294 x_ptr[0] = 0.0; 295 CHKERRQ(VecRestoreArray(user.mup[1],&x_ptr)); 296 297 CHKERRQ(TSSetCostGradients(ts,2,user.lambda,user.mup)); 298 299 CHKERRQ(TSAdjointSolve(ts)); 300 301 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[y(tf)]/d[y0] d[y(tf)]/d[z0]\n")); 302 CHKERRQ(VecView(user.lambda[0],PETSC_VIEWER_STDOUT_WORLD)); 303 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[z(tf)]/d[y0] d[z(tf)]/d[z0]\n")); 304 CHKERRQ(VecView(user.lambda[1],PETSC_VIEWER_STDOUT_WORLD)); 305 306 CHKERRQ(VecGetArray(user.mup[0],&x_ptr)); 307 CHKERRQ(VecGetArray(user.lambda[0],&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 CHKERRQ(VecRestoreArray(user.mup[0],&x_ptr)); 310 CHKERRQ(VecRestoreArray(user.lambda[0],&y_ptr)); 311 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp))); 312 313 CHKERRQ(VecGetArray(user.mup[1],&x_ptr)); 314 CHKERRQ(VecGetArray(user.lambda[1],&y_ptr)); 315 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]; 316 CHKERRQ(VecRestoreArray(user.mup[1],&x_ptr)); 317 CHKERRQ(VecRestoreArray(user.lambda[1],&y_ptr)); 318 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp))); 319 320 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 321 Free work space. All PETSc objects should be destroyed when they 322 are no longer needed. 323 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 324 CHKERRQ(MatDestroy(&user.A)); 325 CHKERRQ(MatDestroy(&user.Jacp)); 326 CHKERRQ(VecDestroy(&user.U)); 327 CHKERRQ(VecDestroy(&user.lambda[0])); 328 CHKERRQ(VecDestroy(&user.lambda[1])); 329 CHKERRQ(VecDestroy(&user.mup[0])); 330 CHKERRQ(VecDestroy(&user.mup[1])); 331 CHKERRQ(TSDestroy(&ts)); 332 333 ierr = PetscFinalize(); 334 return(ierr); 335 } 336 337 /*TEST 338 339 test: 340 requires: revolve 341 args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -viewer_binary_skip_info -ts_dt 0.001 -mu 100000 342 343 test: 344 suffix: 2 345 args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 346 347 test: 348 suffix: 3 349 args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 0 350 output_file: output/ex20adj_2.out 351 352 test: 353 suffix: 4 354 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 355 output_file: output/ex20adj_2.out 356 357 test: 358 suffix: 5 359 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 360 output_file: output/ex20adj_2.out 361 362 test: 363 suffix: 6 364 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 365 output_file: output/ex20adj_2.out 366 367 test: 368 suffix: 7 369 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 370 output_file: output/ex20adj_2.out 371 372 test: 373 suffix: 8 374 requires: revolve !cams 375 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 376 output_file: output/ex20adj_3.out 377 378 test: 379 suffix: 9 380 requires: revolve !cams 381 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 382 output_file: output/ex20adj_4.out 383 384 test: 385 requires: revolve 386 suffix: 10 387 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 388 output_file: output/ex20adj_2.out 389 390 test: 391 requires: revolve 392 suffix: 11 393 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 394 output_file: output/ex20adj_2.out 395 396 test: 397 suffix: 12 398 requires: revolve 399 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 400 output_file: output/ex20adj_2.out 401 402 test: 403 suffix: 13 404 requires: revolve 405 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 406 output_file: output/ex20adj_2.out 407 408 test: 409 suffix: 14 410 requires: revolve 411 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 412 output_file: output/ex20adj_2.out 413 414 test: 415 suffix: 15 416 requires: revolve 417 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 418 output_file: output/ex20adj_2.out 419 420 test: 421 suffix: 16 422 requires: revolve 423 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 424 output_file: output/ex20adj_2.out 425 426 test: 427 suffix: 17 428 requires: revolve 429 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 430 output_file: output/ex20adj_2.out 431 432 test: 433 suffix: 18 434 requires: revolve 435 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 436 output_file: output/ex20adj_2.out 437 438 test: 439 suffix: 19 440 requires: revolve 441 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 442 output_file: output/ex20adj_2.out 443 444 test: 445 suffix: 20 446 requires: revolve 447 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 448 output_file: output/ex20adj_2.out 449 450 test: 451 suffix: 21 452 requires: revolve 453 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 454 output_file: output/ex20adj_2.out 455 456 test: 457 suffix: 22 458 args: -ts_type beuler -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 459 output_file: output/ex20adj_2.out 460 461 test: 462 suffix: 23 463 requires: cams 464 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 465 output_file: output/ex20adj_5.out 466 467 test: 468 suffix: 24 469 requires: cams 470 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 471 output_file: output/ex20adj_6.out 472 473 TEST*/ 474