xref: /petsc/src/ts/event/tsevent.c (revision 9d47de495d3c23378050c1b4a410c12a375cb6c6)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
22dc7a7e3SShri Abhyankar 
36fea3669SShri Abhyankar /*
4ca4445c7SIlya Fursov   Fills array sign[] with signs of array f[]. If abs(f[i]) < vtol[i], the zero sign is taken.
5ca4445c7SIlya Fursov   All arrays should have length 'nev'
6ca4445c7SIlya Fursov */
TSEventCalcSigns(PetscInt nev,const PetscReal * f,const PetscReal * vtol,PetscInt * sign)7ca4445c7SIlya Fursov static inline void TSEventCalcSigns(PetscInt nev, const PetscReal *f, const PetscReal *vtol, PetscInt *sign)
8ca4445c7SIlya Fursov {
9ca4445c7SIlya Fursov   for (PetscInt i = 0; i < nev; i++) {
10ca4445c7SIlya Fursov     if (PetscAbsReal(f[i]) < vtol[i]) sign[i] = 0;
11ca4445c7SIlya Fursov     else sign[i] = PetscSign(f[i]);
12ca4445c7SIlya Fursov   }
13ca4445c7SIlya Fursov }
14ca4445c7SIlya Fursov 
15ca4445c7SIlya Fursov /*
166427ac75SLisandro Dalcin   TSEventInitialize - Initializes TSEvent for TSSolve
176fea3669SShri Abhyankar */
TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U)18d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventInitialize(TSEvent event, TS ts, PetscReal t, Vec U)
19d71ae5a4SJacob Faibussowitsch {
206fea3669SShri Abhyankar   PetscFunctionBegin;
213ba16761SJacob Faibussowitsch   if (!event) PetscFunctionReturn(PETSC_SUCCESS);
224f572ea9SToby Isaac   PetscAssertPointer(event, 1);
236427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
246427ac75SLisandro Dalcin   PetscValidHeaderSpecific(U, VEC_CLASSID, 4);
256fea3669SShri Abhyankar   event->ptime_prev    = t;
2638bf2713SShri Abhyankar   event->iterctr       = 0;
27ca4445c7SIlya Fursov   event->processing    = PETSC_FALSE;
28ca4445c7SIlya Fursov   event->revisit_right = PETSC_FALSE;
29ca4445c7SIlya Fursov   PetscCallBack("TSEvent indicator", (*event->indicator)(ts, t, U, event->fvalue_prev, event->ctx));
30ca4445c7SIlya Fursov   TSEventCalcSigns(event->nevents, event->fvalue_prev, event->vtol, event->fsign_prev); // by this time event->vtol should have been defined
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
326fea3669SShri Abhyankar }
336fea3669SShri Abhyankar 
TSEventDestroy(TSEvent * event)34d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventDestroy(TSEvent *event)
35d71ae5a4SJacob Faibussowitsch {
367dbe0728SLisandro Dalcin   PetscFunctionBegin;
374f572ea9SToby Isaac   PetscAssertPointer(event, 1);
383ba16761SJacob Faibussowitsch   if (!*event) PetscFunctionReturn(PETSC_SUCCESS);
399371c9d4SSatish Balay   if (--(*event)->refct > 0) {
409371c9d4SSatish Balay     *event = NULL;
413ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
429371c9d4SSatish Balay   }
43e7069c78SShri 
449566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue_prev));
45ca4445c7SIlya Fursov   PetscCall(PetscFree((*event)->fvalue));
469566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue_right));
47ca4445c7SIlya Fursov   PetscCall(PetscFree((*event)->fsign_prev));
48ca4445c7SIlya Fursov   PetscCall(PetscFree((*event)->fsign));
49ca4445c7SIlya Fursov   PetscCall(PetscFree((*event)->fsign_right));
509566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->side));
51ca4445c7SIlya Fursov   PetscCall(PetscFree((*event)->side_prev));
52ca4445c7SIlya Fursov   PetscCall(PetscFree((*event)->justrefined_AB));
53ca4445c7SIlya Fursov   PetscCall(PetscFree((*event)->gamma_AB));
549566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->direction));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->terminate));
569566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->events_zero));
579566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->vtol));
58a4ffd976SShri Abhyankar 
59ca4445c7SIlya Fursov   for (PetscInt i = 0; i < (*event)->recsize; i++) PetscCall(PetscFree((*event)->recorder.eventidx[i]));
609566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.eventidx));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.nevents));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.stepnum));
639566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.time));
64a4ffd976SShri Abhyankar 
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*event)->monitor));
669566063dSJacob Faibussowitsch   PetscCall(PetscFree(*event));
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
687dbe0728SLisandro Dalcin }
697dbe0728SLisandro Dalcin 
70fe4ad979SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-unknown
71e3005195SShri Abhyankar /*@
72fe4ad979SIlya Fursov   TSSetPostEventStep - Set the first time step to use after the event
73ca4445c7SIlya Fursov 
74ca4445c7SIlya Fursov   Logically Collective
75ca4445c7SIlya Fursov 
76ca4445c7SIlya Fursov   Input Parameters:
77ca4445c7SIlya Fursov + ts  - time integration context
78fe4ad979SIlya Fursov - dt1 - first post event step
79ca4445c7SIlya Fursov 
80ca4445c7SIlya Fursov   Options Database Key:
81fe4ad979SIlya Fursov . -ts_event_post_event_step <dt1> - first time step after the event
82ca4445c7SIlya Fursov 
83ca4445c7SIlya Fursov   Level: advanced
84ca4445c7SIlya Fursov 
85ca4445c7SIlya Fursov   Notes:
86fe4ad979SIlya Fursov   `TSSetPostEventStep()` allows one to set a time step to use immediately following an event.
87fe4ad979SIlya Fursov   Note, if `TSAdapt` is allowed to interfere and reject steps, a large 'dt1' set by `TSSetPostEventStep()` may get truncated,
88fe4ad979SIlya Fursov   resulting in a smaller actual post-event step. See also the warning below regarding the `TSAdapt`.
89ca4445c7SIlya Fursov 
90fe4ad979SIlya Fursov   The post-event time steps should be selected based on the post-event dynamics.
91ca4445c7SIlya Fursov   If the dynamics are stiff, or a significant jump in the equations or the state vector has taken place at the event,
92fe4ad979SIlya Fursov   conservative (small) steps should be employed. If not, then larger time steps may be appropriate.
93ca4445c7SIlya Fursov 
94fe4ad979SIlya Fursov   This function accepts either a numerical value for `dt1`, or `PETSC_DECIDE`. The special value `PETSC_DECIDE` signals the event handler to follow
95fe4ad979SIlya Fursov   the originally planned trajectory, and is assumed by default.
96fe4ad979SIlya Fursov 
97fe4ad979SIlya Fursov   To describe the way `PETSC_DECIDE` affects the post-event steps, consider a trajectory of time points t1 -> t2 -> t3 -> t4.
9826a11704SBarry Smith   Suppose the `TS` has reached and calculated the solution at point t3, and has planned the next move: t3 -> t4.
99fe4ad979SIlya Fursov   At this moment, an event between t2 and t3 is detected, and after a few iterations it is resolved at point `te`, t2 < te < t3.
100fe4ad979SIlya Fursov   After event `te`, two post-event steps can be specified: the first one dt1 (`TSSetPostEventStep()`),
101fe4ad979SIlya Fursov   and the second one dt2 (`TSSetPostEventSecondStep()`). Both post-event steps can be either `PETSC_DECIDE`, or a number.
102fe4ad979SIlya Fursov   Four different combinations are possible\:
103fe4ad979SIlya Fursov 
10426a11704SBarry Smith   1. dt1 = `PETSC_DECIDE`, dt2 = `PETSC_DECIDE`. Then, after `te` `TS` goes to t3, and then to t4. This is the all-default behaviour.
105fe4ad979SIlya Fursov 
10626a11704SBarry Smith   2. dt1 = `PETSC_DECIDE`, dt2 = x2 (numerical). Then, after `te` `TS` goes to t3, and then to t3+x2.
107fe4ad979SIlya Fursov 
10826a11704SBarry Smith   3. dt1 = x1 (numerical), dt2 = x2 (numerical). Then, after `te` `TS` goes to te+x1, and then to te+x1+x2.
109fe4ad979SIlya Fursov 
11026a11704SBarry Smith   4. dt1 = x1 (numerical), dt2 = `PETSC_DECIDE`. Then, after `te` `TS` goes to te+x1, and event handler does not interfere to the subsequent steps.
111fe4ad979SIlya Fursov 
112fe4ad979SIlya Fursov   In the special case when `te` == t3 with a good precision, the post-event step te -> t3 is not performed, so behaviour of (1) and (2) becomes\:
113fe4ad979SIlya Fursov 
11426a11704SBarry Smith   1a. After `te` `TS` goes to t4, and event handler does not interfere to the subsequent steps.
115fe4ad979SIlya Fursov 
11626a11704SBarry Smith   2a. After `te` `TS` goes to t4, and then to t4+x2.
117fe4ad979SIlya Fursov 
11826a11704SBarry Smith   Warning! When the second post-event step (either `PETSC_DECIDE` or a numerical value) is managed by the event handler, i.e. in cases 1, 2, 3 and 2a,
119fe4ad979SIlya Fursov   `TSAdapt` will never analyse (and never do a reasonable rejection of) the first post-event step. The first post-event step will always be accepted.
120fe4ad979SIlya Fursov   In this situation, it is the user's responsibility to make sure the step size is appropriate!
121fe4ad979SIlya Fursov   In cases 4 and 1a, however, `TSAdapt` will analyse the first post-event step, and is allowed to reject it.
122ca4445c7SIlya Fursov 
123ca4445c7SIlya Fursov   This function can be called not only in the initial setup, but also inside the `postevent()` callback set with `TSSetEventHandler()`,
124fe4ad979SIlya Fursov   affecting the post-event steps for the current event, and the subsequent ones.
125fe4ad979SIlya Fursov   Thus, the strategy of the post-event time step definition can be adjusted on the fly.
126fe4ad979SIlya Fursov   In case several events are triggered in the given time point, only a single postevent handler is invoked,
127ca4445c7SIlya Fursov   and the user is to determine what post-event time step is more appropriate in this situation.
128ca4445c7SIlya Fursov 
129fe4ad979SIlya Fursov   The default value is `PETSC_DECIDE`.
130ca4445c7SIlya Fursov 
131fe4ad979SIlya Fursov   Developer Notes:
132fe4ad979SIlya Fursov   Event processing starts after visiting point t3, which means ts->adapt->dt_span_cached has been set to whatever value is required
133fe4ad979SIlya Fursov   when planning the step t3 -> t4.
134fe4ad979SIlya Fursov 
135fe4ad979SIlya Fursov .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()`, `TSSetPostEventSecondStep()`
136ca4445c7SIlya Fursov @*/
TSSetPostEventStep(TS ts,PetscReal dt1)137fe4ad979SIlya Fursov PetscErrorCode TSSetPostEventStep(TS ts, PetscReal dt1)
138ca4445c7SIlya Fursov {
139ca4445c7SIlya Fursov   PetscFunctionBegin;
140fe4ad979SIlya Fursov   ts->event->timestep_postevent = dt1;
141ca4445c7SIlya Fursov   PetscFunctionReturn(PETSC_SUCCESS);
142ca4445c7SIlya Fursov }
143ca4445c7SIlya Fursov 
144fe4ad979SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-unknown
145ca4445c7SIlya Fursov /*@
146fe4ad979SIlya Fursov   TSSetPostEventSecondStep - Set the second time step to use after the event
147458122a4SShri Abhyankar 
148458122a4SShri Abhyankar   Logically Collective
149458122a4SShri Abhyankar 
1504165533cSJose E. Roman   Input Parameters:
151458122a4SShri Abhyankar + ts  - time integration context
152fe4ad979SIlya Fursov - dt2 - second post event step
153458122a4SShri Abhyankar 
154bcf0153eSBarry Smith   Options Database Key:
155fe4ad979SIlya Fursov . -ts_event_post_event_second_step <dt2> - second time step after the event
156458122a4SShri Abhyankar 
157bcf0153eSBarry Smith   Level: advanced
158458122a4SShri Abhyankar 
159bcf0153eSBarry Smith   Notes:
160fe4ad979SIlya Fursov   `TSSetPostEventSecondStep()` allows one to set the second time step after the event.
161458122a4SShri Abhyankar 
162fe4ad979SIlya Fursov   The post-event time steps should be selected based on the post-event dynamics.
163fe4ad979SIlya Fursov   If the dynamics are stiff, or a significant jump in the equations or the state vector has taken place at the event,
164fe4ad979SIlya Fursov   conservative (small) steps should be employed. If not, then larger time steps may be appropriate.
165fe4ad979SIlya Fursov 
166fe4ad979SIlya Fursov   This function accepts either a numerical value for `dt2`, or `PETSC_DECIDE` (default).
167fe4ad979SIlya Fursov 
168fe4ad979SIlya Fursov   To describe the way `PETSC_DECIDE` affects the post-event steps, consider a trajectory of time points t1 -> t2 -> t3 -> t4.
16926a11704SBarry Smith   Suppose the `TS` has reached and calculated the solution at point t3, and has planned the next move: t3 -> t4.
170fe4ad979SIlya Fursov   At this moment, an event between t2 and t3 is detected, and after a few iterations it is resolved at point `te`, t2 < te < t3.
171fe4ad979SIlya Fursov   After event `te`, two post-event steps can be specified: the first one dt1 (`TSSetPostEventStep()`),
172fe4ad979SIlya Fursov   and the second one dt2 (`TSSetPostEventSecondStep()`). Both post-event steps can be either `PETSC_DECIDE`, or a number.
173fe4ad979SIlya Fursov   Four different combinations are possible\:
174fe4ad979SIlya Fursov 
17526a11704SBarry Smith   1. dt1 = `PETSC_DECIDE`, dt2 = `PETSC_DECIDE`. Then, after `te` `TS` goes to t3, and then to t4. This is the all-default behaviour.
176fe4ad979SIlya Fursov 
17726a11704SBarry Smith   2. dt1 = `PETSC_DECIDE`, dt2 = x2 (numerical). Then, after `te` `TS` goes to t3, and then to t3+x2.
178fe4ad979SIlya Fursov 
17926a11704SBarry Smith   3. dt1 = x1 (numerical), dt2 = x2 (numerical). Then, after `te` `TS` goes to te+x1, and then to te+x1+x2.
180fe4ad979SIlya Fursov 
18126a11704SBarry Smith   4. dt1 = x1 (numerical), dt2 = `PETSC_DECIDE`. Then, after `te` `TS` goes to te+x1, and event handler does not interfere to the subsequent steps.
182fe4ad979SIlya Fursov 
183fe4ad979SIlya Fursov   In the special case when `te` == t3 with a good precision, the post-event step te -> t3 is not performed, so behaviour of (1) and (2) becomes\:
184fe4ad979SIlya Fursov 
18526a11704SBarry Smith   1a. After `te` `TS` goes to t4, and event handler does not interfere to the subsequent steps.
186fe4ad979SIlya Fursov 
18726a11704SBarry Smith   2a. After `te` `TS` goes to t4, and then to t4+x2.
188fe4ad979SIlya Fursov 
18926a11704SBarry Smith   Warning! When the second post-event step (either `PETSC_DECIDE` or a numerical value) is managed by the event handler, i.e. in cases 1, 2, 3 and 2a,
190fe4ad979SIlya Fursov   `TSAdapt` will never analyse (and never do a reasonable rejection of) the first post-event step. The first post-event step will always be accepted.
191fe4ad979SIlya Fursov   In this situation, it is the user's responsibility to make sure the step size is appropriate!
192fe4ad979SIlya Fursov   In cases 4 and 1a, however, `TSAdapt` will analyse the first post-event step, and is allowed to reject it.
193fe4ad979SIlya Fursov 
194fe4ad979SIlya Fursov   This function can be called not only in the initial setup, but also inside the `postevent()` callback set with `TSSetEventHandler()`,
195fe4ad979SIlya Fursov   affecting the post-event steps for the current event, and the subsequent ones.
196fe4ad979SIlya Fursov 
197fe4ad979SIlya Fursov   The default value is `PETSC_DECIDE`.
198fe4ad979SIlya Fursov 
199fe4ad979SIlya Fursov .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()`, `TSSetPostEventStep()`
200458122a4SShri Abhyankar @*/
TSSetPostEventSecondStep(TS ts,PetscReal dt2)201fe4ad979SIlya Fursov PetscErrorCode TSSetPostEventSecondStep(TS ts, PetscReal dt2)
202d71ae5a4SJacob Faibussowitsch {
203458122a4SShri Abhyankar   PetscFunctionBegin;
204fe4ad979SIlya Fursov   ts->event->timestep_2nd_postevent = dt2;
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206458122a4SShri Abhyankar }
207458122a4SShri Abhyankar 
208458122a4SShri Abhyankar /*@
209ca4445c7SIlya Fursov   TSSetEventTolerances - Set tolerances for event (indicator function) zero crossings
210e3005195SShri Abhyankar 
211e3005195SShri Abhyankar   Logically Collective
212e3005195SShri Abhyankar 
2134165533cSJose E. Roman   Input Parameters:
214e3005195SShri Abhyankar + ts   - time integration context
21509cb0f53SBarry Smith . tol  - tolerance, `PETSC_CURRENT` to leave the current value
2169b3d3759SBarry Smith - vtol - array of tolerances or `NULL`, used in preference to `tol` if present
217e3005195SShri Abhyankar 
218bcf0153eSBarry Smith   Options Database Key:
219ca4445c7SIlya Fursov . -ts_event_tol <tol> - tolerance for event (indicator function) zero crossing
220e3005195SShri Abhyankar 
221e3005195SShri Abhyankar   Level: beginner
222e3005195SShri Abhyankar 
223bcf0153eSBarry Smith   Notes:
224ca4445c7SIlya Fursov   One must call `TSSetEventHandler()` before setting the tolerances.
225bcf0153eSBarry Smith 
226ca4445c7SIlya Fursov   The size of `vtol` should be equal to the number of events on the given process.
22720f4b53cSBarry Smith 
228ca4445c7SIlya Fursov   This function can be also called from the `postevent()` callback set with `TSSetEventHandler()`,
229ca4445c7SIlya Fursov   to adjust the tolerances on the fly.
230bcf0153eSBarry Smith 
2311cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()`
232e3005195SShri Abhyankar @*/
TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])233d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetEventTolerances(TS ts, PetscReal tol, PetscReal vtol[])
234d71ae5a4SJacob Faibussowitsch {
2356427ac75SLisandro Dalcin   TSEvent event;
236e3005195SShri Abhyankar 
237e3005195SShri Abhyankar   PetscFunctionBegin;
2386427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2394f572ea9SToby Isaac   if (vtol) PetscAssertPointer(vtol, 3);
2403c633725SBarry Smith   PetscCheck(ts->event, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set the events first by calling TSSetEventHandler()");
2416427ac75SLisandro Dalcin 
2426427ac75SLisandro Dalcin   event = ts->event;
243e3005195SShri Abhyankar   if (vtol) {
244ca4445c7SIlya Fursov     for (PetscInt i = 0; i < event->nevents; i++) event->vtol[i] = vtol[i];
245e3005195SShri Abhyankar   } else {
24609cb0f53SBarry Smith     if (tol != (PetscReal)PETSC_CURRENT) {
247ca4445c7SIlya Fursov       for (PetscInt i = 0; i < event->nevents; i++) event->vtol[i] = tol;
248e3005195SShri Abhyankar     }
249e3005195SShri Abhyankar   }
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251e3005195SShri Abhyankar }
252e3005195SShri Abhyankar 
2532dc7a7e3SShri Abhyankar /*@C
254ca4445c7SIlya Fursov   TSSetEventHandler - Sets functions and parameters used for indicating events and handling them
2552dc7a7e3SShri Abhyankar 
256c3339decSBarry Smith   Logically Collective
2572dc7a7e3SShri Abhyankar 
2582dc7a7e3SShri Abhyankar   Input Parameters:
259bcf0153eSBarry Smith + ts        - the `TS` context obtained from `TSCreate()`
260ca4445c7SIlya Fursov . nevents   - number of local events (i.e. managed by the given MPI process)
261ca4445c7SIlya Fursov . direction - direction of zero crossing to be detected (one for each local event).
262ca4445c7SIlya Fursov               `-1` => zero crossing in negative direction,
263ca4445c7SIlya Fursov               `+1` => zero crossing in positive direction, `0` => both ways
264ca4445c7SIlya Fursov . terminate - flag to indicate whether time stepping should be terminated after
265ca4445c7SIlya Fursov               an event is detected (one for each local event)
266ca4445c7SIlya Fursov . indicator - callback defininig the user indicator functions whose sign changes (see `direction`) mark presence of the events
267ca4445c7SIlya Fursov . postevent - [optional] user post-event callback; it can change the solution, ODE etc at the time of the event
268ca4445c7SIlya Fursov - ctx       - [optional] user-defined context for private data for the
269ca4445c7SIlya Fursov               `indicator()` and `postevent()` routines (use `NULL` if no
270ca4445c7SIlya Fursov               context is desired)
2712dc7a7e3SShri Abhyankar 
2729b3d3759SBarry Smith   Calling sequence of `indicator`:
2732fe279fdSBarry Smith + ts     - the `TS` context
2742dc7a7e3SShri Abhyankar . t      - current time
275ca4445c7SIlya Fursov . U      - current solution
276ca4445c7SIlya Fursov . fvalue - output array with values of local indicator functions (length == `nevents`) for time t and state-vector U
2779b3d3759SBarry Smith - ctx    - the context passed as the final argument to `TSSetEventHandler()`
2782dc7a7e3SShri Abhyankar 
27920f4b53cSBarry Smith   Calling sequence of `postevent`:
2802fe279fdSBarry Smith + ts           - the `TS` context
281ca4445c7SIlya Fursov . nevents_zero - number of triggered local events (whose indicator function is marked as crossing zero, and direction is appropriate)
282ca4445c7SIlya Fursov . events_zero  - indices of the triggered local events
2832dc7a7e3SShri Abhyankar . t            - current time
2842dc7a7e3SShri Abhyankar . U            - current solution
285ca4445c7SIlya Fursov . forwardsolve - flag to indicate whether `TS` is doing a forward solve (`PETSC_TRUE`) or adjoint solve (`PETSC_FALSE`)
2869b3d3759SBarry Smith - ctx          - the context passed as the final argument to `TSSetEventHandler()`
2872dc7a7e3SShri Abhyankar 
288ca4445c7SIlya Fursov   Options Database Keys:
289ca4445c7SIlya Fursov + -ts_event_tol <tol>                       - tolerance for zero crossing check of indicator functions
290ca4445c7SIlya Fursov . -ts_event_monitor                         - print choices made by event handler
291ca4445c7SIlya Fursov . -ts_event_recorder_initial_size <recsize> - initial size of event recorder
292fe4ad979SIlya Fursov . -ts_event_post_event_step <dt1>           - first time step after event
293fe4ad979SIlya Fursov . -ts_event_post_event_second_step <dt2>    - second time step after event
294ca4445c7SIlya Fursov - -ts_event_dt_min <dt>                     - minimum time step considered for TSEvent
295ca4445c7SIlya Fursov 
2962dc7a7e3SShri Abhyankar   Level: intermediate
2972dc7a7e3SShri Abhyankar 
298ca4445c7SIlya Fursov   Notes:
299ca4445c7SIlya Fursov   The indicator functions should be defined in the `indicator` callback using the components of solution `U` and/or time `t`.
300ca4445c7SIlya Fursov   Note that `U` is `PetscScalar`-valued, and the indicator functions are `PetscReal`-valued. It is the user's responsibility to
301ca4445c7SIlya Fursov   properly handle this difference, e.g. by applying `PetscRealPart()` or other appropriate conversion means.
302ca4445c7SIlya Fursov 
303ca4445c7SIlya Fursov   The full set of events is distributed (by the user design) across MPI processes, with each process defining its own local sub-set of events.
304ca4445c7SIlya Fursov   However, the `postevent()` callback invocation is performed synchronously on all processes, including
305ca4445c7SIlya Fursov   those processes which have not currently triggered any events.
306ca4445c7SIlya Fursov 
3071cc06b55SBarry Smith .seealso: [](ch_ts), `TSEvent`, `TSCreate()`, `TSSetTimeStep()`, `TSSetConvergedReason()`
3082dc7a7e3SShri Abhyankar @*/
TSSetEventHandler(TS ts,PetscInt nevents,PetscInt direction[],PetscBool terminate[],PetscErrorCode (* indicator)(TS ts,PetscReal t,Vec U,PetscReal fvalue[],PetscCtx ctx),PetscErrorCode (* postevent)(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,PetscCtx ctx),PetscCtx ctx)3092a8381b2SBarry Smith PetscErrorCode TSSetEventHandler(TS ts, PetscInt nevents, PetscInt direction[], PetscBool terminate[], PetscErrorCode (*indicator)(TS ts, PetscReal t, Vec U, PetscReal fvalue[], PetscCtx ctx), PetscErrorCode (*postevent)(TS ts, PetscInt nevents_zero, PetscInt events_zero[], PetscReal t, Vec U, PetscBool forwardsolve, PetscCtx ctx), PetscCtx ctx)
310d71ae5a4SJacob Faibussowitsch {
311d294eb03SHong Zhang   TSAdapt   adapt;
312d294eb03SHong Zhang   PetscReal hmin;
3132dc7a7e3SShri Abhyankar   TSEvent   event;
314006e6a18SShri Abhyankar   PetscBool flg;
315*beceaeb6SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
316a6c783d2SShri Abhyankar   PetscReal tol = 1e-4;
317a6c783d2SShri Abhyankar #else
318d569cc17SSatish Balay   PetscReal tol = 1e-6;
319a6c783d2SShri Abhyankar #endif
3202dc7a7e3SShri Abhyankar 
3212dc7a7e3SShri Abhyankar   PetscFunctionBegin;
3226427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
3230a82b154SShri   if (nevents) {
3244f572ea9SToby Isaac     PetscAssertPointer(direction, 3);
3254f572ea9SToby Isaac     PetscAssertPointer(terminate, 4);
3260a82b154SShri   }
3274dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&event));
3289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->fvalue_prev));
329ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->fvalue));
3309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->fvalue_right));
331ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->fsign_prev));
332ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->fsign));
333ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->fsign_right));
3349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->side));
335ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->side_prev));
336ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->justrefined_AB));
337ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->gamma_AB));
3389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->direction));
3399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->terminate));
340ca4445c7SIlya Fursov   PetscCall(PetscMalloc1(nevents, &event->events_zero));
3419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->vtol));
342ca4445c7SIlya Fursov   for (PetscInt i = 0; i < nevents; i++) {
343d94325d3SShri Abhyankar     event->direction[i]      = direction[i];
344d94325d3SShri Abhyankar     event->terminate[i]      = terminate[i];
345ca4445c7SIlya Fursov     event->justrefined_AB[i] = PETSC_FALSE;
346ca4445c7SIlya Fursov     event->gamma_AB[i]       = 1;
347ca4445c7SIlya Fursov     event->side[i]           = 2;
348ca4445c7SIlya Fursov     event->side_prev[i]      = 0;
349d94325d3SShri Abhyankar   }
350ca4445c7SIlya Fursov   event->iterctr                = 0;
351ca4445c7SIlya Fursov   event->processing             = PETSC_FALSE;
352ca4445c7SIlya Fursov   event->revisit_right          = PETSC_FALSE;
3532589fa24SLisandro Dalcin   event->nevents                = nevents;
3549b3d3759SBarry Smith   event->indicator              = indicator;
3552dc7a7e3SShri Abhyankar   event->postevent              = postevent;
3566427ac75SLisandro Dalcin   event->ctx                    = ctx;
357fe4ad979SIlya Fursov   event->timestep_postevent     = PETSC_DECIDE;
358fe4ad979SIlya Fursov   event->timestep_2nd_postevent = PETSC_DECIDE;
3599566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
3609566063dSJacob Faibussowitsch   PetscCall(TSAdaptGetStepLimits(adapt, &hmin, NULL));
361d294eb03SHong Zhang   event->timestep_min = hmin;
3622dc7a7e3SShri Abhyankar 
36302749585SLisandro Dalcin   event->recsize = 8; /* Initial size of the recorder */
364d0609cedSBarry Smith   PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "TS Event options", "TS");
3652dc7a7e3SShri Abhyankar   {
366ca4445c7SIlya Fursov     PetscCall(PetscOptionsReal("-ts_event_tol", "Tolerance for zero crossing check of indicator functions", "TSSetEventTolerances", tol, &tol, NULL));
3679566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-ts_event_monitor", "Print choices made by event handler", "", &flg));
3689566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size", "Initial size of event recorder", "", event->recsize, &event->recsize, NULL));
369fe4ad979SIlya Fursov     PetscCall(PetscOptionsDeprecated("-ts_event_post_eventinterval_step", "-ts_event_post_event_second_step", "3.21", NULL));
370fe4ad979SIlya Fursov     PetscCall(PetscOptionsReal("-ts_event_post_event_step", "First time step after event", "", event->timestep_postevent, &event->timestep_postevent, NULL));
371fe4ad979SIlya Fursov     PetscCall(PetscOptionsReal("-ts_event_post_event_second_step", "Second time step after event", "", event->timestep_2nd_postevent, &event->timestep_2nd_postevent, NULL));
3729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_dt_min", "Minimum time step considered for TSEvent", "", event->timestep_min, &event->timestep_min, NULL));
3732dc7a7e3SShri Abhyankar   }
374d0609cedSBarry Smith   PetscOptionsEnd();
375a4ffd976SShri Abhyankar 
3769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.time));
3779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.stepnum));
3789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.nevents));
3799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.eventidx));
380ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &event->recorder.eventidx[i]));
381e7069c78SShri   /* Initialize the event recorder */
382e7069c78SShri   event->recorder.ctr = 0;
383a4ffd976SShri Abhyankar 
384ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents; i++) event->vtol[i] = tol;
3856c23a314SBarry Smith   if (flg) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts), "stdout", &event->monitor));
3869e12be75SShri Abhyankar 
3879566063dSJacob Faibussowitsch   PetscCall(TSEventDestroy(&ts->event));
388d94325d3SShri Abhyankar   ts->event        = event;
389e7069c78SShri   ts->event->refct = 1;
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3912dc7a7e3SShri Abhyankar }
3922dc7a7e3SShri Abhyankar 
393a4ffd976SShri Abhyankar /*
394a4ffd976SShri Abhyankar   TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
395a4ffd976SShri Abhyankar                           is reached.
396a4ffd976SShri Abhyankar */
TSEventRecorderResize(TSEvent event)397d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventRecorderResize(TSEvent event)
398d71ae5a4SJacob Faibussowitsch {
399a4ffd976SShri Abhyankar   PetscReal *time;
400ca4445c7SIlya Fursov   PetscInt  *stepnum, *nevents;
401a4ffd976SShri Abhyankar   PetscInt **eventidx;
402ca4445c7SIlya Fursov   PetscInt   fact = 2;
403a4ffd976SShri Abhyankar 
404a4ffd976SShri Abhyankar   PetscFunctionBegin;
405ca4445c7SIlya Fursov   /* Create larger arrays */
4069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact * event->recsize, &time));
4079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact * event->recsize, &stepnum));
4089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact * event->recsize, &nevents));
4099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact * event->recsize, &eventidx));
410ca4445c7SIlya Fursov   for (PetscInt i = 0; i < fact * event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &eventidx[i]));
411a4ffd976SShri Abhyankar 
412a4ffd976SShri Abhyankar   /* Copy over data */
4139566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(time, event->recorder.time, event->recsize));
4149566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(stepnum, event->recorder.stepnum, event->recsize));
4159566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(nevents, event->recorder.nevents, event->recsize));
416ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->recsize; i++) PetscCall(PetscArraycpy(eventidx[i], event->recorder.eventidx[i], event->recorder.nevents[i]));
417a4ffd976SShri Abhyankar 
418a4ffd976SShri Abhyankar   /* Destroy old arrays */
419ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->recsize; i++) PetscCall(PetscFree(event->recorder.eventidx[i]));
4209566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.eventidx));
4219566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.nevents));
4229566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.stepnum));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.time));
424a4ffd976SShri Abhyankar 
425a4ffd976SShri Abhyankar   /* Set pointers */
426a4ffd976SShri Abhyankar   event->recorder.time     = time;
427a4ffd976SShri Abhyankar   event->recorder.stepnum  = stepnum;
428a4ffd976SShri Abhyankar   event->recorder.nevents  = nevents;
429a4ffd976SShri Abhyankar   event->recorder.eventidx = eventidx;
430a4ffd976SShri Abhyankar 
431ca4445c7SIlya Fursov   /* Update the size */
432a4ffd976SShri Abhyankar   event->recsize *= fact;
4333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
434a4ffd976SShri Abhyankar }
435a4ffd976SShri Abhyankar 
436031fbad4SShri Abhyankar /*
437ca4445c7SIlya Fursov   Helper routine to handle user postevents and recording
438031fbad4SShri Abhyankar */
TSPostEvent(TS ts,PetscReal t,Vec U)439d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPostEvent(TS ts, PetscReal t, Vec U)
440d71ae5a4SJacob Faibussowitsch {
4412dc7a7e3SShri Abhyankar   TSEvent   event        = ts->event;
44228d5b5d6SLisandro Dalcin   PetscBool restart      = PETSC_FALSE;
443ca4445c7SIlya Fursov   PetscBool terminate    = PETSC_FALSE;
444ca4445c7SIlya Fursov   PetscBool statechanged = PETSC_FALSE;
445ca4445c7SIlya Fursov   PetscInt  ctr, stepnum;
446ca4445c7SIlya Fursov   PetscBool inflag[3], outflag[3];
447ca4445c7SIlya Fursov   PetscBool forwardsolve = PETSC_TRUE; // Flag indicating that TS is doing a forward solve
4482dc7a7e3SShri Abhyankar 
4492dc7a7e3SShri Abhyankar   PetscFunctionBegin;
4502dc7a7e3SShri Abhyankar   if (event->postevent) {
45128d5b5d6SLisandro Dalcin     PetscObjectState state_prev, state_post;
4529566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateGet((PetscObject)U, &state_prev));
453ca4445c7SIlya Fursov     PetscCallBack("TSEvent post-event processing", (*event->postevent)(ts, event->nevents_zero, event->events_zero, t, U, forwardsolve, event->ctx)); // TODO update 'restart' here?
4549566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateGet((PetscObject)U, &state_post));
455ca4445c7SIlya Fursov     if (state_prev != state_post) {
456ca4445c7SIlya Fursov       restart      = PETSC_TRUE;
457ca4445c7SIlya Fursov       statechanged = PETSC_TRUE;
458ca4445c7SIlya Fursov     }
4592dc7a7e3SShri Abhyankar   }
4604597913aSLisandro Dalcin 
461ca4445c7SIlya Fursov   // Handle termination events and step restart
462ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents_zero; i++)
4639371c9d4SSatish Balay     if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
4649371c9d4SSatish Balay   inflag[0] = restart;
4659371c9d4SSatish Balay   inflag[1] = terminate;
466ca4445c7SIlya Fursov   inflag[2] = statechanged;
4675440e5dcSBarry Smith   PetscCallMPI(MPIU_Allreduce(inflag, outflag, 3, MPI_C_BOOL, MPI_LOR, ((PetscObject)ts)->comm));
4689371c9d4SSatish Balay   restart      = outflag[0];
4699371c9d4SSatish Balay   terminate    = outflag[1];
470ca4445c7SIlya Fursov   statechanged = outflag[2];
4719566063dSJacob Faibussowitsch   if (restart) PetscCall(TSRestartStep(ts));
4729566063dSJacob Faibussowitsch   if (terminate) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_EVENT));
473f7aea88cSShri Abhyankar 
474ca4445c7SIlya Fursov   /*
475ca4445c7SIlya Fursov     Recalculate the indicator functions and signs if the state has been changed by the user postevent callback.
476ca4445c7SIlya Fursov     Note! If the state HAS NOT changed, the existing event->fsign (equal to zero) is kept, which:
477ca4445c7SIlya Fursov     - might have been defined using the previous (now-possibly-overridden) event->vtol,
478ca4445c7SIlya Fursov     - might have been set to zero on reaching a small time step rather than using the vtol criterion.
479ca4445c7SIlya Fursov     This will enforce keeping event->fsign = 0 where the zero-crossings were actually marked,
480ca4445c7SIlya Fursov     resulting in a more consistent behaviour of fsign's.
481ca4445c7SIlya Fursov   */
482ca4445c7SIlya Fursov   if (statechanged) {
483ca4445c7SIlya Fursov     if (event->monitor) PetscCall(PetscPrintf(((PetscObject)ts)->comm, "TSEvent: at time %g the vector state has been changed by PostEvent, recalculating fvalues and signs\n", (double)t));
4849566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(U));
485ca4445c7SIlya Fursov     PetscCallBack("TSEvent indicator", (*event->indicator)(ts, t, U, event->fvalue, event->ctx));
4869566063dSJacob Faibussowitsch     PetscCall(VecLockReadPop(U));
487ca4445c7SIlya Fursov     TSEventCalcSigns(event->nevents, event->fvalue, event->vtol, event->fsign); // note, event->vtol might have been changed by the postevent()
488f443add6SLisandro Dalcin   }
489f443add6SLisandro Dalcin 
490ca4445c7SIlya Fursov   // Record the event in the event recorder
4919566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &stepnum));
492f7aea88cSShri Abhyankar   ctr = event->recorder.ctr;
49348a46eb9SPierre Jolivet   if (ctr == event->recsize) PetscCall(TSEventRecorderResize(event));
494f7aea88cSShri Abhyankar   event->recorder.time[ctr]    = t;
495d0578d90SShri Abhyankar   event->recorder.stepnum[ctr] = stepnum;
4964597913aSLisandro Dalcin   event->recorder.nevents[ctr] = event->nevents_zero;
497ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
498f7aea88cSShri Abhyankar   event->recorder.ctr++;
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5002dc7a7e3SShri Abhyankar }
5012dc7a7e3SShri Abhyankar 
502ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-sowing-chars
503ca4445c7SIlya Fursov /*
504ca4445c7SIlya Fursov   (modified) Anderson-Bjorck variant of regula falsi method, refines [tleft, t] or [t, tright] based on 'side' (-1 or +1).
505ca4445c7SIlya Fursov   The scaling parameter is defined based on the 'justrefined' flag, the history of repeats of 'side', and the threshold.
506ca4445c7SIlya Fursov   To escape certain failure modes, the algorithm may drift towards the bisection rule.
507ca4445c7SIlya Fursov   The value pointed to by 'side_prev' gets updated.
508ca4445c7SIlya Fursov   This function returns the new time step.
509ca4445c7SIlya Fursov 
510ca4445c7SIlya Fursov   The underlying pure Anderson-Bjorck algorithm was taken as described in
511ca4445c7SIlya Fursov   J.M. Fernandez-Diaz, C.O. Menendez-Perez "A common framework for modified Regula Falsi methods and new methods of this kind", 2023.
512ca4445c7SIlya Fursov   The modifications subsequently introduced have little effect on the behaviour for simple cases requiring only a few iterations
513ca4445c7SIlya Fursov   (some minor convergence slowdown may take place though), but the effect becomes more pronounced for tough cases requiring many iterations.
514ca4445c7SIlya Fursov   For the latter situation the speed-up may be order(s) of magnitude compared to the classical Anderson-Bjorck.
515ca4445c7SIlya Fursov   The modifications (the threshold trick, and the drift towards bisection) were tested and tweaked
516ca4445c7SIlya Fursov   based on a number of test functions from the mentioned paper.
517ca4445c7SIlya Fursov */
RefineAndersonBjorck(PetscReal tleft,PetscReal t,PetscReal tright,PetscReal fleft,PetscReal f,PetscReal fright,PetscInt side,PetscInt * side_prev,PetscBool justrefined,PetscReal * gamma)518ca4445c7SIlya Fursov static inline PetscReal RefineAndersonBjorck(PetscReal tleft, PetscReal t, PetscReal tright, PetscReal fleft, PetscReal f, PetscReal fright, PetscInt side, PetscInt *side_prev, PetscBool justrefined, PetscReal *gamma)
519d71ae5a4SJacob Faibussowitsch {
520ca4445c7SIlya Fursov   PetscReal      new_dt, scal = 1.0, scalB = 1.0, threshold = 0.0, power;
521ca4445c7SIlya Fursov   PetscInt       reps     = 0;
522ca4445c7SIlya Fursov   const PetscInt REPS_CAP = 8; // an upper bound to be imposed on 'reps' (set to 8, somewhat arbitrary number, found after some tweaking)
523ca4445c7SIlya Fursov 
524ca4445c7SIlya Fursov   // Preparations
525ca4445c7SIlya Fursov   if (justrefined) {
526ca4445c7SIlya Fursov     if (*side_prev * side > 0) *side_prev += side;     // the side keeps repeating -> increment the side counter (-ve or +ve)
527ca4445c7SIlya Fursov     else *side_prev = side;                            // reset the counter
528ca4445c7SIlya Fursov     reps      = PetscMin(*side_prev * side, REPS_CAP); // the length of the recent side-repeat series, including the current 'side'
529ca4445c7SIlya Fursov     threshold = PetscPowReal(0.5, reps) * 0.1;         // ad-hoc strategy for threshold calculation (involved some tweaking)
530ca4445c7SIlya Fursov   } else *side_prev = side;                            // initial reset of the counter
531ca4445c7SIlya Fursov 
532ca4445c7SIlya Fursov   // Time step calculation
533e2cdd850SShri Abhyankar   if (side == -1) {
534ca4445c7SIlya Fursov     if (justrefined && fright != 0.0 && fleft != 0.0) {
535ca4445c7SIlya Fursov       scal  = (fright - f) / fright;
536ca4445c7SIlya Fursov       scalB = -f / fleft;
537e2cdd850SShri Abhyankar     }
538ca4445c7SIlya Fursov   } else { // must be side == +1
539ca4445c7SIlya Fursov     if (justrefined && fleft != 0.0 && fright != 0.0) {
540ca4445c7SIlya Fursov       scal  = (fleft - f) / fleft;
541ca4445c7SIlya Fursov       scalB = -f / fright;
542e2cdd850SShri Abhyankar     }
54338bf2713SShri Abhyankar   }
544e2cdd850SShri Abhyankar 
545ca4445c7SIlya Fursov   if (scal < threshold) scal = 0.5;
546ca4445c7SIlya Fursov   if (reps > 1) *gamma *= scal; // side did not switch since the last time, accumulate gamma
547ca4445c7SIlya Fursov   else *gamma = 1.0;            // side switched -> reset gamma
548ca4445c7SIlya Fursov   power = PetscMax(0.0, (reps - 2.0) / (REPS_CAP - 2.0));
549ca4445c7SIlya Fursov   scal  = PetscPowReal(scalB / *gamma, power) * (*gamma); // mix the Anderson-Bjorck scaling and Bisection scaling
550ca4445c7SIlya Fursov 
551ca4445c7SIlya Fursov   if (side == -1) new_dt = scal * fleft / (scal * fleft - f) * (t - tleft);
552ca4445c7SIlya Fursov   else new_dt = f / (f - scal * fright) * (tright - t);
553ca4445c7SIlya Fursov   /*
554ca4445c7SIlya Fursov     In tough cases (e.g. a polynomial of high order), there is a failure mode for the standard Anderson-Bjorck,
555ca4445c7SIlya Fursov     when the new proposed point jumps from one end-point of the bracket to the other, however the bracket
556ca4445c7SIlya Fursov     is contracting very slowly. A larger threshold for 'scal' prevents entering this mode.
557ca4445c7SIlya Fursov     On the other hand, if the iteration gets stuck near one end-point of the bracket, and the 'side' does not switch for a while,
558ca4445c7SIlya Fursov     the 'scal' drifts towards the bisection approach (via scalB), ensuring stable convergence.
559ca4445c7SIlya Fursov   */
560ca4445c7SIlya Fursov   return new_dt;
561ca4445c7SIlya Fursov }
562ca4445c7SIlya Fursov 
563ca4445c7SIlya Fursov /*
564ca4445c7SIlya Fursov   Checks if the current point (t) is the zero-crossing location, based on the indicator function signs and direction[]:
565ca4445c7SIlya Fursov   - using the dt_min criterion,
566ca4445c7SIlya Fursov   - using the vtol criterion.
567ca4445c7SIlya Fursov   The situation (fsign_prev, fsign) = (0, 0) is treated as staying in the near-zero-zone of the previous zero-crossing,
568ca4445c7SIlya Fursov   and is not marked as a new zero-crossing.
569ca4445c7SIlya Fursov   This function may update event->side[].
570ca4445c7SIlya Fursov */
TSEventTestZero(TS ts,PetscReal t)571ca4445c7SIlya Fursov static PetscErrorCode TSEventTestZero(TS ts, PetscReal t)
572d71ae5a4SJacob Faibussowitsch {
573d294eb03SHong Zhang   TSEvent event = ts->event;
574d294eb03SHong Zhang 
575d294eb03SHong Zhang   PetscFunctionBegin;
576ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents; i++) {
577ca4445c7SIlya Fursov     const PetscBool bracket_is_left = (event->fsign_prev[i] * event->fsign[i] < 0 && event->fsign[i] * event->direction[i] >= 0) ? PETSC_TRUE : PETSC_FALSE;
578d294eb03SHong Zhang 
579ca4445c7SIlya Fursov     if (bracket_is_left && ((t - event->ptime_prev <= event->timestep_min) || event->revisit_right)) event->side[i] = 0;          // mark zero-crossing from dt_min; 'bracket_is_left' accounts for direction
580ca4445c7SIlya Fursov     if (event->fsign[i] == 0 && event->fsign_prev[i] != 0 && event->fsign_prev[i] * event->direction[i] <= 0) event->side[i] = 0; // mark zero-crossing from vtol
581d294eb03SHong Zhang   }
5823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
583d294eb03SHong Zhang }
58438bf2713SShri Abhyankar 
585ca4445c7SIlya Fursov /*
586ca4445c7SIlya Fursov   Checks if [fleft, f] or [f, fright] are 'brackets', i.e. intervals with the sign change, satisfying the 'direction'.
587ca4445c7SIlya Fursov   The right interval is only checked if iterctr > 0 (i.e. Anderson-Bjorck refinement has started).
588ca4445c7SIlya Fursov   The intervals like [0, x] and [x, 0] are not counted as brackets, i.e. intervals with the sign change.
589ca4445c7SIlya Fursov   The function returns the 'side' value: -1 (left, or both are brackets), +1 (only right one), +2 (neither).
590ca4445c7SIlya Fursov */
TSEventTestBracket(PetscInt fsign_left,PetscInt fsign,PetscInt fsign_right,PetscInt direction,PetscInt iterctr)591ca4445c7SIlya Fursov static inline PetscInt TSEventTestBracket(PetscInt fsign_left, PetscInt fsign, PetscInt fsign_right, PetscInt direction, PetscInt iterctr)
592ca4445c7SIlya Fursov {
593ca4445c7SIlya Fursov   PetscInt side = 2;
594ca4445c7SIlya Fursov   if (fsign_left * fsign < 0 && fsign * direction >= 0) side = -1;
595ca4445c7SIlya Fursov   if (side != -1 && iterctr > 0 && fsign * fsign_right < 0 && fsign_right * direction >= 0) side = 1;
596ca4445c7SIlya Fursov   return side;
597ca4445c7SIlya Fursov }
598ca4445c7SIlya Fursov 
599ca4445c7SIlya Fursov /*
600136cf249SJames Wright   Caps the time steps, accounting for evaluation time points.
601136cf249SJames Wright   It uses 'event->timestep_cache' as a time step to calculate the tolerance for eval_times points detection. This
602ca4445c7SIlya Fursov   is done since the event resolution may result in significant time step refinement, and we don't use these small steps for tolerances.
603136cf249SJames Wright   To enhance the consistency of eval_times points detection, tolerance 'eval_times->worktol' is reused later in the TSSolve iteration.
604ca4445c7SIlya Fursov   If a user-defined step is cut by this function, the input uncut step is saved to adapt->dt_span_cached.
605ca4445c7SIlya Fursov   Flag 'user_dt' indicates if the step was defined by user.
606ca4445c7SIlya Fursov */
TSEvent_dt_cap(TS ts,PetscReal t,PetscReal dt,PetscBool user_dt)607ca4445c7SIlya Fursov static inline PetscReal TSEvent_dt_cap(TS ts, PetscReal t, PetscReal dt, PetscBool user_dt)
608ca4445c7SIlya Fursov {
609ca4445c7SIlya Fursov   PetscReal res = dt;
610ca4445c7SIlya Fursov   if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
611136cf249SJames Wright     PetscReal maxdt    = ts->max_time - t; // this may be overridden by eval_times
612ca4445c7SIlya Fursov     PetscBool cut_made = PETSC_FALSE;
613ca4445c7SIlya Fursov     PetscReal eps      = 10 * PETSC_MACHINE_EPSILON;
614136cf249SJames Wright     if (ts->eval_times) {
615136cf249SJames Wright       PetscInt   idx = ts->eval_times->time_point_idx;
616136cf249SJames Wright       PetscInt   Ns  = ts->eval_times->num_time_points;
617136cf249SJames Wright       PetscReal *st  = ts->eval_times->time_points;
618ca4445c7SIlya Fursov 
619136cf249SJames Wright       if (ts->eval_times->worktol == 0) ts->eval_times->worktol = ts->eval_times->reltol * ts->event->timestep_cache + ts->eval_times->abstol; // in case TSAdaptChoose() has not defined it
620136cf249SJames Wright       if (idx < Ns && PetscIsCloseAtTol(t, st[idx], ts->eval_times->worktol, 0)) {                                                             // just hit a evaluation time point
621136cf249SJames Wright         if (idx + 1 < Ns) maxdt = st[idx + 1] - t;                                                                                             // ok to use the next evaluation time point
622136cf249SJames Wright         else maxdt = ts->max_time - t;                                                                                                         // can't use the next evaluation time point: they have finished
623136cf249SJames Wright       } else if (idx < Ns) maxdt = st[idx] - t;                                                                                                // haven't hit a evaluation time point, use the nearest one
624ca4445c7SIlya Fursov     }
625ca4445c7SIlya Fursov     maxdt = PetscMin(maxdt, ts->max_time - t);
6268343f784SJames Wright     PetscCheck((maxdt > eps) || (PetscAbsReal(maxdt) <= eps && PetscIsCloseAtTol(t, ts->max_time, eps, 0)), PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Unexpected state: bad maxdt (%g) in TSEvent_dt_cap()", (double)maxdt);
627ca4445c7SIlya Fursov 
628ca4445c7SIlya Fursov     if (PetscIsCloseAtTol(dt, maxdt, eps, 0)) res = maxdt; // no cut
629ca4445c7SIlya Fursov     else {
630ca4445c7SIlya Fursov       if (dt > maxdt) {
631ca4445c7SIlya Fursov         res      = maxdt; // yes cut
632ca4445c7SIlya Fursov         cut_made = PETSC_TRUE;
633ca4445c7SIlya Fursov       } else res = dt; // no cut
634ca4445c7SIlya Fursov     }
635ca4445c7SIlya Fursov     if (ts->adapt && user_dt) { // only update dt_span_cached for the user-defined step
636136cf249SJames Wright       if (cut_made) ts->adapt->dt_eval_times_cached = dt;
637136cf249SJames Wright       else ts->adapt->dt_eval_times_cached = 0;
638ca4445c7SIlya Fursov     }
639ca4445c7SIlya Fursov   }
640ca4445c7SIlya Fursov   return res;
641ca4445c7SIlya Fursov }
642ca4445c7SIlya Fursov 
643ca4445c7SIlya Fursov /*
644ca4445c7SIlya Fursov   Updates the left-end values
645ca4445c7SIlya Fursov */
TSEvent_update_left(TSEvent event,PetscReal t)646ca4445c7SIlya Fursov static inline void TSEvent_update_left(TSEvent event, PetscReal t)
647ca4445c7SIlya Fursov {
648ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents; i++) {
649ca4445c7SIlya Fursov     event->fvalue_prev[i] = event->fvalue[i];
650ca4445c7SIlya Fursov     event->fsign_prev[i]  = event->fsign[i];
651ca4445c7SIlya Fursov   }
652ca4445c7SIlya Fursov   event->ptime_prev = t;
653ca4445c7SIlya Fursov }
654ca4445c7SIlya Fursov 
655ca4445c7SIlya Fursov /*
656ca4445c7SIlya Fursov   Updates the right-end values
657ca4445c7SIlya Fursov */
TSEvent_update_right(TSEvent event,PetscReal t)658ca4445c7SIlya Fursov static inline void TSEvent_update_right(TSEvent event, PetscReal t)
659ca4445c7SIlya Fursov {
660ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents; i++) {
661ca4445c7SIlya Fursov     event->fvalue_right[i] = event->fvalue[i];
662ca4445c7SIlya Fursov     event->fsign_right[i]  = event->fsign[i];
663ca4445c7SIlya Fursov   }
664ca4445c7SIlya Fursov   event->ptime_right = t;
665ca4445c7SIlya Fursov }
666ca4445c7SIlya Fursov 
667ca4445c7SIlya Fursov /*
668ca4445c7SIlya Fursov   Updates the current values from the right-end values
669ca4445c7SIlya Fursov */
TSEvent_update_from_right(TSEvent event)670ca4445c7SIlya Fursov static inline PetscReal TSEvent_update_from_right(TSEvent event)
671ca4445c7SIlya Fursov {
672ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents; i++) {
673ca4445c7SIlya Fursov     event->fvalue[i] = event->fvalue_right[i];
674ca4445c7SIlya Fursov     event->fsign[i]  = event->fsign_right[i];
675ca4445c7SIlya Fursov   }
676ca4445c7SIlya Fursov   return event->ptime_right;
677ca4445c7SIlya Fursov }
678ca4445c7SIlya Fursov 
Not_PETSC_DECIDE(PetscReal dt)679fe4ad979SIlya Fursov static inline PetscBool Not_PETSC_DECIDE(PetscReal dt)
680fe4ad979SIlya Fursov {
6814ad8454bSPierre Jolivet   return dt == PETSC_DECIDE ? PETSC_FALSE : PETSC_TRUE;
682fe4ad979SIlya Fursov }
683fe4ad979SIlya Fursov 
684ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-spacing
685ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-unknown
686ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-spelling
687ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-missing
688ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-fishy-header
689ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
690ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-synopsis-missing-description
691ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-sowing-chars
692ca4445c7SIlya Fursov /*
693ca4445c7SIlya Fursov   TSEventHandler - the main function to perform a single iteration of event detection.
694ca4445c7SIlya Fursov 
695ca4445c7SIlya Fursov   Developer notes:
696fe4ad979SIlya Fursov   A) The 'event->iterctr > 0' is used as an indicator that Anderson-Bjorck refinement has started.
697fe4ad979SIlya Fursov   B) If event->iterctr == 0, then justrefined_AB[i] is always false.
698fe4ad979SIlya Fursov   C) The right-end quantities: ptime_right, fvalue_right[i] and fsign_right[i] are only guaranteed to be valid
699ca4445c7SIlya Fursov   for event->iterctr > 0.
700fe4ad979SIlya Fursov   D) If event->iterctr > 0, then event->processing is PETSC_TRUE; the opposite may not hold.
701fe4ad979SIlya Fursov   E) When event->processing == PETSC_TRUE and event->iterctr == 0, the event handler iterations are complete, but
702fe4ad979SIlya Fursov   the event handler continues managing the 1st and 2nd post-event steps. In this case the 1st post-event step
703fe4ad979SIlya Fursov   proposed by the event handler is not checked by TSAdapt, and is always accepted (beware!).
704fe4ad979SIlya Fursov   However, if the 2nd post-event step is not managed by the event handler (e.g. 1st = numerical, 2nd = PETSC_DECIDE),
705fe4ad979SIlya Fursov   condition "E" does not hold, and TSAdapt may reject/adjust the 1st post-event step.
706fe4ad979SIlya Fursov   F) event->side[i] may take values: 0 <=> point t is a zero-crossing for indicator function i (via vtol/dt_min criterion);
707ca4445c7SIlya Fursov   -1/+1 <=> detected a bracket to the left/right of t for indicator function i; +2 <=> no brackets/zero-crossings.
708fe4ad979SIlya Fursov   G) The signs event->fsign[i] (with values 0/-1/+1) are calculated for each new point. Zero sign is set if the function value is
709ca4445c7SIlya Fursov   smaller than the tolerance. Besides, zero sign is enforced after marking a zero-crossing due to small bracket size criterion.
710ca4445c7SIlya Fursov 
711ca4445c7SIlya Fursov   The intervals with the indicator function sign change (i.e. containing the potential zero-crossings) are called 'brackets'.
712ca4445c7SIlya Fursov   To find a zero-crossing, the algorithm first locates a bracket, and then sequentially subdivides it, generating a sequence
713ca4445c7SIlya Fursov   of brackets whose length tends to zero. The bracket subdivision involves the (modified) Anderson-Bjorck method.
714ca4445c7SIlya Fursov 
715ca4445c7SIlya Fursov   Apart from the comments scattered throughout the code to clarify different lines and blocks,
716ca4445c7SIlya Fursov   a few tricky aspects of the algorithm (and the underlying reasoning) are discussed in detail below:
717ca4445c7SIlya Fursov 
718ca4445c7SIlya Fursov   =Sign tracking=
719ca4445c7SIlya Fursov   When a zero-crossing is found, the sign variable (event->fsign[i]) is set to zero for the current point t.
720ca4445c7SIlya Fursov   This happens both for zero-crossings triggered via the vtol criterion, and those triggered via the dt_min
721ca4445c7SIlya Fursov   criterion. After the event, as the TS steps forward, the current sign values are handed over to event->fsign_prev[i].
722ca4445c7SIlya Fursov   The recalculation of signs is avoided if possible: e.g. if a 'vtol' criterion resulted in a zero-crossing at point t,
723ca4445c7SIlya Fursov   but the subsequent call to postevent() handler decreased 'vtol', making the indicator function no longer "close to zero"
724ca4445c7SIlya Fursov   at point t, the fsign[i] will still consistently keep the zero value. This allows avoiding the erroneous duplication
725ca4445c7SIlya Fursov   of events:
726ca4445c7SIlya Fursov 
727ca4445c7SIlya Fursov   E.g. consider a bracket [t0, t2], where f0 < 0, f2 > 0, which resulted in a zero-crossing t1 with f1 < 0, abs(f1) < vtol.
728ca4445c7SIlya Fursov   Suppose the postevent() handler changes vtol to vtol*, such that abs(f1) > vtol*. The TS makes a step t1 -> t3, where
729ca4445c7SIlya Fursov   again f1 < 0, f3 > 0, and the event handler will find a new event near t1, which is actually a duplication of the
730ca4445c7SIlya Fursov   original event at t1. The duplications are avoided by NOT counting the sign progressions 0 -> +1, or 0 -> -1
731ca4445c7SIlya Fursov   as brackets. Tracking (instead of recalculating) the sign values makes this procedure work more consistently.
732ca4445c7SIlya Fursov 
733ca4445c7SIlya Fursov   The sign values are however recalculated if the postevent() callback has changed the current solution vector U
734ca4445c7SIlya Fursov   (such a change resets everything).
735ca4445c7SIlya Fursov   The sign value is also set to zero if the dt_min criterion has triggered the event. This allows the algorithm to
736ca4445c7SIlya Fursov   work consistently, irrespective of the type of criterion involved (vtol/dt_min).
737ca4445c7SIlya Fursov 
738ca4445c7SIlya Fursov   =Event from min bracket=
739ca4445c7SIlya Fursov   When the event handler ends up with a bracket [t0, t1] with size <= dt_min, a zero crossing is reported at t1,
740ca4445c7SIlya Fursov   and never at t0. If such a bracket is discovered when TS is staying at t0, one more step forward (to t1) is necessary
741ca4445c7SIlya Fursov   to mark the found event. This is the situation of revisiting t1, which is described below (see Revisiting).
742ca4445c7SIlya Fursov 
743ca4445c7SIlya Fursov   Why t0 is not reported as event location? Suppose it is, and let f0 < 0, f1 > 0. Also suppose that the
744ca4445c7SIlya Fursov   postevent() handler has slightly changed the solution U, so the sign at t0 is recalculated: it equals -1. As the TS steps
745ca4445c7SIlya Fursov   further: t0 -> t2, with sign0 == -1, and sign2 == +1, the event handler will locate the bracket [t0, t2], eventually
746ca4445c7SIlya Fursov   resolving a new event near t1, i.e. finding a duplicate event.
747ca4445c7SIlya Fursov   This situation is avoided by reporting the event at t1 in the first place.
748ca4445c7SIlya Fursov 
749ca4445c7SIlya Fursov   =Revisiting=
750ca4445c7SIlya Fursov   When handling the situation with small bracket size, the TS solver may happen to visit the same point twice,
751ca4445c7SIlya Fursov   but with different results.
752ca4445c7SIlya Fursov 
753ca4445c7SIlya Fursov   E.g. originally it discovered a bracket with sign change [t0, t10], and started resolving the zero-crossing,
754ca4445c7SIlya Fursov   visiting the points t1,...,t9 : t0 < t1 < ... < t9 < t10. Suppose that at t9 the algorithm discovers
755ca4445c7SIlya Fursov   that [t9, t10] is a bracket with the sign change it was looking for, and that |t10 - t9| is too small.
756ca4445c7SIlya Fursov   So point t10 should be revisited and marked as the zero crossing (by the minimum bracket size criterion).
757ca4445c7SIlya Fursov   On re-visiting t10, via the refined sequence of steps t0,...,t10, the TS solver may arrive at a solution U*
758ca4445c7SIlya Fursov   different from the solution U it found at t10 originally. Hence, the indicator functions at t10 may become different,
759ca4445c7SIlya Fursov   and the condition of the sign change, which existed originally, may disappear, breaking the logic of the algorithm.
760ca4445c7SIlya Fursov 
761ca4445c7SIlya Fursov   To handle such (-=unlikely=-, but possible) situations, two strategies can be considered:
762ca4445c7SIlya Fursov   1) [not used here] Allow the brackets with sign change to disappear during iterations. The algorithm should be able
763ca4445c7SIlya Fursov   to cleanly exit the iteration and leave all the objects/variables/caches involved in a valid state.
764ca4445c7SIlya Fursov   2) [ADOPTED HERE!] On revisiting t10, the event handler reuses the indicator functions previously calculated for the
765ca4445c7SIlya Fursov   original solution U. This U may be less precise than U*, but this trick does not allow the algorithm logic to break down.
766ca4445c7SIlya Fursov   HOWEVER, the original U is not stored anywhere, it is essentially lost since the TS performed the rollback from it.
767ca4445c7SIlya Fursov   On revisiting t10, the updated solution U* will inevitably be found and used everywhere EXCEPT the current
768ca4445c7SIlya Fursov   indicator functions calculation, e.g. U* will be used in the postevent() handler call. Since t10 is the event location,
769ca4445c7SIlya Fursov   the appropriate indicator-function-signs will be enforced to be 0 (regardless if the solution was U or U*).
770ca4445c7SIlya Fursov   If the solution is then changed by the postevent(), the indicator-function-signs will be recalculated.
771ca4445c7SIlya Fursov 
772ca4445c7SIlya Fursov   Whether the algorithm is revisiting a point in the current TSEventHandler() call is flagged by 'event->revisit_right'.
773ca4445c7SIlya Fursov */
TSEventHandler(TS ts)774d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventHandler(TS ts)
775d71ae5a4SJacob Faibussowitsch {
7766427ac75SLisandro Dalcin   TSEvent   event;
777ca4445c7SIlya Fursov   PetscReal t, dt_next = 0.0;
7782dc7a7e3SShri Abhyankar   Vec       U;
779ca4445c7SIlya Fursov   PetscInt  minsidein = 2, minsideout = 2; // minsideout is sync on all ranks
780ca4445c7SIlya Fursov   PetscBool finished = PETSC_FALSE;        // should stay sync on all ranks
781ca4445c7SIlya Fursov   PetscBool revisit_right_cache;           // [sync] flag for inner consistency checks
7822dc7a7e3SShri Abhyankar 
7832dc7a7e3SShri Abhyankar   PetscFunctionBegin;
7846427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
785ca4445c7SIlya Fursov 
7863ba16761SJacob Faibussowitsch   if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS);
7876427ac75SLisandro Dalcin   event               = ts->event;
788ca4445c7SIlya Fursov   event->nevents_zero = 0;
789ca4445c7SIlya Fursov   revisit_right_cache = event->revisit_right;
790ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents; i++) event->side[i] = 2; // side's are reset on each new iteration
791ca4445c7SIlya Fursov   if (event->iterctr == 0)
792ca4445c7SIlya Fursov     for (PetscInt i = 0; i < event->nevents; i++) event->justrefined_AB[i] = PETSC_FALSE;
7932dc7a7e3SShri Abhyankar 
7949566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
795ca4445c7SIlya Fursov   if (!event->processing) { // update the caches
796ca4445c7SIlya Fursov     PetscReal dt;
7979566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
798ca4445c7SIlya Fursov     event->ptime_cache    = t;
799ca4445c7SIlya Fursov     event->timestep_cache = dt; // the next TS move is planned to be: t -> t+dt
8002dc7a7e3SShri Abhyankar   }
801fe4ad979SIlya Fursov   if (event->processing && event->iterctr == 0 && Not_PETSC_DECIDE(event->timestep_2nd_postevent)) { // update the caches while processing the post-event steps
802fe4ad979SIlya Fursov     event->ptime_cache    = t;
803fe4ad979SIlya Fursov     event->timestep_cache = event->timestep_2nd_postevent;
804fe4ad979SIlya Fursov   }
8052dc7a7e3SShri Abhyankar 
806ca4445c7SIlya Fursov   PetscCall(TSGetSolution(ts, &U)); // if revisiting, this will be the updated U* (see discussion on "Revisiting" in the Developer notes above)
807ca4445c7SIlya Fursov   if (event->revisit_right) {
808ca4445c7SIlya Fursov     PetscReal tr = TSEvent_update_from_right(event);
809ca4445c7SIlya Fursov     PetscCheck(PetscAbsReal(tr - t) < PETSC_SMALL, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Inconsistent time value when performing 'revisiting' in TSEventHandler()");
8102dc7a7e3SShri Abhyankar   } else {
811ca4445c7SIlya Fursov     PetscCall(VecLockReadPush(U));
812ca4445c7SIlya Fursov     PetscCallBack("TSEvent indicator", (*event->indicator)(ts, t, U, event->fvalue, event->ctx)); // fill fvalue's at point 't'
813ca4445c7SIlya Fursov     PetscCall(VecLockReadPop(U));
814ca4445c7SIlya Fursov     TSEventCalcSigns(event->nevents, event->fvalue, event->vtol, event->fsign); // fill fvalue signs
815ca4445c7SIlya Fursov   }
816ca4445c7SIlya Fursov   PetscCall(TSEventTestZero(ts, t)); // check if the current point 't' is the event location; event->side[] may get updated
817ca4445c7SIlya Fursov 
818ca4445c7SIlya Fursov   for (PetscInt i = 0; i < event->nevents; i++) { // check for brackets on the left/right of 't'
819ca4445c7SIlya Fursov     if (event->side[i] != 0) event->side[i] = TSEventTestBracket(event->fsign_prev[i], event->fsign[i], event->fsign_right[i], event->direction[i], event->iterctr);
820ca4445c7SIlya Fursov     minsidein = PetscMin(minsidein, event->side[i]);
821ca4445c7SIlya Fursov   }
822462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(&minsidein, &minsideout, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)ts)));
823ca4445c7SIlya Fursov   /*
824ca4445c7SIlya Fursov     minsideout (sync on all ranks) indicates the minimum of the following states:
825ca4445c7SIlya Fursov     -1 : [ptime_prev, t] is a bracket for some indicator-function-i
826ca4445c7SIlya Fursov     +1 : [t, ptime_right] is a bracket for some indicator-function-i
827ca4445c7SIlya Fursov      0 : t is a zero-crossing for some indicator-function-i
828ca4445c7SIlya Fursov      2 : none of the above
829ca4445c7SIlya Fursov   */
830ca4445c7SIlya Fursov   PetscCheck(!event->revisit_right || minsideout == 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "minsideout != 0 when performing 'revisiting' in TSEventHandler()");
831ca4445c7SIlya Fursov 
832ca4445c7SIlya Fursov   if (minsideout == -1 || minsideout == +1) {                                                           // this if-branch will refine the left/right bracket
833ca4445c7SIlya Fursov     const PetscReal bracket_size = (minsideout == -1) ? t - event->ptime_prev : event->ptime_right - t; // sync on all ranks
834ca4445c7SIlya Fursov 
835ca4445c7SIlya Fursov     if (minsideout == +1 && bracket_size <= event->timestep_min) { // check if the bracket (right) is small
836ca4445c7SIlya Fursov       // [--------------------|-]
837ca4445c7SIlya Fursov       dt_next              = bracket_size; // need one more step to get to event->ptime_right
838ca4445c7SIlya Fursov       event->revisit_right = PETSC_TRUE;
839ca4445c7SIlya Fursov       TSEvent_update_left(event, t);
840ca4445c7SIlya Fursov       if (event->monitor)
841ca4445c7SIlya Fursov         PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - reached too small bracket [%g - %g], next stepping to its right end %g (revisiting)\n", PetscGlobalRank, event->iterctr, (double)event->ptime_prev,
842ca4445c7SIlya Fursov                                          (double)event->ptime_right, (double)(event->ptime_prev + dt_next)));
843ca4445c7SIlya Fursov     } else { // the bracket is not very small -> refine it
844ca4445c7SIlya Fursov       // [--------|-------------]
845ca4445c7SIlya Fursov       if (bracket_size <= 2 * event->timestep_min) dt_next = bracket_size / 2; // the bracket is almost small -> bisect it
846ca4445c7SIlya Fursov       else {                                                                   // the bracket is not small -> use Anderson-Bjorck
847ca4445c7SIlya Fursov         PetscReal dti_min = PETSC_MAX_REAL;
848ca4445c7SIlya Fursov         for (PetscInt i = 0; i < event->nevents; i++) {
849ca4445c7SIlya Fursov           if (event->side[i] == minsideout) { // only refine the appropriate brackets
850ca4445c7SIlya Fursov             PetscReal dti = RefineAndersonBjorck(event->ptime_prev, t, event->ptime_right, event->fvalue_prev[i], event->fvalue[i], event->fvalue_right[i], event->side[i], &event->side_prev[i], event->justrefined_AB[i], &event->gamma_AB[i]);
851ca4445c7SIlya Fursov             dti_min       = PetscMin(dti_min, dti);
852ca4445c7SIlya Fursov           }
853ca4445c7SIlya Fursov         }
854462c564dSBarry Smith         PetscCallMPI(MPIU_Allreduce(&dti_min, &dt_next, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
855ca4445c7SIlya Fursov         if (dt_next < event->timestep_min) dt_next = event->timestep_min;
856ca4445c7SIlya Fursov         if (bracket_size - dt_next < event->timestep_min) dt_next = bracket_size - event->timestep_min;
857ca4445c7SIlya Fursov       }
858ca4445c7SIlya Fursov 
859ca4445c7SIlya Fursov       if (minsideout == -1) { // minsideout == -1, update the right-end values, retain the left-end values
860ca4445c7SIlya Fursov         TSEvent_update_right(event, t);
861ca4445c7SIlya Fursov         PetscCall(TSRollBack(ts));
862ca4445c7SIlya Fursov         PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); // e.g. to override TS_CONVERGED_TIME on reaching ts->max_time
863ca4445c7SIlya Fursov       } else TSEvent_update_left(event, t);                          // minsideout == +1, update the left-end values, retain the right-end values
864ca4445c7SIlya Fursov 
865ca4445c7SIlya Fursov       for (PetscInt i = 0; i < event->nevents; i++) { // update the "Anderson-Bjorck" flags
866ca4445c7SIlya Fursov         if (event->side[i] == minsideout) {
867ca4445c7SIlya Fursov           event->justrefined_AB[i] = PETSC_TRUE; // only for these i's Anderson-Bjorck was invoked
868ca4445c7SIlya Fursov           if (event->monitor)
869ca4445c7SIlya Fursov             PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " refining the bracket with sign change [%g - %g], next stepping to %g\n", PetscGlobalRank, event->iterctr, i, (double)event->ptime_prev,
870ca4445c7SIlya Fursov                                              (double)event->ptime_right, (double)(event->ptime_prev + dt_next)));
871ca4445c7SIlya Fursov         } else event->justrefined_AB[i] = PETSC_FALSE; // for these i's Anderson-Bjorck was not invoked
872ca4445c7SIlya Fursov       }
873ca4445c7SIlya Fursov     }
874ca4445c7SIlya Fursov     event->iterctr++;
875ca4445c7SIlya Fursov     event->processing = PETSC_TRUE;
876ca4445c7SIlya Fursov   } else if (minsideout == 0) { // found the appropriate zero-crossing (and no brackets to the left), finishing!
877ca4445c7SIlya Fursov     // [--------0-------------]
878ca4445c7SIlya Fursov     finished             = PETSC_TRUE;
879ca4445c7SIlya Fursov     event->revisit_right = PETSC_FALSE;
880ca4445c7SIlya Fursov     for (PetscInt i = 0; i < event->nevents; i++)
881ca4445c7SIlya Fursov       if (event->side[i] == minsideout) {
882ca4445c7SIlya Fursov         event->events_zero[event->nevents_zero++] = i;
883ca4445c7SIlya Fursov         if (event->fsign[i] == 0) { // vtol was engaged
884ca4445c7SIlya Fursov           if (event->monitor)
885ca4445c7SIlya Fursov             PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " zero crossing located at time %g (tol=%g)\n", PetscGlobalRank, event->iterctr, i, (double)t, (double)event->vtol[i]));
886ca4445c7SIlya Fursov         } else {               // dt_min was engaged
887ca4445c7SIlya Fursov           event->fsign[i] = 0; // sign = 0 is enforced further
888ca4445c7SIlya Fursov           if (event->monitor)
889ca4445c7SIlya Fursov             PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " accepting time %g as event location, due to reaching too small bracket [%g - %g]\n", PetscGlobalRank, event->iterctr, i, (double)t,
890ca4445c7SIlya Fursov                                              (double)event->ptime_prev, (double)t));
891ca4445c7SIlya Fursov         }
892ca4445c7SIlya Fursov       }
893ca4445c7SIlya Fursov     event->iterctr++;
894ca4445c7SIlya Fursov     event->processing = PETSC_TRUE;
895ca4445c7SIlya Fursov   } else { // minsideout == 2: no brackets, no zero-crossings
896ca4445c7SIlya Fursov     // [----------------------]
897ca4445c7SIlya Fursov     PetscCheck(event->iterctr == 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Unexpected state (event->iterctr != 0) in TSEventHandler()");
898fe4ad979SIlya Fursov     if (event->processing) {
899fe4ad979SIlya Fursov       PetscReal dt2;
900fe4ad979SIlya Fursov       if (event->timestep_2nd_postevent == PETSC_DECIDE) dt2 = event->timestep_cache;                            // (1)
901fe4ad979SIlya Fursov       else dt2 = event->timestep_2nd_postevent;                                                                  // (2), (2a), (3)
902fe4ad979SIlya Fursov       PetscCall(TSSetTimeStep(ts, TSEvent_dt_cap(ts, t, dt2, Not_PETSC_DECIDE(event->timestep_2nd_postevent)))); // set the second post-event step
903fe4ad979SIlya Fursov     }
904ca4445c7SIlya Fursov     event->processing = PETSC_FALSE;
905ca4445c7SIlya Fursov   }
906ca4445c7SIlya Fursov 
907ca4445c7SIlya Fursov   // if 'revisit_right' was flagged before the current iteration started, the iteration is expected to finish
908ca4445c7SIlya Fursov   PetscCheck(!revisit_right_cache || finished, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Unexpected state of 'revisit_right_cache' in TSEventHandler()");
909ca4445c7SIlya Fursov 
910ca4445c7SIlya Fursov   if (finished) { // finished handling the current event
911ca4445c7SIlya Fursov     PetscCall(TSPostEvent(ts, t, U));
912ca4445c7SIlya Fursov 
913fe4ad979SIlya Fursov     PetscReal dt1;
914fe4ad979SIlya Fursov     if (event->timestep_postevent == PETSC_DECIDE) { // (1), (2)
915fe4ad979SIlya Fursov       dt1               = event->ptime_cache - t;
916ca4445c7SIlya Fursov       event->processing = PETSC_TRUE;
917fe4ad979SIlya Fursov       if (PetscAbsReal(dt1) < PETSC_SMALL) { // (1a), (2a): the cached post-event point == event point
918fe4ad979SIlya Fursov         dt1               = event->timestep_cache;
919fe4ad979SIlya Fursov         event->processing = Not_PETSC_DECIDE(event->timestep_2nd_postevent);
920ca4445c7SIlya Fursov       }
921fe4ad979SIlya Fursov     } else {                                         // (3), (4)
922fe4ad979SIlya Fursov       dt1               = event->timestep_postevent; // 1st post-event dt = user-provided value
923fe4ad979SIlya Fursov       event->processing = Not_PETSC_DECIDE(event->timestep_2nd_postevent);
924ca4445c7SIlya Fursov     }
925fe4ad979SIlya Fursov 
926fe4ad979SIlya Fursov     PetscCall(TSSetTimeStep(ts, TSEvent_dt_cap(ts, t, dt1, Not_PETSC_DECIDE(event->timestep_postevent)))); // set the first post-event step
927ca4445c7SIlya Fursov     event->iterctr = 0;
928ca4445c7SIlya Fursov   } // if-finished
929ca4445c7SIlya Fursov 
930ca4445c7SIlya Fursov   if (event->iterctr == 0) TSEvent_update_left(event, t); // not found an event, or finished the event
931ca4445c7SIlya Fursov   else {
932ca4445c7SIlya Fursov     PetscCall(TSGetTime(ts, &t));                                              // update 't' to account for potential rollback
933ca4445c7SIlya Fursov     PetscCall(TSSetTimeStep(ts, TSEvent_dt_cap(ts, t, dt_next, PETSC_FALSE))); // continue resolving the event
93438bf2713SShri Abhyankar   }
9353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9362dc7a7e3SShri Abhyankar }
9372dc7a7e3SShri Abhyankar 
TSAdjointEventHandler(TS ts)938d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointEventHandler(TS ts)
939d71ae5a4SJacob Faibussowitsch {
9406427ac75SLisandro Dalcin   TSEvent   event;
941d0578d90SShri Abhyankar   PetscReal t;
942d0578d90SShri Abhyankar   Vec       U;
943d0578d90SShri Abhyankar   PetscInt  ctr;
944ca4445c7SIlya Fursov   PetscBool forwardsolve = PETSC_FALSE; // Flag indicating that TS is doing an adjoint solve
945d0578d90SShri Abhyankar 
946d0578d90SShri Abhyankar   PetscFunctionBegin;
9476427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
9483ba16761SJacob Faibussowitsch   if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS);
9496427ac75SLisandro Dalcin   event = ts->event;
950d0578d90SShri Abhyankar 
9519566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
9529566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
953d0578d90SShri Abhyankar 
954d0578d90SShri Abhyankar   ctr = event->recorder.ctr - 1;
955bcbf8bb3SShri Abhyankar   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
956ca4445c7SIlya Fursov     // Call the user post-event function
957d0578d90SShri Abhyankar     if (event->postevent) {
958ca4445c7SIlya Fursov       PetscCallBack("TSEvent post-event processing", (*event->postevent)(ts, event->recorder.nevents[ctr], event->recorder.eventidx[ctr], t, U, forwardsolve, event->ctx));
959d0578d90SShri Abhyankar       event->recorder.ctr--;
960d0578d90SShri Abhyankar     }
961d0578d90SShri Abhyankar   }
9623ba16761SJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)ts));
9633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
964d0578d90SShri Abhyankar }
9651ea83e56SMiguel 
9661ea83e56SMiguel /*@
967ca4445c7SIlya Fursov   TSGetNumEvents - Get the number of events defined on the given MPI process
9681ea83e56SMiguel 
9691ea83e56SMiguel   Logically Collective
9701ea83e56SMiguel 
9714165533cSJose E. Roman   Input Parameter:
972bcf0153eSBarry Smith . ts - the `TS` context
9731ea83e56SMiguel 
9744165533cSJose E. Roman   Output Parameter:
975ca4445c7SIlya Fursov . nevents - the number of local events on each MPI process
9761ea83e56SMiguel 
9771ea83e56SMiguel   Level: intermediate
9781ea83e56SMiguel 
9791cc06b55SBarry Smith .seealso: [](ch_ts), `TSEvent`, `TSSetEventHandler()`
9801ea83e56SMiguel @*/
TSGetNumEvents(TS ts,PetscInt * nevents)981d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetNumEvents(TS ts, PetscInt *nevents)
982d71ae5a4SJacob Faibussowitsch {
9831ea83e56SMiguel   PetscFunctionBegin;
9841ea83e56SMiguel   *nevents = ts->event->nevents;
9853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9861ea83e56SMiguel }
987