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 20*9371c9d4SSatish Balay PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) { 21c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx; 22c4762a1bSJed Brown const PetscScalar *u; 23c4762a1bSJed Brown 247510d9b0SBarry Smith PetscFunctionBeginUser; 25c4762a1bSJed Brown /* Event for ball height */ 269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 27c4762a1bSJed Brown fvalue[0] = u[0]; 28c4762a1bSJed Brown /* Event for number of bounces */ 29c4762a1bSJed Brown fvalue[1] = app->maxbounces - app->nbounces; 309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34*9371c9d4SSatish Balay PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) { 35c4762a1bSJed Brown AppCtx *app = (AppCtx *)ctx; 36c4762a1bSJed Brown PetscScalar *u; 37c4762a1bSJed Brown 387510d9b0SBarry Smith PetscFunctionBeginUser; 39c4762a1bSJed Brown if (event_list[0] == 0) { 409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball hit the ground at t = %5.2f seconds\n", (double)t)); 41c4762a1bSJed Brown /* Set new initial conditions with .9 attenuation */ 429566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 43c4762a1bSJed Brown u[0] = 0.0; 44c4762a1bSJed Brown u[1] = -0.9 * u[1]; 459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 46c4762a1bSJed Brown } else if (event_list[0] == 1) { 4763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Ball bounced %" PetscInt_FMT " times\n", app->nbounces)); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown app->nbounces++; 50c4762a1bSJed Brown PetscFunctionReturn(0); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* 54c4762a1bSJed Brown Defines the ODE passed to the ODE solver in explicit form: U_t = F(U) 55c4762a1bSJed Brown */ 56*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) { 57c4762a1bSJed Brown PetscScalar *f; 58c4762a1bSJed Brown const PetscScalar *u; 59c4762a1bSJed Brown 607510d9b0SBarry Smith PetscFunctionBeginUser; 61c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 639566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown f[0] = u[1]; 66c4762a1bSJed Brown f[1] = -9.8; 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 70c4762a1bSJed Brown PetscFunctionReturn(0); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* 74c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of the Jacobian. 75c4762a1bSJed Brown */ 76*9371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) { 77c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 78c4762a1bSJed Brown PetscScalar J[2][2]; 79c4762a1bSJed Brown const PetscScalar *u; 80c4762a1bSJed Brown 817510d9b0SBarry Smith PetscFunctionBeginUser; 829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 83c4762a1bSJed Brown 84*9371c9d4SSatish Balay J[0][0] = 0.0; 85*9371c9d4SSatish Balay J[0][1] = 1.0; 86*9371c9d4SSatish Balay J[1][0] = 0.0; 87*9371c9d4SSatish Balay J[1][1] = 0.0; 889566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 94c4762a1bSJed Brown if (A != B) { 959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 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 */ 104*9371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) { 105c4762a1bSJed Brown PetscScalar *f; 106c4762a1bSJed Brown const PetscScalar *u, *udot; 107c4762a1bSJed Brown 1087510d9b0SBarry Smith PetscFunctionBeginUser; 109c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown f[0] = udot[0] - u[1]; 115c4762a1bSJed Brown f[1] = udot[1] + 9.8; 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 120c4762a1bSJed Brown PetscFunctionReturn(0); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* 124c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 125c4762a1bSJed Brown */ 126*9371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) { 127c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 128c4762a1bSJed Brown PetscScalar J[2][2]; 129c4762a1bSJed Brown const PetscScalar *u, *udot; 130c4762a1bSJed Brown 1317510d9b0SBarry Smith PetscFunctionBeginUser; 1329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 134c4762a1bSJed Brown 135*9371c9d4SSatish Balay J[0][0] = a; 136*9371c9d4SSatish Balay J[0][1] = -1.0; 137*9371c9d4SSatish Balay J[1][0] = 0.0; 138*9371c9d4SSatish Balay J[1][1] = a; 1399566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 143c4762a1bSJed Brown 1449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 146c4762a1bSJed Brown if (A != B) { 1479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown PetscFunctionReturn(0); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153*9371c9d4SSatish Balay int main(int argc, char **argv) { 154c4762a1bSJed Brown TS ts; /* ODE integrator */ 155c4762a1bSJed Brown Vec U; /* solution will be stored here */ 156c4762a1bSJed Brown PetscMPIInt size; 157c4762a1bSJed Brown PetscInt n = 2; 158c4762a1bSJed Brown PetscScalar *u; 159c4762a1bSJed Brown AppCtx app; 160c4762a1bSJed Brown PetscInt direction[2]; 161c4762a1bSJed Brown PetscBool terminate[2]; 162c4762a1bSJed Brown PetscBool rhs_form = PETSC_FALSE, hist = PETSC_TRUE; 163c4762a1bSJed Brown TSAdapt adapt; 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 166c4762a1bSJed Brown Initialize program 167c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 168327415f7SBarry Smith PetscFunctionBeginUser; 1699566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1713c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 172c4762a1bSJed Brown 173c4762a1bSJed Brown app.nbounces = 0; 174c4762a1bSJed Brown app.maxbounces = 10; 175d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex40 options", ""); 1769566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-maxbounces", "", "", app.maxbounces, &app.maxbounces, NULL)); 1779566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_adapthistory", "", "", hist, &hist, NULL)); 178d0609cedSBarry Smith PetscOptionsEnd(); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Create timestepping solver context 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1839566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1849566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 187c4762a1bSJed Brown Set ODE routines 188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1899566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 190c4762a1bSJed Brown /* Users are advised against the following branching and code duplication. 191c4762a1bSJed Brown For problems without a mass matrix like the one at hand, the RHSFunction 192c4762a1bSJed Brown (and companion RHSJacobian) interface is enough to support both explicit 193c4762a1bSJed Brown and implicit timesteppers. This tutorial example also deals with the 194c4762a1bSJed Brown IFunction/IJacobian interface for demonstration and testing purposes. */ 1959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-rhs-form", &rhs_form, NULL)); 196c4762a1bSJed Brown if (rhs_form) { 1979566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, NULL)); 1989566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, NULL)); 199c4762a1bSJed Brown } else { 200c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 2019566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2029566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 2039566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 2049566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2059566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 2069566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, NULL)); 2079566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, NULL)); 2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 212c4762a1bSJed Brown Set initial conditions 213c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2149566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &U)); 2159566063dSJacob Faibussowitsch PetscCall(VecSetSizes(U, n, PETSC_DETERMINE)); 2169566063dSJacob Faibussowitsch PetscCall(VecSetUp(U)); 2179566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 218c4762a1bSJed Brown u[0] = 0.0; 219c4762a1bSJed Brown u[1] = 20.0; 2209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 2219566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 224c4762a1bSJed Brown Set solver options 225c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2269566063dSJacob Faibussowitsch if (hist) PetscCall(TSSetSaveTrajectory(ts)); 2279566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 30.0)); 2289566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 2299566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.1)); 230a5b23f4aSJose E. Roman /* The adaptive time step controller could take very large timesteps resulting in 231a5b23f4aSJose E. Roman the same event occurring multiple times in the same interval. A maximum step size 232c4762a1bSJed Brown limit is enforced here to avoid this issue. */ 2339566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 2349566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 237*9371c9d4SSatish Balay direction[0] = -1; 238*9371c9d4SSatish Balay direction[1] = -1; 239*9371c9d4SSatish Balay terminate[0] = PETSC_FALSE; 240*9371c9d4SSatish Balay terminate[1] = PETSC_TRUE; 2419566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)&app)); 242c4762a1bSJed Brown 2439566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 244c4762a1bSJed Brown 245c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 246c4762a1bSJed Brown Run timestepping solver 247c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2489566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 249c4762a1bSJed Brown 250c4762a1bSJed Brown if (hist) { /* replay following history */ 251c4762a1bSJed Brown TSTrajectory tj; 252c4762a1bSJed Brown PetscReal tf, t0, dt; 253c4762a1bSJed Brown 254c4762a1bSJed Brown app.nbounces = 0; 2559566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tf)); 2569566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, tf)); 2579566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 2589566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts)); 2599566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2609566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2619566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 2629566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, TSADAPTHISTORY)); 2639566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj)); 2649566063dSJacob Faibussowitsch PetscCall(TSAdaptHistorySetTrajectory(adapt, tj, PETSC_FALSE)); 2659566063dSJacob Faibussowitsch PetscCall(TSAdaptHistoryGetStep(adapt, 0, &t0, &dt)); 266c4762a1bSJed Brown /* this example fails with single (or smaller) precision */ 267c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) || defined(PETSC_USE_REAL__FP16) 2689566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, TSADAPTBASIC)); 2699566063dSJacob Faibussowitsch PetscCall(TSAdaptSetStepLimits(adapt, 0.0, 0.5)); 2709566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 271c4762a1bSJed Brown #endif 2729566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, t0)); 2739566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 2749566063dSJacob Faibussowitsch PetscCall(TSResetTrajectory(ts)); 2759566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 276c4762a1bSJed Brown u[0] = 0.0; 277c4762a1bSJed Brown u[1] = 20.0; 2789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 2799566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 282c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 283c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2859566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 286c4762a1bSJed Brown 2879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 288b122ec5aSJacob Faibussowitsch return 0; 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown /*TEST 292c4762a1bSJed Brown 293c4762a1bSJed Brown test: 294c4762a1bSJed Brown suffix: a 295c4762a1bSJed Brown args: -snes_stol 1e-4 -ts_trajectory_dirname ex40_a_dir 296c4762a1bSJed Brown output_file: output/ex40.out 297c4762a1bSJed Brown 298c4762a1bSJed Brown test: 299c4762a1bSJed Brown suffix: b 300c4762a1bSJed Brown args: -ts_type arkimex -snes_stol 1e-4 -ts_trajectory_dirname ex40_b_dir 301c4762a1bSJed Brown output_file: output/ex40.out 302c4762a1bSJed Brown 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: c 305c4762a1bSJed Brown args: -ts_type theta -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_c_dir 306c4762a1bSJed Brown output_file: output/ex40.out 307c4762a1bSJed Brown 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: d 310c4762a1bSJed Brown args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -ts_trajectory_dirname ex40_d_dir 311c4762a1bSJed Brown output_file: output/ex40.out 312c4762a1bSJed Brown 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown suffix: e 315c4762a1bSJed Brown args: -ts_type bdf -ts_adapt_dt_max 0.025 -ts_max_steps 1500 -ts_trajectory_dirname ex40_e_dir 316c4762a1bSJed Brown output_file: output/ex40.out 317c4762a1bSJed Brown 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: f 320c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 3bs -ts_trajectory_dirname ex40_f_dir 321c4762a1bSJed Brown output_file: output/ex40.out 322c4762a1bSJed Brown 323c4762a1bSJed Brown test: 324c4762a1bSJed Brown suffix: g 325c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 5bs -ts_trajectory_dirname ex40_g_dir 326c4762a1bSJed Brown output_file: output/ex40.out 327c4762a1bSJed Brown 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: h 330c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 6vr -ts_trajectory_dirname ex40_h_dir 331c4762a1bSJed Brown output_file: output/ex40.out 332c4762a1bSJed Brown 333c4762a1bSJed Brown test: 334c4762a1bSJed Brown suffix: i 335c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 7vr -ts_trajectory_dirname ex40_i_dir 336c4762a1bSJed Brown output_file: output/ex40.out 337c4762a1bSJed Brown 338c4762a1bSJed Brown test: 339c4762a1bSJed Brown suffix: j 340c4762a1bSJed Brown args: -rhs-form -ts_type rk -ts_rk_type 8vr -ts_trajectory_dirname ex40_j_dir 3415385558fSLisandro Dalcin output_file: output/ex40.out 3425385558fSLisandro Dalcin 3435385558fSLisandro Dalcin test: 3445385558fSLisandro Dalcin suffix: k 3455385558fSLisandro Dalcin args: -ts_type theta -ts_adapt_type dsp -ts_trajectory_dirname ex40_k_dir 3465385558fSLisandro Dalcin output_file: output/ex40.out 3475385558fSLisandro Dalcin 3485385558fSLisandro Dalcin test: 3495385558fSLisandro Dalcin suffix: l 3505385558fSLisandro Dalcin args: -rhs-form -ts_type rk -ts_rk_type 2a -ts_trajectory_dirname ex40_l_dir 3515385558fSLisandro Dalcin args: -ts_adapt_type dsp -ts_adapt_always_accept {{false true}} -ts_adapt_dt_min 0.01 352c4762a1bSJed Brown output_file: output/ex40.out 353c4762a1bSJed Brown 35408c4bdbbSLisandro Dalcin test: 35508c4bdbbSLisandro Dalcin suffix: m 35608c4bdbbSLisandro Dalcin args: -ts_type alpha -ts_adapt_type basic -ts_atol 1e-1 -snes_stol 1e-4 -test_adapthistory false 35708c4bdbbSLisandro Dalcin args: -ts_max_time 10 -ts_exact_final_time {{STEPOVER MATCHSTEP INTERPOLATE}} 35808c4bdbbSLisandro Dalcin 35908c4bdbbSLisandro Dalcin test: 36008c4bdbbSLisandro Dalcin requires: !single 36108c4bdbbSLisandro Dalcin suffix: n 36208c4bdbbSLisandro Dalcin args: -test_adapthistory false 36308c4bdbbSLisandro Dalcin args: -ts_type alpha -ts_alpha_radius 1.0 -ts_view 36408c4bdbbSLisandro Dalcin args: -ts_dt 0.25 -ts_adapt_type basic -ts_adapt_wnormtype INFINITY -ts_adapt_monitor 36508c4bdbbSLisandro Dalcin args: -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false 36608c4bdbbSLisandro Dalcin 367e7685601SHong Zhang test: 368e7685601SHong Zhang requires: !single 369e7685601SHong Zhang suffix: o 370e7685601SHong Zhang args: -rhs-form -ts_type rk -ts_rk_type 2b -ts_trajectory_dirname ex40_o_dir 371e7685601SHong Zhang output_file: output/ex40.out 372c4762a1bSJed Brown TEST*/ 373