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