xref: /petsc/src/ts/event/tests/ex2.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 two event functions = piecewise-polynomials, zeros = 1 (rank-0), 9 (last rank)\n"
10*ca4445c7SIlya Fursov                      "Options:\n"
11*ca4445c7SIlya Fursov                      "-dir    d : zero-crossing direction for events\n"
12*ca4445c7SIlya Fursov                      "-flg      : additional output in Postevent\n"
13*ca4445c7SIlya Fursov                      "-restart  : flag for TSRestartStep() in PostEvent\n"
14*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"
15*ca4445c7SIlya Fursov                      "            if x == 0, nothing happens\n";
16*ca4445c7SIlya Fursov 
17*ca4445c7SIlya Fursov #define MAX_NFUNC 100  // max event functions per rank
18*ca4445c7SIlya Fursov #define MAX_NEV   5000 // max zero crossings for each rank
19*ca4445c7SIlya Fursov 
20*ca4445c7SIlya Fursov typedef struct {
21*ca4445c7SIlya Fursov   PetscMPIInt rank, size;
22*ca4445c7SIlya Fursov   PetscReal   pi;
23*ca4445c7SIlya Fursov   PetscReal   fvals[MAX_NFUNC]; // helper array for reporting the residuals
24*ca4445c7SIlya Fursov   PetscReal   evres[MAX_NEV];   // times of found zero-crossings
25*ca4445c7SIlya Fursov   PetscInt    evnum[MAX_NEV];   // number of zero-crossings at each time
26*ca4445c7SIlya Fursov   PetscInt    cnt;              // counter
27*ca4445c7SIlya Fursov   PetscBool   flg;              // flag for additional print in PostEvent
28*ca4445c7SIlya Fursov   PetscBool   restart;          // flag for TSRestartStep() in PostEvent
29*ca4445c7SIlya Fursov   PetscReal   dtpost;           // post-event step
30*ca4445c7SIlya Fursov   PetscInt    postcnt;          // counter for PostEvent calls
31*ca4445c7SIlya Fursov } AppCtx;
32*ca4445c7SIlya Fursov 
33*ca4445c7SIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal gval[], void *ctx);
34*ca4445c7SIlya Fursov PetscErrorCode Postevent(TS ts, PetscInt nev_zero, PetscInt evs_zero[], PetscReal t, Vec U, PetscBool fwd, void *ctx);
35*ca4445c7SIlya Fursov 
36*ca4445c7SIlya Fursov int main(int argc, char **argv)
37*ca4445c7SIlya Fursov {
38*ca4445c7SIlya Fursov   TS           ts;
39*ca4445c7SIlya Fursov   Mat          A;
40*ca4445c7SIlya Fursov   Vec          sol;
41*ca4445c7SIlya Fursov   PetscInt     n, dir0, m = 0;
42*ca4445c7SIlya Fursov   PetscInt     dir[MAX_NFUNC], inds[2];
43*ca4445c7SIlya Fursov   PetscBool    term[MAX_NFUNC];
44*ca4445c7SIlya Fursov   PetscScalar *x, vals[4];
45*ca4445c7SIlya Fursov   AppCtx       ctx;
46*ca4445c7SIlya Fursov 
47*ca4445c7SIlya Fursov   PetscFunctionBeginUser;
48*ca4445c7SIlya Fursov   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
49*ca4445c7SIlya Fursov   setbuf(stdout, NULL);
50*ca4445c7SIlya Fursov   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &ctx.rank));
51*ca4445c7SIlya Fursov   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &ctx.size));
52*ca4445c7SIlya Fursov   ctx.pi      = PetscAcosReal(-1.0);
53*ca4445c7SIlya Fursov   ctx.cnt     = 0;
54*ca4445c7SIlya Fursov   ctx.flg     = PETSC_FALSE;
55*ca4445c7SIlya Fursov   ctx.restart = PETSC_FALSE;
56*ca4445c7SIlya Fursov   ctx.dtpost  = 0;
57*ca4445c7SIlya Fursov   ctx.postcnt = 0;
58*ca4445c7SIlya Fursov 
59*ca4445c7SIlya Fursov   // The linear problem has a 2*2 matrix. The matrix is constant
60*ca4445c7SIlya Fursov   if (ctx.rank == 0) m = 2;
61*ca4445c7SIlya Fursov   inds[0] = 0;
62*ca4445c7SIlya Fursov   inds[1] = 1;
63*ca4445c7SIlya Fursov   vals[0] = 0;
64*ca4445c7SIlya Fursov   vals[1] = 0.2;
65*ca4445c7SIlya Fursov   vals[2] = -0.2;
66*ca4445c7SIlya Fursov   vals[3] = 0;
67*ca4445c7SIlya Fursov   PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, m, m, PETSC_DETERMINE, PETSC_DETERMINE, 2, NULL, 0, NULL, &A));
68*ca4445c7SIlya Fursov   PetscCall(MatSetValues(A, m, inds, m, inds, vals, INSERT_VALUES));
69*ca4445c7SIlya Fursov   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
70*ca4445c7SIlya Fursov   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
71*ca4445c7SIlya Fursov   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
72*ca4445c7SIlya Fursov 
73*ca4445c7SIlya Fursov   PetscCall(MatCreateVecs(A, &sol, NULL));
74*ca4445c7SIlya Fursov   PetscCall(VecGetArray(sol, &x));
75*ca4445c7SIlya Fursov   if (ctx.rank == 0) { // initial conditions
76*ca4445c7SIlya Fursov     x[0] = 0;          // sin(0)
77*ca4445c7SIlya Fursov     x[1] = 1;          // cos(0)
78*ca4445c7SIlya Fursov   }
79*ca4445c7SIlya Fursov   PetscCall(VecRestoreArray(sol, &x));
80*ca4445c7SIlya Fursov 
81*ca4445c7SIlya Fursov   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
82*ca4445c7SIlya Fursov   PetscCall(TSSetProblemType(ts, TS_LINEAR));
83*ca4445c7SIlya Fursov 
84*ca4445c7SIlya Fursov   PetscCall(TSSetRHSFunction(ts, NULL, TSComputeRHSFunctionLinear, NULL));
85*ca4445c7SIlya Fursov   PetscCall(TSSetRHSJacobian(ts, A, A, TSComputeRHSJacobianConstant, NULL));
86*ca4445c7SIlya Fursov 
87*ca4445c7SIlya Fursov   PetscCall(TSSetTimeStep(ts, 0.1));
88*ca4445c7SIlya Fursov   PetscCall(TSSetType(ts, TSBEULER));
89*ca4445c7SIlya Fursov   PetscCall(TSSetMaxSteps(ts, 10000));
90*ca4445c7SIlya Fursov   PetscCall(TSSetMaxTime(ts, 10.0));
91*ca4445c7SIlya Fursov   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
92*ca4445c7SIlya Fursov   PetscCall(TSSetFromOptions(ts));
93*ca4445c7SIlya Fursov 
94*ca4445c7SIlya Fursov   // Set the event handling
95*ca4445c7SIlya Fursov   dir0 = 0;
96*ca4445c7SIlya Fursov   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dir", &dir0, NULL));             // desired zero-crossing direction
97*ca4445c7SIlya Fursov   PetscCall(PetscOptionsHasName(NULL, NULL, "-flg", &ctx.flg));               // flag for additional output
98*ca4445c7SIlya Fursov   PetscCall(PetscOptionsGetBool(NULL, NULL, "-restart", &ctx.restart, NULL)); // flag for TSRestartStep()
99*ca4445c7SIlya Fursov   PetscCall(PetscOptionsGetReal(NULL, NULL, "-dtpost", &ctx.dtpost, NULL));   // post-event step
100*ca4445c7SIlya Fursov 
101*ca4445c7SIlya Fursov   n = 0;               // event counter
102*ca4445c7SIlya Fursov   if (ctx.rank == 0) { // first event -- on rank-0
103*ca4445c7SIlya Fursov     dir[n]    = dir0;
104*ca4445c7SIlya Fursov     term[n++] = PETSC_FALSE;
105*ca4445c7SIlya Fursov   }
106*ca4445c7SIlya Fursov   if (ctx.rank == ctx.size - 1) { // second event -- on last rank
107*ca4445c7SIlya Fursov     dir[n]    = dir0;
108*ca4445c7SIlya Fursov     term[n++] = PETSC_FALSE;
109*ca4445c7SIlya Fursov   }
110*ca4445c7SIlya Fursov   PetscCall(TSSetEventHandler(ts, n, dir, term, EventFunction, Postevent, &ctx));
111*ca4445c7SIlya Fursov 
112*ca4445c7SIlya Fursov   // Solution
113*ca4445c7SIlya Fursov   PetscCall(TSSolve(ts, sol));
114*ca4445c7SIlya Fursov 
115*ca4445c7SIlya Fursov   // The 3 columns printed are: [RANK] [num. of events at the given time] [time of event]
116*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]));
117*ca4445c7SIlya Fursov   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
118*ca4445c7SIlya Fursov 
119*ca4445c7SIlya Fursov   PetscCall(MatDestroy(&A));
120*ca4445c7SIlya Fursov   PetscCall(TSDestroy(&ts));
121*ca4445c7SIlya Fursov   PetscCall(VecDestroy(&sol));
122*ca4445c7SIlya Fursov 
123*ca4445c7SIlya Fursov   PetscCall(PetscFinalize());
124*ca4445c7SIlya Fursov   return 0;
125*ca4445c7SIlya Fursov }
126*ca4445c7SIlya Fursov 
127*ca4445c7SIlya Fursov /*
128*ca4445c7SIlya Fursov   User callback for defining the event-functions
129*ca4445c7SIlya Fursov */
130*ca4445c7SIlya Fursov PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal gval[], void *ctx)
131*ca4445c7SIlya Fursov {
132*ca4445c7SIlya Fursov   PetscInt n   = 0;
133*ca4445c7SIlya Fursov   AppCtx  *Ctx = (AppCtx *)ctx;
134*ca4445c7SIlya Fursov 
135*ca4445c7SIlya Fursov   PetscFunctionBeginUser;
136*ca4445c7SIlya Fursov   // for the test purposes, event-functions are defined based on t
137*ca4445c7SIlya Fursov   // first event -- on rank-0
138*ca4445c7SIlya Fursov   if (Ctx->rank == 0) {
139*ca4445c7SIlya Fursov     if (t < 2.0) gval[n++] = 0.5 * (1 - PetscPowReal(t - 2.0, 12));
140*ca4445c7SIlya Fursov     else gval[n++] = 0.5;
141*ca4445c7SIlya Fursov   }
142*ca4445c7SIlya Fursov 
143*ca4445c7SIlya Fursov   // second event -- on last rank
144*ca4445c7SIlya Fursov   if (Ctx->rank == Ctx->size - 1) {
145*ca4445c7SIlya Fursov     if (t > 8.0) gval[n++] = 0.25 * (1 - PetscPowReal(t - 8.0, 12));
146*ca4445c7SIlya Fursov     else gval[n++] = 0.25;
147*ca4445c7SIlya Fursov   }
148*ca4445c7SIlya Fursov   PetscFunctionReturn(PETSC_SUCCESS);
149*ca4445c7SIlya Fursov }
150*ca4445c7SIlya Fursov 
151*ca4445c7SIlya Fursov /*
152*ca4445c7SIlya Fursov   User callback for the post-event stuff
153*ca4445c7SIlya Fursov */
154*ca4445c7SIlya Fursov PetscErrorCode Postevent(TS ts, PetscInt nev_zero, PetscInt evs_zero[], PetscReal t, Vec U, PetscBool fwd, void *ctx)
155*ca4445c7SIlya Fursov {
156*ca4445c7SIlya Fursov   AppCtx *Ctx = (AppCtx *)ctx;
157*ca4445c7SIlya Fursov 
158*ca4445c7SIlya Fursov   PetscFunctionBeginUser;
159*ca4445c7SIlya Fursov   if (Ctx->flg) {
160*ca4445c7SIlya Fursov     PetscCallBack("EventFunction", EventFunction(ts, t, U, Ctx->fvals, ctx));
161*ca4445c7SIlya Fursov     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] At t = %20.16g : %" PetscInt_FMT " events triggered, fvalues =", Ctx->rank, (double)t, nev_zero));
162*ca4445c7SIlya Fursov     for (PetscInt j = 0; j < nev_zero; j++) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "\t%g", (double)Ctx->fvals[evs_zero[j]]));
163*ca4445c7SIlya Fursov     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "\n"));
164*ca4445c7SIlya Fursov     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
165*ca4445c7SIlya Fursov   }
166*ca4445c7SIlya Fursov 
167*ca4445c7SIlya Fursov   if (Ctx->cnt < MAX_NEV && nev_zero > 0) {
168*ca4445c7SIlya Fursov     Ctx->evres[Ctx->cnt]   = t;
169*ca4445c7SIlya Fursov     Ctx->evnum[Ctx->cnt++] = nev_zero;
170*ca4445c7SIlya Fursov   }
171*ca4445c7SIlya Fursov 
172*ca4445c7SIlya Fursov #ifdef NEW_VERSION
173*ca4445c7SIlya Fursov   Ctx->postcnt++; // sync
174*ca4445c7SIlya Fursov   if (Ctx->dtpost > 0) {
175*ca4445c7SIlya Fursov     if (Ctx->postcnt % 2 == 0) PetscCall(TSSetPostEventStep(ts, Ctx->dtpost));
176*ca4445c7SIlya Fursov     else PetscCall(TSSetPostEventStep(ts, 0));
177*ca4445c7SIlya Fursov   }
178*ca4445c7SIlya Fursov #endif
179*ca4445c7SIlya Fursov 
180*ca4445c7SIlya Fursov   if (Ctx->restart) PetscCall(TSRestartStep(ts));
181*ca4445c7SIlya Fursov   PetscFunctionReturn(PETSC_SUCCESS);
182*ca4445c7SIlya Fursov }
183*ca4445c7SIlya Fursov /*---------------------------------------------------------------------------------------------*/
184*ca4445c7SIlya Fursov /*TEST
185*ca4445c7SIlya Fursov   test:
186*ca4445c7SIlya Fursov     suffix: 0s1
187*ca4445c7SIlya Fursov     output_file: output/ex2_0s1.out
188*ca4445c7SIlya Fursov     args: -dir 0
189*ca4445c7SIlya Fursov     args: -restart {{0 1}}
190*ca4445c7SIlya Fursov     args: -dtpost {{0 0.25}}
191*ca4445c7SIlya Fursov     args: -ts_event_post_event_step {{0 0.31}}
192*ca4445c7SIlya Fursov     args: -ts_type {{beuler rk}}
193*ca4445c7SIlya Fursov     args: -ts_adapt_type {{none basic}}
194*ca4445c7SIlya Fursov     nsize: 1
195*ca4445c7SIlya Fursov 
196*ca4445c7SIlya Fursov   test:
197*ca4445c7SIlya Fursov     suffix: 0s4
198*ca4445c7SIlya Fursov     output_file: output/ex2_0s4.out
199*ca4445c7SIlya Fursov     args: -dir 0
200*ca4445c7SIlya Fursov     args: -restart {{0 1}}
201*ca4445c7SIlya Fursov     args: -dtpost {{0 0.25}}
202*ca4445c7SIlya Fursov     args: -ts_event_post_event_step {{0 0.31}}
203*ca4445c7SIlya Fursov     args: -ts_type {{beuler rk}}
204*ca4445c7SIlya Fursov     args: -ts_adapt_type {{none basic}}
205*ca4445c7SIlya Fursov     nsize: 4
206*ca4445c7SIlya Fursov     filter: sort
207*ca4445c7SIlya Fursov     filter_output: sort
208*ca4445c7SIlya Fursov 
209*ca4445c7SIlya Fursov   test:
210*ca4445c7SIlya Fursov     suffix: pos
211*ca4445c7SIlya Fursov     output_file: output/ex2_pos.out
212*ca4445c7SIlya Fursov     args: -dir 1
213*ca4445c7SIlya Fursov     args: -restart {{0 1}}
214*ca4445c7SIlya Fursov     args: -dtpost {{0 0.25}}
215*ca4445c7SIlya Fursov     args: -ts_event_post_event_step {{0 0.31}}
216*ca4445c7SIlya Fursov     args: -ts_type {{beuler rk}}
217*ca4445c7SIlya Fursov     args: -ts_adapt_type {{none basic}}
218*ca4445c7SIlya Fursov     nsize: {{1 4}}
219*ca4445c7SIlya Fursov     filter: sort
220*ca4445c7SIlya Fursov     filter_output: sort
221*ca4445c7SIlya Fursov 
222*ca4445c7SIlya Fursov   test:
223*ca4445c7SIlya Fursov     suffix: ns1
224*ca4445c7SIlya Fursov     output_file: output/ex2_ns1.out
225*ca4445c7SIlya Fursov     args: -dir -1
226*ca4445c7SIlya Fursov     args: -restart {{0 1}}
227*ca4445c7SIlya Fursov     args: -dtpost {{0 0.25}}
228*ca4445c7SIlya Fursov     args: -ts_event_post_event_step {{0 0.31}}
229*ca4445c7SIlya Fursov     args: -ts_type {{beuler rk}}
230*ca4445c7SIlya Fursov     args: -ts_adapt_type {{none basic}}
231*ca4445c7SIlya Fursov     nsize: 1
232*ca4445c7SIlya Fursov 
233*ca4445c7SIlya Fursov   test:
234*ca4445c7SIlya Fursov     suffix: ns4
235*ca4445c7SIlya Fursov     output_file: output/ex2_ns4.out
236*ca4445c7SIlya Fursov     args: -dir -1
237*ca4445c7SIlya Fursov     args: -restart {{0 1}}
238*ca4445c7SIlya Fursov     args: -dtpost {{0 0.25}}
239*ca4445c7SIlya Fursov     args: -ts_event_post_event_step {{0 0.31}}
240*ca4445c7SIlya Fursov     args: -ts_type {{beuler rk}}
241*ca4445c7SIlya Fursov     args: -ts_adapt_type {{none basic}}
242*ca4445c7SIlya Fursov     nsize: 4
243*ca4445c7SIlya Fursov     filter: sort
244*ca4445c7SIlya Fursov     filter_output: sort
245*ca4445c7SIlya Fursov TEST*/
246