1c4762a1bSJed Brown static char help[] = "Serial 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 There are two events set in this example. The first one checks for the ball hitting the 9c4762a1bSJed Brown ground (u1 = 0). Every time the ball hits the ground, its velocity u2 is attenuated by 10c4762a1bSJed Brown a factor of 0.9. The second event sets a limit on the number of ball bounces. 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscts.h> 14c4762a1bSJed Brown 15c4762a1bSJed Brown typedef struct { 16c4762a1bSJed Brown PetscInt maxbounces; 17c4762a1bSJed Brown PetscInt nbounces; 18c4762a1bSJed Brown } AppCtx; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx) 21c4762a1bSJed Brown { 22c4762a1bSJed Brown AppCtx *app=(AppCtx*)ctx; 23c4762a1bSJed Brown const PetscScalar *u; 24c4762a1bSJed Brown 25c4762a1bSJed Brown PetscFunctionBegin; 26c4762a1bSJed Brown /* Event for ball height */ 279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 28c4762a1bSJed Brown fvalue[0] = u[0]; 29c4762a1bSJed Brown /* Event for number of bounces */ 30c4762a1bSJed Brown fvalue[1] = app->maxbounces - app->nbounces; 319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 32c4762a1bSJed Brown PetscFunctionReturn(0); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 36c4762a1bSJed Brown { 37c4762a1bSJed Brown AppCtx *app=(AppCtx*)ctx; 38c4762a1bSJed Brown PetscScalar *u; 39c4762a1bSJed Brown 40c4762a1bSJed Brown PetscFunctionBegin; 41c4762a1bSJed Brown if (event_list[0] == 0) { 429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Ball hit the ground at t = %5.2f seconds\n",(double)t)); 43c4762a1bSJed Brown /* Set new initial conditions with .9 attenuation */ 449566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 45c4762a1bSJed Brown u[0] = 0.0; 46c4762a1bSJed Brown u[1] = -0.9*u[1]; 479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 48c4762a1bSJed Brown } else if (event_list[0] == 1) { 4963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Ball bounced %" PetscInt_FMT " times\n",app->nbounces)); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown app->nbounces++; 52c4762a1bSJed Brown PetscFunctionReturn(0); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 56c4762a1bSJed Brown Defines the ODE passed to the ODE solver in explicit form: U_t = F(U) 57c4762a1bSJed Brown */ 58c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 59c4762a1bSJed Brown { 60c4762a1bSJed Brown PetscScalar *f; 61c4762a1bSJed Brown const PetscScalar *u; 62c4762a1bSJed Brown 63c4762a1bSJed Brown PetscFunctionBegin; 64c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 669566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown f[0] = u[1]; 69c4762a1bSJed Brown f[1] = - 9.8; 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 73c4762a1bSJed Brown PetscFunctionReturn(0); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of the Jacobian. 78c4762a1bSJed Brown */ 79c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 82c4762a1bSJed Brown PetscScalar J[2][2]; 83c4762a1bSJed Brown const PetscScalar *u; 84c4762a1bSJed Brown 85c4762a1bSJed Brown PetscFunctionBegin; 869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown J[0][0] = 0.0; J[0][1] = 1.0; 89c4762a1bSJed Brown J[1][0] = 0.0; J[1][1] = 0.0; 909566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 96c4762a1bSJed Brown if (A != B) { 979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown PetscFunctionReturn(0); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* 104c4762a1bSJed Brown Defines the ODE passed to the ODE solver in implicit form: F(U_t,U) = 0 105c4762a1bSJed Brown */ 106c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 107c4762a1bSJed Brown { 108c4762a1bSJed Brown PetscScalar *f; 109c4762a1bSJed Brown const PetscScalar *u,*udot; 110c4762a1bSJed Brown 111c4762a1bSJed Brown PetscFunctionBegin; 112c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 1149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 1159566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown f[0] = udot[0] - u[1]; 118c4762a1bSJed Brown f[1] = udot[1] + 9.8; 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 123c4762a1bSJed Brown PetscFunctionReturn(0); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* 127c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 128c4762a1bSJed Brown */ 129c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 130c4762a1bSJed Brown { 131c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 132c4762a1bSJed Brown PetscScalar J[2][2]; 133c4762a1bSJed Brown const PetscScalar *u,*udot; 134c4762a1bSJed Brown 135c4762a1bSJed Brown PetscFunctionBegin; 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 140c4762a1bSJed Brown J[1][0] = 0.0; J[1][1] = a; 1419566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1479566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 148c4762a1bSJed Brown if (A != B) { 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown PetscFunctionReturn(0); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown int main(int argc,char **argv) 156c4762a1bSJed Brown { 157c4762a1bSJed Brown TS ts; /* ODE integrator */ 158c4762a1bSJed Brown Vec U; /* solution will be stored here */ 159c4762a1bSJed Brown PetscMPIInt size; 160c4762a1bSJed Brown PetscInt n = 2; 161c4762a1bSJed Brown PetscScalar *u; 162c4762a1bSJed Brown AppCtx app; 163c4762a1bSJed Brown PetscInt direction[2]; 164c4762a1bSJed Brown PetscBool terminate[2]; 165c4762a1bSJed Brown PetscBool rhs_form=PETSC_FALSE,hist=PETSC_TRUE; 166c4762a1bSJed Brown TSAdapt adapt; 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 169c4762a1bSJed Brown Initialize program 170c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 171*327415f7SBarry Smith PetscFunctionBeginUser; 1729566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 1739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1743c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 175c4762a1bSJed Brown 176c4762a1bSJed Brown app.nbounces = 0; 177c4762a1bSJed Brown app.maxbounces = 10; 178d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex40 options",""); 1799566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-maxbounces","","",app.maxbounces,&app.maxbounces,NULL)); 1809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_adapthistory","","",hist,&hist,NULL)); 181d0609cedSBarry Smith PetscOptionsEnd(); 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184c4762a1bSJed Brown Create timestepping solver context 185c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1869566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 1879566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSROSW)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 190c4762a1bSJed Brown Set ODE routines 191c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1929566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 193c4762a1bSJed Brown /* Users are advised against the following branching and code duplication. 194c4762a1bSJed Brown For problems without a mass matrix like the one at hand, the RHSFunction 195c4762a1bSJed Brown (and companion RHSJacobian) interface is enough to support both explicit 196c4762a1bSJed Brown and implicit timesteppers. This tutorial example also deals with the 197c4762a1bSJed Brown IFunction/IJacobian interface for demonstration and testing purposes. */ 1989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-rhs-form",&rhs_form,NULL)); 199c4762a1bSJed Brown if (rhs_form) { 2009566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,NULL)); 2019566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,NULL)); 202c4762a1bSJed Brown } else { 203c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 2049566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 2059566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 2069566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATDENSE)); 2079566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2089566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 2099566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,IFunction,NULL)); 2109566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,A,A,IJacobian,NULL)); 2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 215c4762a1bSJed Brown Set initial conditions 216c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2179566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&U)); 2189566063dSJacob Faibussowitsch PetscCall(VecSetSizes(U,n,PETSC_DETERMINE)); 2199566063dSJacob Faibussowitsch PetscCall(VecSetUp(U)); 2209566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 221c4762a1bSJed Brown u[0] = 0.0; 222c4762a1bSJed Brown u[1] = 20.0; 2239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 2249566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,U)); 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 227c4762a1bSJed Brown Set solver options 228c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2299566063dSJacob Faibussowitsch if (hist) PetscCall(TSSetSaveTrajectory(ts)); 2309566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,30.0)); 2319566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 2329566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,0.1)); 233a5b23f4aSJose E. Roman /* The adaptive time step controller could take very large timesteps resulting in 234a5b23f4aSJose E. Roman the same event occurring multiple times in the same interval. A maximum step size 235c4762a1bSJed Brown limit is enforced here to avoid this issue. */ 2369566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 2379566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt,0.0,0.5)); 238c4762a1bSJed Brown 239c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 240c4762a1bSJed Brown direction[0] = -1; direction[1] = -1; 241c4762a1bSJed Brown terminate[0] = PETSC_FALSE; terminate[1] = PETSC_TRUE; 2429566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&app)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 245c4762a1bSJed Brown 246c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 247c4762a1bSJed Brown Run timestepping solver 248c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2499566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 250c4762a1bSJed Brown 251c4762a1bSJed Brown if (hist) { /* replay following history */ 252c4762a1bSJed Brown TSTrajectory tj; 253c4762a1bSJed Brown PetscReal tf,t0,dt; 254c4762a1bSJed Brown 255c4762a1bSJed Brown app.nbounces = 0; 2569566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&tf)); 2579566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,tf)); 2589566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts,0)); 2599566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts)); 2609566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 2619566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2629566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 2639566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt,TSADAPTHISTORY)); 2649566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts,&tj)); 2659566063dSJacob Faibussowitsch PetscCall(TSAdaptHistorySetTrajectory(adapt,tj,PETSC_FALSE)); 2669566063dSJacob Faibussowitsch PetscCall(TSAdaptHistoryGetStep(adapt,0,&t0,&dt)); 267c4762a1bSJed Brown /* this example fails with single (or smaller) precision */ 268c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL__FP16) 2699566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt,TSADAPTBASIC)); 2709566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt,0.0,0.5)); 2719566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 272c4762a1bSJed Brown #endif 2739566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,t0)); 2749566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 2759566063dSJacob Faibussowitsch PetscCall(TSResetTrajectory(ts)); 2769566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 277c4762a1bSJed Brown u[0] = 0.0; 278c4762a1bSJed Brown u[1] = 20.0; 2799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 2809566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 281c4762a1bSJed Brown } 282c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 283c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 284c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2869566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 287c4762a1bSJed Brown 2889566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 289b122ec5aSJacob Faibussowitsch return 0; 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown /*TEST 293c4762a1bSJed Brown 294c4762a1bSJed Brown test: 295c4762a1bSJed Brown suffix: a 296c4762a1bSJed Brown args: -snes_stol 1e-4 -ts_trajectory_dirname ex40_a_dir 297c4762a1bSJed Brown output_file: output/ex40.out 298c4762a1bSJed Brown 299c4762a1bSJed Brown test: 300c4762a1bSJed Brown suffix: b 301c4762a1bSJed Brown args: -ts_type arkimex -snes_stol 1e-4 -ts_trajectory_dirname ex40_b_dir 302c4762a1bSJed Brown output_file: output/ex40.out 303c4762a1bSJed Brown 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown suffix: c 306c4762a1bSJed Brown args: -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_c_dir 307c4762a1bSJed Brown output_file: output/ex40.out 308c4762a1bSJed Brown 309c4762a1bSJed Brown test: 310c4762a1bSJed Brown suffix: d 311c4762a1bSJed Brown args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_d_dir 312c4762a1bSJed Brown output_file: output/ex40.out 313c4762a1bSJed Brown 314c4762a1bSJed Brown test: 315c4762a1bSJed Brown suffix: e 316c4762a1bSJed Brown args: -ts_type bdf -ts_adapt_dt_max 0.025 -ts_max_steps 1500 -ts_trajectory_dirname ex40_e_dir 317c4762a1bSJed Brown output_file: output/ex40.out 318c4762a1bSJed Brown 319c4762a1bSJed Brown test: 320c4762a1bSJed Brown suffix: f 321c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 3bs -ts_trajectory_dirname ex40_f_dir 322c4762a1bSJed Brown output_file: output/ex40.out 323c4762a1bSJed Brown 324c4762a1bSJed Brown test: 325c4762a1bSJed Brown suffix: g 326c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 5bs -ts_trajectory_dirname ex40_g_dir 327c4762a1bSJed Brown output_file: output/ex40.out 328c4762a1bSJed Brown 329c4762a1bSJed Brown test: 330c4762a1bSJed Brown suffix: h 331c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 6vr -ts_trajectory_dirname ex40_h_dir 332c4762a1bSJed Brown output_file: output/ex40.out 333c4762a1bSJed Brown 334c4762a1bSJed Brown test: 335c4762a1bSJed Brown suffix: i 336c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 7vr -ts_trajectory_dirname ex40_i_dir 337c4762a1bSJed Brown output_file: output/ex40.out 338c4762a1bSJed Brown 339c4762a1bSJed Brown test: 340c4762a1bSJed Brown suffix: j 341c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 8vr -ts_trajectory_dirname ex40_j_dir 3425385558fSLisandro Dalcin output_file: output/ex40.out 3435385558fSLisandro Dalcin 3445385558fSLisandro Dalcin test: 3455385558fSLisandro Dalcin suffix: k 3465385558fSLisandro Dalcin args: -ts_type theta -ts_adapt_type dsp -ts_trajectory_dirname ex40_k_dir 3475385558fSLisandro Dalcin output_file: output/ex40.out 3485385558fSLisandro Dalcin 3495385558fSLisandro Dalcin test: 3505385558fSLisandro Dalcin suffix: l 3515385558fSLisandro Dalcin args: -rhs-form -ts_type rk -ts_rk_type 2a -ts_trajectory_dirname ex40_l_dir 3525385558fSLisandro Dalcin args: -ts_adapt_type dsp -ts_adapt_always_accept {{false true}} -ts_adapt_dt_min 0.01 353c4762a1bSJed Brown output_file: output/ex40.out 354c4762a1bSJed Brown 35508c4bdbbSLisandro Dalcin test: 35608c4bdbbSLisandro Dalcin suffix: m 35708c4bdbbSLisandro Dalcin args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -test_adapthistory false 35808c4bdbbSLisandro Dalcin args: -ts_max_time 10 -ts_exact_final_time {{STEPOVER MATCHSTEP INTERPOLATE}} 35908c4bdbbSLisandro Dalcin 36008c4bdbbSLisandro Dalcin test: 36108c4bdbbSLisandro Dalcin requires: !single 36208c4bdbbSLisandro Dalcin suffix: n 36308c4bdbbSLisandro Dalcin args: -test_adapthistory false 36408c4bdbbSLisandro Dalcin args: -ts_type alpha -ts_alpha_radius 1.0 -ts_view 36508c4bdbbSLisandro Dalcin args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor 36608c4bdbbSLisandro Dalcin args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false 36708c4bdbbSLisandro Dalcin 368e7685601SHong Zhang test: 369e7685601SHong Zhang requires: !single 370e7685601SHong Zhang suffix: o 371e7685601SHong Zhang args: -rhs-form -ts_type rk -ts_rk_type 2b -ts_trajectory_dirname ex40_o_dir 372e7685601SHong Zhang output_file: output/ex40.out 373c4762a1bSJed Brown TEST*/ 374