1ca4445c7SIlya Fursov #include <petscts.h> 2ca4445c7SIlya Fursov #include <stdio.h> 3ca4445c7SIlya Fursov 4*fe4ad979SIlya Fursov #define NEW_VERSION // Applicable for the new features; avoid this for the older PETSc versions (without TSSetPostEventStep()) 5ca4445c7SIlya Fursov 6ca4445c7SIlya Fursov static char help[] = "Simple linear problem with events\n" 7ca4445c7SIlya Fursov "x_dot = 0.2*y\n" 8ca4445c7SIlya Fursov "y_dot = -0.2*x\n" 9ca4445c7SIlya Fursov "Using 16 event functions:\n" 10ca4445c7SIlya Fursov "7 polynomials (dir=+1) with zeros: 1+2^i, i=-3,...3, on ranks=(i+3)%size\n" 11ca4445c7SIlya Fursov "7 polynomials (dir=-1) with zeros: 1+(8-2^i), i=-3,...3, on ranks=(i+3)%size\n" 12ca4445c7SIlya Fursov "(t-5)^2 * sin(pi*t), with zeros = 1,2,...10, on rank-0\n" 13ca4445c7SIlya Fursov " 0.5 * cos(pi*t), with zeros = 0.5,1.5,...9.5, on last-rank\n" 14ca4445c7SIlya Fursov "Options:\n" 15ca4445c7SIlya Fursov "-dir d : zero-crossing direction for events\n" 16ca4445c7SIlya Fursov "-flg : additional output in Postevent\n" 17*fe4ad979SIlya Fursov "-errtol e : error tolerance, for printing 'pass/fail' for located events (1e-5 by default)\n" 18ca4445c7SIlya Fursov "-restart : flag for TSRestartStep() in PostEvent\n" 19*fe4ad979SIlya Fursov "-dtpost x : if x > 0, then on even PostEvent calls 1st-post-event-step = x is set,\n" 20*fe4ad979SIlya Fursov " on odd PostEvent calls 1st-post-event-step = PETSC_DECIDE is set,\n" 21ca4445c7SIlya Fursov " if x == 0, nothing happens\n"; 22ca4445c7SIlya Fursov 23ca4445c7SIlya Fursov #define MAX_NFUNC 100 // max event functions per rank 24ca4445c7SIlya Fursov #define MAX_NEV 5000 // max zero crossings for each rank 25ca4445c7SIlya Fursov 26ca4445c7SIlya Fursov typedef struct { 27ca4445c7SIlya Fursov PetscMPIInt rank, size; 28ca4445c7SIlya Fursov PetscReal pi; 29ca4445c7SIlya Fursov PetscReal fvals[MAX_NFUNC]; // helper array for reporting the residuals 30ca4445c7SIlya Fursov PetscReal evres[MAX_NEV]; // times of found zero-crossings 31*fe4ad979SIlya Fursov PetscReal ref[MAX_NEV]; // reference times of zero-crossings, for checking 32ca4445c7SIlya Fursov PetscInt cnt; // counter 33*fe4ad979SIlya Fursov PetscInt cntref; // actual length of 'ref' on the given rank 34ca4445c7SIlya Fursov PetscBool flg; // flag for additional print in PostEvent 35*fe4ad979SIlya Fursov PetscReal errtol; // error tolerance, for printing 'pass/fail' for located events (1e-5 by default) 36ca4445c7SIlya Fursov PetscBool restart; // flag for TSRestartStep() in PostEvent 37ca4445c7SIlya Fursov PetscReal dtpost; // post-event step 38ca4445c7SIlya Fursov PetscInt postcnt; // counter for PostEvent calls 39ca4445c7SIlya Fursov PetscReal vtol[MAX_NFUNC]; // vtol array, with extra storage 40ca4445c7SIlya Fursov PetscInt dir0; // desired zero-crossing direction 41ca4445c7SIlya Fursov } AppCtx; 42ca4445c7SIlya Fursov 43ca4445c7SIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal gval[], void *ctx); 44ca4445c7SIlya Fursov PetscErrorCode Postevent(TS ts, PetscInt nev_zero, PetscInt evs_zero[], PetscReal t, Vec U, PetscBool fwd, void *ctx); 45ca4445c7SIlya Fursov static inline void SetVtols(PetscMPIInt rank, PetscMPIInt size, PetscReal tol0, PetscReal tolsin, PetscReal *vtol); // helper function to fill vtol[] 46ca4445c7SIlya Fursov 47ca4445c7SIlya Fursov int main(int argc, char **argv) 48ca4445c7SIlya Fursov { 49ca4445c7SIlya Fursov TS ts; 50ca4445c7SIlya Fursov Mat A; 51ca4445c7SIlya Fursov Vec sol; 52ca4445c7SIlya Fursov PetscInt n, m = 0; 53ca4445c7SIlya Fursov PetscInt dir[MAX_NFUNC], inds[2]; 54ca4445c7SIlya Fursov PetscBool term[MAX_NFUNC]; 55ca4445c7SIlya Fursov PetscScalar *x, vals[4]; 56*fe4ad979SIlya Fursov PetscReal aux; 57ca4445c7SIlya Fursov AppCtx ctx; 58ca4445c7SIlya Fursov 59ca4445c7SIlya Fursov PetscFunctionBeginUser; 60ca4445c7SIlya Fursov PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 61ca4445c7SIlya Fursov setbuf(stdout, NULL); 62ca4445c7SIlya Fursov PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &ctx.rank)); 63ca4445c7SIlya Fursov PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &ctx.size)); 64ca4445c7SIlya Fursov ctx.pi = PetscAcosReal(-1.0); 65ca4445c7SIlya Fursov ctx.cnt = 0; 66*fe4ad979SIlya Fursov ctx.cntref = 0; 67ca4445c7SIlya Fursov ctx.flg = PETSC_FALSE; 68*fe4ad979SIlya Fursov ctx.errtol = 1e-5; 69ca4445c7SIlya Fursov ctx.restart = PETSC_FALSE; 70ca4445c7SIlya Fursov ctx.dtpost = 0; 71ca4445c7SIlya Fursov ctx.postcnt = 0; 72ca4445c7SIlya Fursov 73ca4445c7SIlya Fursov // The linear problem has a 2*2 matrix. The matrix is constant 74ca4445c7SIlya Fursov if (ctx.rank == 0) m = 2; 75ca4445c7SIlya Fursov inds[0] = 0; 76ca4445c7SIlya Fursov inds[1] = 1; 77ca4445c7SIlya Fursov vals[0] = 0; 78ca4445c7SIlya Fursov vals[1] = 0.2; 79ca4445c7SIlya Fursov vals[2] = -0.2; 80ca4445c7SIlya Fursov vals[3] = 0; 81ca4445c7SIlya Fursov PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, m, m, PETSC_DETERMINE, PETSC_DETERMINE, 2, NULL, 0, NULL, &A)); 82ca4445c7SIlya Fursov PetscCall(MatSetValues(A, m, inds, m, inds, vals, INSERT_VALUES)); 83ca4445c7SIlya Fursov PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 84ca4445c7SIlya Fursov PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 85ca4445c7SIlya Fursov PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 86ca4445c7SIlya Fursov 87ca4445c7SIlya Fursov PetscCall(MatCreateVecs(A, &sol, NULL)); 88ca4445c7SIlya Fursov PetscCall(VecGetArray(sol, &x)); 89ca4445c7SIlya Fursov if (ctx.rank == 0) { // initial conditions 90ca4445c7SIlya Fursov x[0] = 0; // sin(0) 91ca4445c7SIlya Fursov x[1] = 1; // cos(0) 92ca4445c7SIlya Fursov } 93ca4445c7SIlya Fursov PetscCall(VecRestoreArray(sol, &x)); 94ca4445c7SIlya Fursov 95ca4445c7SIlya Fursov PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 96ca4445c7SIlya Fursov PetscCall(TSSetProblemType(ts, TS_LINEAR)); 97ca4445c7SIlya Fursov 98ca4445c7SIlya Fursov PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, NULL)); 99ca4445c7SIlya Fursov PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, NULL)); 100ca4445c7SIlya Fursov 101ca4445c7SIlya Fursov PetscCall(TSSetTimeStep(ts, 0.1)); 102ca4445c7SIlya Fursov PetscCall(TSSetType(ts, TSBEULER)); 103ca4445c7SIlya Fursov PetscCall(TSSetMaxSteps(ts, 10000)); 104ca4445c7SIlya Fursov PetscCall(TSSetMaxTime(ts, 10.0)); 105ca4445c7SIlya Fursov PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 106ca4445c7SIlya Fursov PetscCall(TSSetFromOptions(ts)); 107ca4445c7SIlya Fursov 108ca4445c7SIlya Fursov // Set the event handling 109ca4445c7SIlya Fursov ctx.dir0 = 0; 110ca4445c7SIlya Fursov PetscCall(PetscOptionsGetInt(NULL, NULL, "-dir", &ctx.dir0, NULL)); // desired zero-crossing direction 111ca4445c7SIlya Fursov PetscCall(PetscOptionsHasName(NULL, NULL, "-flg", &ctx.flg)); // flag for additional output 112*fe4ad979SIlya Fursov PetscCall(PetscOptionsGetReal(NULL, NULL, "-errtol", &ctx.errtol, NULL)); // error tolerance for located events 113ca4445c7SIlya Fursov PetscCall(PetscOptionsGetBool(NULL, NULL, "-restart", &ctx.restart, NULL)); // flag for TSRestartStep() 114ca4445c7SIlya Fursov PetscCall(PetscOptionsGetReal(NULL, NULL, "-dtpost", &ctx.dtpost, NULL)); // post-event step 115ca4445c7SIlya Fursov 116ca4445c7SIlya Fursov n = 0; // event counter 117*fe4ad979SIlya Fursov aux = 1.0 / 8.0; 118ca4445c7SIlya Fursov for (PetscInt i = -3; i <= 3; i++) { // pos-polynomials 119ca4445c7SIlya Fursov if (ctx.rank == (i + 3) % ctx.size) { 120ca4445c7SIlya Fursov dir[n] = ctx.dir0; 121ca4445c7SIlya Fursov term[n++] = PETSC_FALSE; 122*fe4ad979SIlya Fursov if (ctx.dir0 >= 0) ctx.ref[ctx.cntref++] = 1.0 + aux; 123ca4445c7SIlya Fursov } 124*fe4ad979SIlya Fursov aux *= 2; 125ca4445c7SIlya Fursov } 126*fe4ad979SIlya Fursov aux = 1.0 / 8.0; 127ca4445c7SIlya Fursov for (PetscInt i = -3; i <= 3; i++) { // neg-polynomials 128ca4445c7SIlya Fursov if (ctx.rank == (i + 3) % ctx.size) { 129ca4445c7SIlya Fursov dir[n] = ctx.dir0; 130ca4445c7SIlya Fursov term[n++] = PETSC_FALSE; 131*fe4ad979SIlya Fursov if (ctx.dir0 <= 0) ctx.ref[ctx.cntref++] = 9.0 - aux; 132ca4445c7SIlya Fursov } 133*fe4ad979SIlya Fursov aux *= 2; 134ca4445c7SIlya Fursov } 135ca4445c7SIlya Fursov if (ctx.rank == 0) { // sin-event -- on rank-0 136ca4445c7SIlya Fursov dir[n] = ctx.dir0; 137ca4445c7SIlya Fursov term[n++] = PETSC_FALSE; 138*fe4ad979SIlya Fursov for (PetscInt i = 1; i < MAX_NEV / 2 - 10; i++) { 139*fe4ad979SIlya Fursov if (i % 2 == 1 && ctx.dir0 <= 0) ctx.ref[ctx.cntref++] = i; 140*fe4ad979SIlya Fursov if (i % 2 == 0 && ctx.dir0 >= 0) ctx.ref[ctx.cntref++] = i; 141*fe4ad979SIlya Fursov } 142ca4445c7SIlya Fursov } 143ca4445c7SIlya Fursov if (ctx.rank == ctx.size - 1) { // cos-event -- on last rank 144ca4445c7SIlya Fursov dir[n] = ctx.dir0; 145ca4445c7SIlya Fursov term[n++] = PETSC_FALSE; 146*fe4ad979SIlya Fursov for (PetscInt i = 1; i < MAX_NEV / 2 - 10; i++) { 147*fe4ad979SIlya Fursov if (i % 2 == 1 && ctx.dir0 <= 0) ctx.ref[ctx.cntref++] = i - 0.5; 148*fe4ad979SIlya Fursov if (i % 2 == 0 && ctx.dir0 >= 0) ctx.ref[ctx.cntref++] = i - 0.5; 149ca4445c7SIlya Fursov } 150*fe4ad979SIlya Fursov } 151*fe4ad979SIlya Fursov if (ctx.cntref > 0) PetscCall(PetscSortReal(ctx.cntref, ctx.ref)); 152ca4445c7SIlya Fursov PetscCall(TSSetEventHandler(ts, n, dir, term, EventFunction, Postevent, &ctx)); 153ca4445c7SIlya Fursov SetVtols(ctx.rank, ctx.size, 1e-8, 1e-8, ctx.vtol); 154ca4445c7SIlya Fursov PetscCall(TSSetEventTolerances(ts, PETSC_DECIDE, ctx.vtol)); 155ca4445c7SIlya Fursov 156ca4445c7SIlya Fursov // Solution 157ca4445c7SIlya Fursov PetscCall(TSSolve(ts, sol)); 158ca4445c7SIlya Fursov 159*fe4ad979SIlya Fursov // The 4 columns printed are: [RANK] [time of event] [error w.r.t. reference] ["pass"/"fail"] 160*fe4ad979SIlya Fursov for (PetscInt j = 0; j < ctx.cnt; j++) { 161*fe4ad979SIlya Fursov PetscReal err = 10.0; 162*fe4ad979SIlya Fursov if (j < ctx.cntref) err = PetscAbsReal(ctx.evres[j] - ctx.ref[j]); 163*fe4ad979SIlya Fursov PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%d\t%g\t%g\t%s\n", ctx.rank, (double)ctx.evres[j], (double)err, err < ctx.errtol ? "pass" : "fail")); 164*fe4ad979SIlya Fursov } 165ca4445c7SIlya Fursov PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 166ca4445c7SIlya Fursov 167ca4445c7SIlya Fursov PetscCall(MatDestroy(&A)); 168ca4445c7SIlya Fursov PetscCall(TSDestroy(&ts)); 169ca4445c7SIlya Fursov PetscCall(VecDestroy(&sol)); 170ca4445c7SIlya Fursov 171ca4445c7SIlya Fursov PetscCall(PetscFinalize()); 172ca4445c7SIlya Fursov return 0; 173ca4445c7SIlya Fursov } 174ca4445c7SIlya Fursov 175ca4445c7SIlya Fursov /* 176ca4445c7SIlya Fursov User callback for defining the event-functions 177ca4445c7SIlya Fursov */ 178ca4445c7SIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal gval[], void *ctx) 179ca4445c7SIlya Fursov { 180ca4445c7SIlya Fursov PetscInt n = 0; 181ca4445c7SIlya Fursov AppCtx *Ctx = (AppCtx *)ctx; 182ca4445c7SIlya Fursov PetscReal P; 183ca4445c7SIlya Fursov 184ca4445c7SIlya Fursov PetscFunctionBeginUser; 185ca4445c7SIlya Fursov // for the test purposes, event-functions are defined based on t 186ca4445c7SIlya Fursov for (PetscInt i = -3; i <= 3; i++) { // pos-polynomials 187ca4445c7SIlya Fursov if (Ctx->rank == (i + 3) % Ctx->size) { 188ca4445c7SIlya Fursov P = PetscPowReal(2.0, i); 189ca4445c7SIlya Fursov if (t < 2 + P) gval[n++] = 1 - PetscPowReal(2 + P - t, i + 5); 190ca4445c7SIlya Fursov else gval[n++] = 1; 191ca4445c7SIlya Fursov } 192ca4445c7SIlya Fursov } 193ca4445c7SIlya Fursov for (PetscInt i = -3; i <= 3; i++) { // neg-polynomials 194ca4445c7SIlya Fursov if (Ctx->rank == (i + 3) % Ctx->size) { 195ca4445c7SIlya Fursov P = PetscPowReal(2.0, i); 196ca4445c7SIlya Fursov if (t > 8 - P) gval[n++] = 1 - PetscPowReal(t - 8 + P, i + 5); 197ca4445c7SIlya Fursov else gval[n++] = 1; 198ca4445c7SIlya Fursov } 199ca4445c7SIlya Fursov } 200ca4445c7SIlya Fursov if (Ctx->rank == 0) gval[n++] = (t - 5) * (t - 5) * PetscSinReal(Ctx->pi * t); // sin-event -- on rank-0 201ca4445c7SIlya Fursov if (Ctx->rank == Ctx->size - 1) gval[n++] = 0.5 * PetscCosReal(Ctx->pi * t); // cos-event -- on last rank 202ca4445c7SIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 203ca4445c7SIlya Fursov } 204ca4445c7SIlya Fursov 205ca4445c7SIlya Fursov /* 206ca4445c7SIlya Fursov User callback for the post-event stuff 207ca4445c7SIlya Fursov */ 208ca4445c7SIlya Fursov PetscErrorCode Postevent(TS ts, PetscInt nev_zero, PetscInt evs_zero[], PetscReal t, Vec U, PetscBool fwd, void *ctx) 209ca4445c7SIlya Fursov { 210ca4445c7SIlya Fursov AppCtx *Ctx = (AppCtx *)ctx; 211ca4445c7SIlya Fursov 212ca4445c7SIlya Fursov PetscFunctionBeginUser; 213ca4445c7SIlya Fursov if (Ctx->flg) { 214ca4445c7SIlya Fursov PetscCallBack("EventFunction", EventFunction(ts, t, U, Ctx->fvals, ctx)); 215ca4445c7SIlya Fursov PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] At t = %20.16g : %" PetscInt_FMT " events triggered, fvalues =", Ctx->rank, (double)t, nev_zero)); 216ca4445c7SIlya Fursov for (PetscInt j = 0; j < nev_zero; j++) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "\t%g", (double)Ctx->fvals[evs_zero[j]])); 217ca4445c7SIlya Fursov PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "\n")); 218ca4445c7SIlya Fursov PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 219ca4445c7SIlya Fursov } 220ca4445c7SIlya Fursov 221*fe4ad979SIlya Fursov if (Ctx->cnt + nev_zero < MAX_NEV) 222*fe4ad979SIlya Fursov for (PetscInt i = 0; i < nev_zero; i++) Ctx->evres[Ctx->cnt++] = t; // save the repeating zeros separately for easier/unified testing 223ca4445c7SIlya Fursov 224ca4445c7SIlya Fursov #ifdef NEW_VERSION 225ca4445c7SIlya Fursov Ctx->postcnt++; // sync 226ca4445c7SIlya Fursov if (Ctx->dtpost > 0) { 227ca4445c7SIlya Fursov if (Ctx->postcnt % 2 == 0) PetscCall(TSSetPostEventStep(ts, Ctx->dtpost)); 228*fe4ad979SIlya Fursov else PetscCall(TSSetPostEventStep(ts, PETSC_DECIDE)); 229ca4445c7SIlya Fursov } 230ca4445c7SIlya Fursov #endif 231ca4445c7SIlya Fursov 232ca4445c7SIlya Fursov if ((Ctx->dir0 == 0 && PetscAbsReal(t - 4.0) < 0.01) || (Ctx->dir0 == -1 && PetscAbsReal(t - 3.0) < 0.01)) { 233ca4445c7SIlya Fursov SetVtols(Ctx->rank, Ctx->size, 1e-8, 1e-26, Ctx->vtol); // for better resolution of sin-event at t=5.0 234ca4445c7SIlya Fursov PetscCall(TSSetEventTolerances(ts, PETSC_DECIDE, Ctx->vtol)); 235ca4445c7SIlya Fursov } 236ca4445c7SIlya Fursov if (PetscAbsReal(t - 5.0) < 0.01) { 237ca4445c7SIlya Fursov SetVtols(Ctx->rank, Ctx->size, 1e-8, 1e-8, Ctx->vtol); // back to normal 238ca4445c7SIlya Fursov PetscCall(TSSetEventTolerances(ts, PETSC_DECIDE, Ctx->vtol)); 239ca4445c7SIlya Fursov } 240ca4445c7SIlya Fursov 241ca4445c7SIlya Fursov if (Ctx->restart) PetscCall(TSRestartStep(ts)); 242ca4445c7SIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 243ca4445c7SIlya Fursov } 244ca4445c7SIlya Fursov 245ca4445c7SIlya Fursov // helper function to fill vtol[] 246ca4445c7SIlya Fursov static inline void SetVtols(PetscMPIInt rank, PetscMPIInt size, PetscReal tol0, PetscReal tolsin, PetscReal *vtol) 247ca4445c7SIlya Fursov { 248ca4445c7SIlya Fursov PetscInt n = 0; 249ca4445c7SIlya Fursov for (PetscInt i = -3; i <= 3; i++) 250ca4445c7SIlya Fursov if (rank == (i + 3) % size) vtol[n++] = tol0; // pos-polynomials 251ca4445c7SIlya Fursov for (PetscInt i = -3; i <= 3; i++) 252ca4445c7SIlya Fursov if (rank == (i + 3) % size) vtol[n++] = tol0; // neg-polynomials 253ca4445c7SIlya Fursov if (rank == 0) vtol[n++] = tolsin; // sin-event -- on rank-0 254ca4445c7SIlya Fursov if (rank == size - 1) vtol[n++] = tol0; // cos-event -- on last rank 255ca4445c7SIlya Fursov } 256ca4445c7SIlya Fursov /*---------------------------------------------------------------------------------------------*/ 257*fe4ad979SIlya Fursov /* 258*fe4ad979SIlya Fursov Note, in the tests below, -ts_event_post_event_step is occasionally set to -1, 259*fe4ad979SIlya Fursov which corresponds to PETSC_DECIDE in the API. It is not a very good practice to 260*fe4ad979SIlya Fursov explicitly specify -1 in this option. Rather, if PETSC_DECIDE behaviour is needed, 261*fe4ad979SIlya Fursov simply remove this option altogether. This will result in using the defaults 262*fe4ad979SIlya Fursov (which is PETSC_DECIDE). 263*fe4ad979SIlya Fursov */ 264ca4445c7SIlya Fursov /*TEST 265ca4445c7SIlya Fursov test: 266ca4445c7SIlya Fursov suffix: pos1 267ca4445c7SIlya Fursov output_file: output/ex5_pos1.out 268ca4445c7SIlya Fursov args: -dir 1 -ts_event_dt_min 1e-6 269*fe4ad979SIlya Fursov args: -restart 1 270ca4445c7SIlya Fursov args: -dtpost {{0 0.25}} 271*fe4ad979SIlya Fursov args: -ts_event_post_event_step 0.31 272*fe4ad979SIlya Fursov args: -ts_type rk 273ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 274ca4445c7SIlya Fursov nsize: 1 275ca4445c7SIlya Fursov 276ca4445c7SIlya Fursov test: 277ca4445c7SIlya Fursov suffix: pos4 278ca4445c7SIlya Fursov output_file: output/ex5_pos4.out 279*fe4ad979SIlya Fursov args: -dir 1 -ts_event_dt_min 1e-6 -ts_dt 0.25 280*fe4ad979SIlya Fursov args: -restart 0 281*fe4ad979SIlya Fursov args: -dtpost 0 282*fe4ad979SIlya Fursov args: -ts_event_post_event_step -1 283ca4445c7SIlya Fursov args: -ts_type {{beuler rk}} 284ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 285ca4445c7SIlya Fursov nsize: 4 286ca4445c7SIlya Fursov filter: sort 287ca4445c7SIlya Fursov filter_output: sort 288ca4445c7SIlya Fursov 289ca4445c7SIlya Fursov test: 290ca4445c7SIlya Fursov suffix: neu1 291ca4445c7SIlya Fursov output_file: output/ex5_neu1.out 292ca4445c7SIlya Fursov args: -dir 0 -ts_event_dt_min 1e-6 293*fe4ad979SIlya Fursov args: -restart 0 294ca4445c7SIlya Fursov args: -dtpost {{0 0.25}} 295*fe4ad979SIlya Fursov args: -ts_event_post_event_step -1 296*fe4ad979SIlya Fursov args: -ts_type rk 297ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 298ca4445c7SIlya Fursov nsize: 1 299ca4445c7SIlya Fursov 300ca4445c7SIlya Fursov test: 301ca4445c7SIlya Fursov suffix: neu4 302ca4445c7SIlya Fursov output_file: output/ex5_neu4.out 303*fe4ad979SIlya Fursov args: -dir 0 -ts_event_dt_min 1e-6 -ts_dt 0.25 304*fe4ad979SIlya Fursov args: -dtpost 0 305*fe4ad979SIlya Fursov args: -ts_event_post_event_step {{-1 0.29}} 306*fe4ad979SIlya Fursov args: -ts_event_post_event_second_step {{-1 0.31}} 307*fe4ad979SIlya Fursov args: -ts_type rk 308ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 309ca4445c7SIlya Fursov nsize: 4 310ca4445c7SIlya Fursov filter: sort 311ca4445c7SIlya Fursov filter_output: sort 312ca4445c7SIlya Fursov 313ca4445c7SIlya Fursov test: 314ca4445c7SIlya Fursov suffix: neg2 315ca4445c7SIlya Fursov output_file: output/ex5_neg2.out 316ca4445c7SIlya Fursov args: -dir -1 -ts_event_dt_min 1e-6 317*fe4ad979SIlya Fursov args: -restart 1 318ca4445c7SIlya Fursov args: -dtpost {{0 0.25}} 319*fe4ad979SIlya Fursov args: -ts_event_post_event_step 0.31 320*fe4ad979SIlya Fursov args: -ts_type beuler 321ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 322ca4445c7SIlya Fursov nsize: 2 323ca4445c7SIlya Fursov filter: sort 324ca4445c7SIlya Fursov filter_output: sort 325ca4445c7SIlya Fursov 326ca4445c7SIlya Fursov test: 327ca4445c7SIlya Fursov suffix: neg4 328ca4445c7SIlya Fursov output_file: output/ex5_neg4.out 329*fe4ad979SIlya Fursov args: -dir -1 -ts_event_dt_min 1e-6 -ts_dt 0.25 330*fe4ad979SIlya Fursov args: -restart 0 331*fe4ad979SIlya Fursov args: -dtpost 0 332*fe4ad979SIlya Fursov args: -ts_event_post_event_step -1 333ca4445c7SIlya Fursov args: -ts_type {{beuler rk}} 334ca4445c7SIlya Fursov args: -ts_adapt_type {{none basic}} 335ca4445c7SIlya Fursov nsize: 4 336ca4445c7SIlya Fursov filter: sort 337ca4445c7SIlya Fursov filter_output: sort 338ca4445c7SIlya Fursov 339ca4445c7SIlya Fursov TEST*/ 340