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