1c4762a1bSJed Brown static char help[] = "Parallel bouncing ball example to test TS event feature.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown The dynamics of the bouncing ball is described by the ODE 5c4762a1bSJed Brown u1_t = u2 6c4762a1bSJed Brown u2_t = -9.8 7c4762a1bSJed Brown 8c4762a1bSJed Brown Each processor is assigned one ball. 9c4762a1bSJed Brown 10c4762a1bSJed Brown The event function routine checks for the ball hitting the 11c4762a1bSJed Brown ground (u1 = 0). Every time the ball hits the ground, its velocity u2 is attenuated by 12c4762a1bSJed Brown a factor of 0.9 and its height set to 1.0*rank. 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscts.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx) 18c4762a1bSJed Brown { 19c4762a1bSJed Brown PetscErrorCode ierr; 20c4762a1bSJed Brown const PetscScalar *u; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBegin; 23c4762a1bSJed Brown /* Event for ball height */ 24c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 25c4762a1bSJed Brown fvalue[0] = u[0]; 26c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 27c4762a1bSJed Brown PetscFunctionReturn(0); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown 30c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 31c4762a1bSJed Brown { 32c4762a1bSJed Brown PetscErrorCode ierr; 33c4762a1bSJed Brown PetscScalar *u; 34c4762a1bSJed Brown PetscMPIInt rank; 35c4762a1bSJed Brown 36c4762a1bSJed Brown PetscFunctionBegin; 37ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 38c4762a1bSJed Brown if (nevents) { 39c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Ball hit the ground at t = %5.2f seconds -> Processor[%d]\n",(double)t,rank);CHKERRQ(ierr); 40c4762a1bSJed Brown /* Set new initial conditions with .9 attenuation */ 41c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 42c4762a1bSJed Brown u[0] = 1.0*rank; 43c4762a1bSJed Brown u[1] = -0.9*u[1]; 44c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown PetscFunctionReturn(0); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* 50c4762a1bSJed Brown Defines the ODE passed to the ODE solver in explicit form: U_t = F(U) 51c4762a1bSJed Brown */ 52c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown PetscErrorCode ierr; 55c4762a1bSJed Brown PetscScalar *f; 56c4762a1bSJed Brown const PetscScalar *u; 57c4762a1bSJed Brown 58c4762a1bSJed Brown PetscFunctionBegin; 59c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 60c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 62c4762a1bSJed Brown 63c4762a1bSJed Brown f[0] = u[1]; 64c4762a1bSJed Brown f[1] = - 9.8; 65c4762a1bSJed Brown 66c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 68c4762a1bSJed Brown PetscFunctionReturn(0); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning the Jacobian. 73c4762a1bSJed Brown */ 74c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 75c4762a1bSJed Brown { 76c4762a1bSJed Brown PetscErrorCode ierr; 77c4762a1bSJed Brown PetscInt rowcol[2],rstart; 78c4762a1bSJed Brown PetscScalar J[2][2]; 79c4762a1bSJed Brown const PetscScalar *u; 80c4762a1bSJed Brown 81c4762a1bSJed Brown PetscFunctionBegin; 82c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 83c4762a1bSJed Brown 84c4762a1bSJed Brown ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 85c4762a1bSJed Brown rowcol[0] = rstart; rowcol[1] = rstart+1; 86c4762a1bSJed Brown 87c4762a1bSJed Brown J[0][0] = 0.0; J[0][1] = 1.0; 88c4762a1bSJed Brown J[1][0] = 0.0; J[1][1] = 0.0; 89c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 90c4762a1bSJed Brown 91c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94c4762a1bSJed Brown if (A != B) { 95c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* 102c4762a1bSJed Brown Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0 103c4762a1bSJed Brown */ 104c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown PetscErrorCode ierr; 107c4762a1bSJed Brown PetscScalar *f; 108c4762a1bSJed Brown const PetscScalar *u,*udot; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBegin; 111c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 112c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown f[0] = udot[0] - u[1]; 117c4762a1bSJed Brown f[1] = udot[1] + 9.8; 118c4762a1bSJed Brown 119c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 122c4762a1bSJed Brown PetscFunctionReturn(0); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 127c4762a1bSJed Brown */ 128c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 129c4762a1bSJed Brown { 130c4762a1bSJed Brown PetscErrorCode ierr; 131c4762a1bSJed Brown PetscInt rowcol[2],rstart; 132c4762a1bSJed Brown PetscScalar J[2][2]; 133c4762a1bSJed Brown const PetscScalar *u,*udot; 134c4762a1bSJed Brown 135c4762a1bSJed Brown PetscFunctionBegin; 136c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 138c4762a1bSJed Brown 139c4762a1bSJed Brown ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr); 140c4762a1bSJed Brown rowcol[0] = rstart; rowcol[1] = rstart+1; 141c4762a1bSJed Brown 142c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 143c4762a1bSJed Brown J[1][0] = 0.0; J[1][1] = a; 144c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 145c4762a1bSJed Brown 146c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 148c4762a1bSJed Brown 149c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151c4762a1bSJed Brown if (A != B) { 152c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown PetscFunctionReturn(0); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown int main(int argc,char **argv) 159c4762a1bSJed Brown { 160c4762a1bSJed Brown TS ts; /* ODE integrator */ 161c4762a1bSJed Brown Vec U; /* solution will be stored here */ 162c4762a1bSJed Brown PetscErrorCode ierr; 163c4762a1bSJed Brown PetscMPIInt rank; 164c4762a1bSJed Brown PetscInt n = 2; 165c4762a1bSJed Brown PetscScalar *u; 166c4762a1bSJed Brown PetscInt direction=-1; 167c4762a1bSJed Brown PetscBool terminate=PETSC_FALSE; 168c4762a1bSJed Brown PetscBool rhs_form=PETSC_FALSE,hist=PETSC_TRUE; 169c4762a1bSJed Brown TSAdapt adapt; 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172c4762a1bSJed Brown Initialize program 173c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 174c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 175ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 178c4762a1bSJed Brown Create timestepping solver context 179c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 180c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr); 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184c4762a1bSJed Brown Set ODE routines 185c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 186c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 187c4762a1bSJed Brown /* Users are advised against the following branching and code duplication. 188c4762a1bSJed Brown For problems without a mass matrix like the one at hand, the RHSFunction 189c4762a1bSJed Brown (and companion RHSJacobian) interface is enough to support both explicit 190c4762a1bSJed Brown and implicit timesteppers. This tutorial example also deals with the 191c4762a1bSJed Brown IFunction/IJacobian interface for demonstration and testing purposes. */ 192c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-rhs-form",&rhs_form,NULL);CHKERRQ(ierr); 193c4762a1bSJed Brown if (rhs_form) { 194c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,NULL);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,NULL);CHKERRQ(ierr); 196c4762a1bSJed Brown } else { 197c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,NULL);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = TSSetIJacobian(ts,NULL,NULL,IJacobian,NULL);CHKERRQ(ierr); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202c4762a1bSJed Brown Set initial conditions 203c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 204c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&U);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = VecSetSizes(U,n,PETSC_DETERMINE);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = VecSetUp(U);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 208c4762a1bSJed Brown u[0] = 1.0*rank; 209c4762a1bSJed Brown u[1] = 20.0; 210c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 214c4762a1bSJed Brown Set solver options 215c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 216c4762a1bSJed Brown ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = TSSetMaxTime(ts,30.0);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.1);CHKERRQ(ierr); 220*a5b23f4aSJose E. Roman /* The adaptive time step controller could take very large timesteps resulting in 221*a5b23f4aSJose E. Roman the same event occurring multiple times in the same interval. A maximum step size 222c4762a1bSJed Brown limit is enforced here to avoid this issue. */ 223c4762a1bSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = TSAdaptSetStepLimits(adapt,0.0,0.5);CHKERRQ(ierr); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* Set direction and terminate flag for the event */ 228c4762a1bSJed Brown ierr = TSSetEventHandler(ts,1,&direction,&terminate,EventFunction,PostEventFunction,NULL);CHKERRQ(ierr); 229c4762a1bSJed Brown 230c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 231c4762a1bSJed Brown 232c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 233c4762a1bSJed Brown Run timestepping solver 234c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 235c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 236c4762a1bSJed Brown 237c4762a1bSJed Brown if (hist) { /* replay following history */ 238c4762a1bSJed Brown TSTrajectory tj; 239c4762a1bSJed Brown PetscReal tf,t0,dt; 240c4762a1bSJed Brown 241c4762a1bSJed Brown ierr = TSGetTime(ts,&tf);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = TSSetMaxTime(ts,tf);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = TSRestartStep(ts);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = TSSetEventHandler(ts,1,&direction,&terminate,EventFunction,PostEventFunction,NULL);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = TSAdaptSetType(adapt,TSADAPTHISTORY);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = TSAdaptHistorySetTrajectory(adapt,tj,PETSC_FALSE);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = TSAdaptHistoryGetStep(adapt,0,&t0,&dt);CHKERRQ(ierr); 253c4762a1bSJed Brown /* this example fails with single (or smaller) precision */ 254c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL__FP16) 255c4762a1bSJed Brown ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = TSAdaptSetStepLimits(adapt,0.0,0.5);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 258c4762a1bSJed Brown #endif 259c4762a1bSJed Brown ierr = TSSetTime(ts,t0);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 261c4762a1bSJed Brown ierr = TSResetTrajectory(ts);CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 263c4762a1bSJed Brown u[0] = 1.0*rank; 264c4762a1bSJed Brown u[1] = 20.0; 265c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 267c4762a1bSJed Brown } 268c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 269c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 270c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 271c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 273c4762a1bSJed Brown 274c4762a1bSJed Brown ierr = PetscFinalize(); 275c4762a1bSJed Brown return ierr; 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 278c4762a1bSJed Brown /*TEST 279c4762a1bSJed Brown 280c4762a1bSJed Brown test: 281c4762a1bSJed Brown suffix: a 282c4762a1bSJed Brown nsize: 2 283c4762a1bSJed Brown args: -ts_trajectory_type memory -snes_stol 1e-4 284c4762a1bSJed Brown filter: sort -b 285c4762a1bSJed Brown 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown suffix: b 288c4762a1bSJed Brown nsize: 2 289c4762a1bSJed Brown args: -ts_trajectory_type memory -ts_type arkimex -snes_stol 1e-4 290c4762a1bSJed Brown filter: sort -b 291c4762a1bSJed Brown 292c4762a1bSJed Brown test: 293c4762a1bSJed Brown suffix: c 294c4762a1bSJed Brown nsize: 2 295c4762a1bSJed Brown args: -ts_trajectory_type memory -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 296c4762a1bSJed Brown filter: sort -b 297c4762a1bSJed Brown 298c4762a1bSJed Brown test: 299c4762a1bSJed Brown suffix: d 300c4762a1bSJed Brown nsize: 2 301c4762a1bSJed Brown args: -ts_trajectory_type memory -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 302c4762a1bSJed Brown filter: sort -b 303c4762a1bSJed Brown 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown suffix: e 306c4762a1bSJed Brown nsize: 2 307c4762a1bSJed Brown args: -ts_trajectory_type memory -ts_type bdf -ts_adapt_dt_max 0.015 -ts_max_steps 3000 308c4762a1bSJed Brown filter: sort -b 309c4762a1bSJed Brown 310c4762a1bSJed Brown test: 311c4762a1bSJed Brown suffix: f 312c4762a1bSJed Brown nsize: 2 313c4762a1bSJed Brown args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 3bs 314c4762a1bSJed Brown filter: sort -b 315c4762a1bSJed Brown 316c4762a1bSJed Brown test: 317c4762a1bSJed Brown suffix: g 318c4762a1bSJed Brown nsize: 2 319c4762a1bSJed Brown args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 5bs 320c4762a1bSJed Brown filter: sort -b 321c4762a1bSJed Brown 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown suffix: h 324c4762a1bSJed Brown nsize: 2 325c4762a1bSJed Brown args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 6vr 326c4762a1bSJed Brown filter: sort -b 327c4762a1bSJed Brown output_file: output/ex41_g.out 328c4762a1bSJed Brown 329c4762a1bSJed Brown test: 330c4762a1bSJed Brown suffix: i 331c4762a1bSJed Brown nsize: 2 332c4762a1bSJed Brown args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 7vr 333c4762a1bSJed Brown filter: sort -b 334c4762a1bSJed Brown output_file: output/ex41_g.out 335c4762a1bSJed Brown 336c4762a1bSJed Brown test: 337c4762a1bSJed Brown suffix: j 338c4762a1bSJed Brown nsize: 2 339c4762a1bSJed Brown args: -ts_trajectory_type memory -rhs-form -ts_type rk -ts_rk_type 8vr 340c4762a1bSJed Brown filter: sort -b 341c4762a1bSJed Brown output_file: output/ex41_g.out 342c4762a1bSJed Brown 343c4762a1bSJed Brown TEST*/ 344