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