xref: /petsc/src/ts/event/tsevent.c (revision 862e4a309d45a165aaa4da0d704ba733429d833a)
1 #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
2 
3 /*
4   TSEventInitialize - Initializes TSEvent for TSSolve
5 */
6 PetscErrorCode TSEventInitialize(TSEvent event, TS ts, PetscReal t, Vec U)
7 {
8   PetscFunctionBegin;
9   if (!event) PetscFunctionReturn(0);
10   PetscValidPointer(event, 1);
11   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
12   PetscValidHeaderSpecific(U, VEC_CLASSID, 4);
13   event->ptime_prev = t;
14   event->iterctr    = 0;
15   PetscCall((*event->eventhandler)(ts, t, U, event->fvalue_prev, event->ctx));
16   PetscFunctionReturn(0);
17 }
18 
19 PetscErrorCode TSEventDestroy(TSEvent *event)
20 {
21   PetscInt i;
22 
23   PetscFunctionBegin;
24   PetscValidPointer(event, 1);
25   if (!*event) PetscFunctionReturn(0);
26   if (--(*event)->refct > 0) {
27     *event = NULL;
28     PetscFunctionReturn(0);
29   }
30 
31   PetscCall(PetscFree((*event)->fvalue));
32   PetscCall(PetscFree((*event)->fvalue_prev));
33   PetscCall(PetscFree((*event)->fvalue_right));
34   PetscCall(PetscFree((*event)->zerocrossing));
35   PetscCall(PetscFree((*event)->side));
36   PetscCall(PetscFree((*event)->direction));
37   PetscCall(PetscFree((*event)->terminate));
38   PetscCall(PetscFree((*event)->events_zero));
39   PetscCall(PetscFree((*event)->vtol));
40 
41   for (i = 0; i < (*event)->recsize; i++) PetscCall(PetscFree((*event)->recorder.eventidx[i]));
42   PetscCall(PetscFree((*event)->recorder.eventidx));
43   PetscCall(PetscFree((*event)->recorder.nevents));
44   PetscCall(PetscFree((*event)->recorder.stepnum));
45   PetscCall(PetscFree((*event)->recorder.time));
46 
47   PetscCall(PetscViewerDestroy(&(*event)->monitor));
48   PetscCall(PetscFree(*event));
49   PetscFunctionReturn(0);
50 }
51 
52 /*@
53   TSSetPostEventIntervalStep - Set the time-step used immediately following the event interval
54 
55   Logically Collective
56 
57   Input Parameters:
58 + ts - time integration context
59 - dt - post event interval step
60 
61   Options Database Keys:
62 . -ts_event_post_eventinterval_step <dt> time-step after event interval
63 
64   Notes:
65   TSSetPostEventIntervalStep allows one to set a time-step that is used immediately following an event interval.
66 
67   This function should be called from the postevent function set with TSSetEventHandler().
68 
69   The post event interval time-step should be selected based on the dynamics following the event.
70   If the dynamics are stiff, a conservative (small) step should be used.
71   If not, then a larger time-step can be used.
72 
73   Level: Advanced
74   .seealso: `TS`, `TSEvent`, `TSSetEventHandler()`
75 @*/
76 PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
77 {
78   PetscFunctionBegin;
79   ts->event->timestep_posteventinterval = dt;
80   PetscFunctionReturn(0);
81 }
82 
83 /*@
84    TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
85 
86    Logically Collective
87 
88    Input Parameters:
89 +  ts - time integration context
90 .  tol - scalar tolerance, PETSC_DECIDE to leave current value
91 -  vtol - array of tolerances or NULL, used in preference to tol if present
92 
93    Options Database Keys:
94 .  -ts_event_tol <tol> - tolerance for event zero crossing
95 
96    Notes:
97    Must call TSSetEventHandler() before setting the tolerances.
98 
99    The size of vtol is equal to the number of events.
100 
101    Level: beginner
102 
103 .seealso: `TS`, `TSEvent`, `TSSetEventHandler()`
104 @*/
105 PetscErrorCode TSSetEventTolerances(TS ts, PetscReal tol, PetscReal vtol[])
106 {
107   TSEvent  event;
108   PetscInt i;
109 
110   PetscFunctionBegin;
111   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
112   if (vtol) PetscValidRealPointer(vtol, 3);
113   PetscCheck(ts->event, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set the events first by calling TSSetEventHandler()");
114 
115   event = ts->event;
116   if (vtol) {
117     for (i = 0; i < event->nevents; i++) event->vtol[i] = vtol[i];
118   } else {
119     if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
120       for (i = 0; i < event->nevents; i++) event->vtol[i] = tol;
121     }
122   }
123   PetscFunctionReturn(0);
124 }
125 
126 /*@C
127    TSSetEventHandler - Sets a function used for detecting events
128 
129    Logically Collective on TS
130 
131    Input Parameters:
132 +  ts - the TS context obtained from TSCreate()
133 .  nevents - number of local events
134 .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
135                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
136 .  terminate - flag to indicate whether time stepping should be terminated after
137                event is detected (one for each event)
138 .  eventhandler - event monitoring routine
139 .  postevent - [optional] post-event function
140 -  ctx       - [optional] user-defined context for private data for the
141                event monitor and post event routine (use NULL if no
142                context is desired)
143 
144    Calling sequence of eventhandler:
145    PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
146 
147    Input Parameters:
148 +  ts  - the TS context
149 .  t   - current time
150 .  U   - current iterate
151 -  ctx - [optional] context passed with eventhandler
152 
153    Output parameters:
154 .  fvalue    - function value of events at time t
155 
156    Calling sequence of postevent:
157    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
158 
159    Input Parameters:
160 +  ts - the TS context
161 .  nevents_zero - number of local events whose event function is zero
162 .  events_zero  - indices of local events which have reached zero
163 .  t            - current time
164 .  U            - current solution
165 .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
166 -  ctx          - the context passed with eventhandler
167 
168    Level: intermediate
169 
170 .seealso: `TSCreate()`, `TSSetTimeStep()`, `TSSetConvergedReason()`
171 @*/
172 PetscErrorCode TSSetEventHandler(TS ts, PetscInt nevents, PetscInt direction[], PetscBool terminate[], PetscErrorCode (*eventhandler)(TS, PetscReal, Vec, PetscScalar[], void *), PetscErrorCode (*postevent)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *ctx)
173 {
174   TSAdapt   adapt;
175   PetscReal hmin;
176   TSEvent   event;
177   PetscInt  i;
178   PetscBool flg;
179 #if defined PETSC_USE_REAL_SINGLE
180   PetscReal tol = 1e-4;
181 #else
182   PetscReal tol = 1e-6;
183 #endif
184 
185   PetscFunctionBegin;
186   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
187   if (nevents) {
188     PetscValidIntPointer(direction, 3);
189     PetscValidBoolPointer(terminate, 4);
190   }
191 
192   PetscCall(PetscNew(&event));
193   PetscCall(PetscMalloc1(nevents, &event->fvalue));
194   PetscCall(PetscMalloc1(nevents, &event->fvalue_prev));
195   PetscCall(PetscMalloc1(nevents, &event->fvalue_right));
196   PetscCall(PetscMalloc1(nevents, &event->zerocrossing));
197   PetscCall(PetscMalloc1(nevents, &event->side));
198   PetscCall(PetscMalloc1(nevents, &event->direction));
199   PetscCall(PetscMalloc1(nevents, &event->terminate));
200   PetscCall(PetscMalloc1(nevents, &event->vtol));
201   for (i = 0; i < nevents; i++) {
202     event->direction[i]    = direction[i];
203     event->terminate[i]    = terminate[i];
204     event->zerocrossing[i] = PETSC_FALSE;
205     event->side[i]         = 0;
206   }
207   PetscCall(PetscMalloc1(nevents, &event->events_zero));
208   event->nevents                    = nevents;
209   event->eventhandler               = eventhandler;
210   event->postevent                  = postevent;
211   event->ctx                        = ctx;
212   event->timestep_posteventinterval = ts->time_step;
213   PetscCall(TSGetAdapt(ts, &adapt));
214   PetscCall(TSAdaptGetStepLimits(adapt, &hmin, NULL));
215   event->timestep_min = hmin;
216 
217   event->recsize = 8; /* Initial size of the recorder */
218   PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "TS Event options", "TS");
219   {
220     PetscCall(PetscOptionsReal("-ts_event_tol", "Scalar event tolerance for zero crossing check", "TSSetEventTolerances", tol, &tol, NULL));
221     PetscCall(PetscOptionsName("-ts_event_monitor", "Print choices made by event handler", "", &flg));
222     PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size", "Initial size of event recorder", "", event->recsize, &event->recsize, NULL));
223     PetscCall(PetscOptionsReal("-ts_event_post_eventinterval_step", "Time step after event interval", "", event->timestep_posteventinterval, &event->timestep_posteventinterval, NULL));
224     PetscCall(PetscOptionsReal("-ts_event_post_event_step", "Time step after event", "", event->timestep_postevent, &event->timestep_postevent, NULL));
225     PetscCall(PetscOptionsReal("-ts_event_dt_min", "Minimum time step considered for TSEvent", "", event->timestep_min, &event->timestep_min, NULL));
226   }
227   PetscOptionsEnd();
228 
229   PetscCall(PetscMalloc1(event->recsize, &event->recorder.time));
230   PetscCall(PetscMalloc1(event->recsize, &event->recorder.stepnum));
231   PetscCall(PetscMalloc1(event->recsize, &event->recorder.nevents));
232   PetscCall(PetscMalloc1(event->recsize, &event->recorder.eventidx));
233   for (i = 0; i < event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &event->recorder.eventidx[i]));
234   /* Initialize the event recorder */
235   event->recorder.ctr = 0;
236 
237   for (i = 0; i < event->nevents; i++) event->vtol[i] = tol;
238   if (flg) PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF, "stdout", &event->monitor));
239 
240   PetscCall(TSEventDestroy(&ts->event));
241   ts->event        = event;
242   ts->event->refct = 1;
243   PetscFunctionReturn(0);
244 }
245 
246 /*
247   TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
248                           is reached.
249 */
250 static PetscErrorCode TSEventRecorderResize(TSEvent event)
251 {
252   PetscReal *time;
253   PetscInt  *stepnum;
254   PetscInt  *nevents;
255   PetscInt **eventidx;
256   PetscInt   i, fact = 2;
257 
258   PetscFunctionBegin;
259 
260   /* Create large arrays */
261   PetscCall(PetscMalloc1(fact * event->recsize, &time));
262   PetscCall(PetscMalloc1(fact * event->recsize, &stepnum));
263   PetscCall(PetscMalloc1(fact * event->recsize, &nevents));
264   PetscCall(PetscMalloc1(fact * event->recsize, &eventidx));
265   for (i = 0; i < fact * event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &eventidx[i]));
266 
267   /* Copy over data */
268   PetscCall(PetscArraycpy(time, event->recorder.time, event->recsize));
269   PetscCall(PetscArraycpy(stepnum, event->recorder.stepnum, event->recsize));
270   PetscCall(PetscArraycpy(nevents, event->recorder.nevents, event->recsize));
271   for (i = 0; i < event->recsize; i++) PetscCall(PetscArraycpy(eventidx[i], event->recorder.eventidx[i], event->recorder.nevents[i]));
272 
273   /* Destroy old arrays */
274   for (i = 0; i < event->recsize; i++) PetscCall(PetscFree(event->recorder.eventidx[i]));
275   PetscCall(PetscFree(event->recorder.eventidx));
276   PetscCall(PetscFree(event->recorder.nevents));
277   PetscCall(PetscFree(event->recorder.stepnum));
278   PetscCall(PetscFree(event->recorder.time));
279 
280   /* Set pointers */
281   event->recorder.time     = time;
282   event->recorder.stepnum  = stepnum;
283   event->recorder.nevents  = nevents;
284   event->recorder.eventidx = eventidx;
285 
286   /* Double size */
287   event->recsize *= fact;
288 
289   PetscFunctionReturn(0);
290 }
291 
292 /*
293    Helper routine to handle user postevents and recording
294 */
295 static PetscErrorCode TSPostEvent(TS ts, PetscReal t, Vec U)
296 {
297   TSEvent   event     = ts->event;
298   PetscBool terminate = PETSC_FALSE;
299   PetscBool restart   = PETSC_FALSE;
300   PetscInt  i, ctr, stepnum;
301   PetscBool inflag[2], outflag[2];
302   PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
303 
304   PetscFunctionBegin;
305   if (event->postevent) {
306     PetscObjectState state_prev, state_post;
307     PetscCall(PetscObjectStateGet((PetscObject)U, &state_prev));
308     PetscCall((*event->postevent)(ts, event->nevents_zero, event->events_zero, t, U, forwardsolve, event->ctx));
309     PetscCall(PetscObjectStateGet((PetscObject)U, &state_post));
310     if (state_prev != state_post) restart = PETSC_TRUE;
311   }
312 
313   /* Handle termination events and step restart */
314   for (i = 0; i < event->nevents_zero; i++)
315     if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
316   inflag[0] = restart;
317   inflag[1] = terminate;
318   PetscCall(MPIU_Allreduce(inflag, outflag, 2, MPIU_BOOL, MPI_LOR, ((PetscObject)ts)->comm));
319   restart   = outflag[0];
320   terminate = outflag[1];
321   if (restart) PetscCall(TSRestartStep(ts));
322   if (terminate) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_EVENT));
323   event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
324 
325   /* Reset event residual functions as states might get changed by the postevent callback */
326   if (event->postevent) {
327     PetscCall(VecLockReadPush(U));
328     PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx));
329     PetscCall(VecLockReadPop(U));
330   }
331 
332   /* Cache current time and event residual functions */
333   event->ptime_prev = t;
334   for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
335 
336   /* Record the event in the event recorder */
337   PetscCall(TSGetStepNumber(ts, &stepnum));
338   ctr = event->recorder.ctr;
339   if (ctr == event->recsize) PetscCall(TSEventRecorderResize(event));
340   event->recorder.time[ctr]    = t;
341   event->recorder.stepnum[ctr] = stepnum;
342   event->recorder.nevents[ctr] = event->nevents_zero;
343   for (i = 0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
344   event->recorder.ctr++;
345   PetscFunctionReturn(0);
346 }
347 
348 /* Uses Anderson-Bjorck variant of regula falsi method */
349 static inline PetscReal TSEventComputeStepSize(PetscReal tleft, PetscReal t, PetscReal tright, PetscScalar fleft, PetscScalar f, PetscScalar fright, PetscInt side, PetscReal dt)
350 {
351   PetscReal new_dt, scal = 1.0;
352   if (PetscRealPart(fleft) * PetscRealPart(f) < 0) {
353     if (side == 1) {
354       scal = (PetscRealPart(fright) - PetscRealPart(f)) / PetscRealPart(fright);
355       if (scal < PETSC_SMALL) scal = 0.5;
356     }
357     new_dt = (scal * PetscRealPart(fleft) * t - PetscRealPart(f) * tleft) / (scal * PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
358   } else {
359     if (side == -1) {
360       scal = (PetscRealPart(fleft) - PetscRealPart(f)) / PetscRealPart(fleft);
361       if (scal < PETSC_SMALL) scal = 0.5;
362     }
363     new_dt = (PetscRealPart(f) * tright - scal * PetscRealPart(fright) * t) / (PetscRealPart(f) - scal * PetscRealPart(fright)) - t;
364   }
365   return PetscMin(dt, new_dt);
366 }
367 
368 static PetscErrorCode TSEventDetection(TS ts)
369 {
370   TSEvent   event = ts->event;
371   PetscReal t;
372   PetscInt  i;
373   PetscInt  fvalue_sign, fvalueprev_sign;
374   PetscInt  in, out;
375 
376   PetscFunctionBegin;
377   PetscCall(TSGetTime(ts, &t));
378   for (i = 0; i < event->nevents; i++) {
379     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
380       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
381       event->status = TSEVENT_LOCATED_INTERVAL;
382       if (event->monitor) {
383         PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " interval detected due to zero value (tol=%g) [%g - %g]\n", event->iterctr, i, (double)event->vtol[i], (double)event->ptime_prev, (double)t));
384       }
385       continue;
386     }
387     if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */
388     fvalue_sign     = PetscSign(PetscRealPart(event->fvalue[i]));
389     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
390     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) {
391       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
392       event->status = TSEVENT_LOCATED_INTERVAL;
393       if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " interval detected due to sign change [%g - %g]\n", event->iterctr, i, (double)event->ptime_prev, (double)t));
394     }
395   }
396   in = (PetscInt)event->status;
397   PetscCall(MPIU_Allreduce(&in, &out, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts)));
398   event->status = (TSEventStatus)out;
399   PetscFunctionReturn(0);
400 }
401 
402 static PetscErrorCode TSEventLocation(TS ts, PetscReal *dt)
403 {
404   TSEvent   event = ts->event;
405   PetscInt  i;
406   PetscReal t;
407   PetscInt  fvalue_sign, fvalueprev_sign;
408   PetscInt  rollback = 0, in[2], out[2];
409 
410   PetscFunctionBegin;
411   PetscCall(TSGetTime(ts, &t));
412   event->nevents_zero = 0;
413   for (i = 0; i < event->nevents; i++) {
414     if (event->zerocrossing[i]) {
415       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i] || *dt < event->timestep_min || PetscAbsReal((*dt) / ((event->ptime_right - event->ptime_prev) / 2)) < event->vtol[i]) { /* stopping criteria */
416         event->status          = TSEVENT_ZERO;
417         event->fvalue_right[i] = event->fvalue[i];
418         continue;
419       }
420       /* Compute new time step */
421       *dt             = TSEventComputeStepSize(event->ptime_prev, t, event->ptime_right, event->fvalue_prev[i], event->fvalue[i], event->fvalue_right[i], event->side[i], *dt);
422       fvalue_sign     = PetscSign(PetscRealPart(event->fvalue[i]));
423       fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
424       switch (event->direction[i]) {
425       case -1:
426         if (fvalue_sign < 0) {
427           rollback               = 1;
428           event->fvalue_right[i] = event->fvalue[i];
429           event->side[i]         = 1;
430         }
431         break;
432       case 1:
433         if (fvalue_sign > 0) {
434           rollback               = 1;
435           event->fvalue_right[i] = event->fvalue[i];
436           event->side[i]         = 1;
437         }
438         break;
439       case 0:
440         if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */
441           rollback               = 1;
442           event->fvalue_right[i] = event->fvalue[i];
443           event->side[i]         = 1;
444         }
445         break;
446       }
447       if (event->status == TSEVENT_PROCESSING) event->side[i] = -1;
448     }
449   }
450   in[0] = (PetscInt)event->status;
451   in[1] = rollback;
452   PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts)));
453   event->status = (TSEventStatus)out[0];
454   rollback      = out[1];
455   /* If rollback is true, the status will be overwritten so that an event at the endtime of current time step will be postponed to guarantee corret order */
456   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
457   if (event->status == TSEVENT_ZERO) {
458     for (i = 0; i < event->nevents; i++) {
459       if (event->zerocrossing[i]) {
460         if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i] || *dt < event->timestep_min || PetscAbsReal((*dt) / ((event->ptime_right - event->ptime_prev) / 2)) < event->vtol[i]) { /* stopping criteria */
461           event->events_zero[event->nevents_zero++] = i;
462           if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " zero crossing located at time %g\n", event->iterctr, i, (double)t));
463           event->zerocrossing[i] = PETSC_FALSE;
464         }
465       }
466       event->side[i] = 0;
467     }
468   }
469   PetscFunctionReturn(0);
470 }
471 
472 PetscErrorCode TSEventHandler(TS ts)
473 {
474   TSEvent   event;
475   PetscReal t;
476   Vec       U;
477   PetscInt  i;
478   PetscReal dt, dt_min, dt_reset = 0.0;
479 
480   PetscFunctionBegin;
481   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
482   if (!ts->event) PetscFunctionReturn(0);
483   event = ts->event;
484 
485   PetscCall(TSGetTime(ts, &t));
486   PetscCall(TSGetTimeStep(ts, &dt));
487   PetscCall(TSGetSolution(ts, &U));
488 
489   if (event->status == TSEVENT_NONE) {
490     event->timestep_prev = dt;
491     event->ptime_end     = t;
492   }
493   if (event->status == TSEVENT_RESET_NEXTSTEP) {
494     /* user has specified a PostEventInterval dt */
495     dt = event->timestep_posteventinterval;
496     if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
497       PetscReal maxdt = ts->max_time - t;
498       dt              = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt);
499     }
500     PetscCall(TSSetTimeStep(ts, dt));
501     event->status = TSEVENT_NONE;
502   }
503 
504   PetscCall(VecLockReadPush(U));
505   PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx));
506   PetscCall(VecLockReadPop(U));
507 
508   /* Detect the events */
509   PetscCall(TSEventDetection(ts));
510 
511   /* Locate the events */
512   if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) {
513     /* Approach the zero crosing by setting a new step size */
514     PetscCall(TSEventLocation(ts, &dt));
515     /* Roll back when new events are detected */
516     if (event->status == TSEVENT_LOCATED_INTERVAL) {
517       PetscCall(TSRollBack(ts));
518       PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING));
519       event->iterctr++;
520     }
521     PetscCall(MPIU_Allreduce(&dt, &dt_min, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
522     if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset;
523     PetscCall(TSSetTimeStep(ts, dt_min));
524     /* Found the zero crossing */
525     if (event->status == TSEVENT_ZERO) {
526       PetscCall(TSPostEvent(ts, t, U));
527 
528       dt = event->ptime_end - t;
529       if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */
530         dt            = event->timestep_prev;
531         event->status = TSEVENT_NONE;
532       }
533       if (event->timestep_postevent) { /* user has specified a PostEvent dt*/
534         dt = event->timestep_postevent;
535       }
536       if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
537         PetscReal maxdt = ts->max_time - t;
538         dt              = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt);
539       }
540       PetscCall(TSSetTimeStep(ts, dt));
541       event->iterctr = 0;
542     }
543     /* Have not found the zero crosing yet */
544     if (event->status == TSEVENT_PROCESSING) {
545       if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Stepping forward as no event detected in interval [%g - %g]\n", event->iterctr, (double)event->ptime_prev, (double)t));
546       event->iterctr++;
547     }
548   }
549   if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */
550     event->status      = TSEVENT_PROCESSING;
551     event->ptime_right = t;
552   } else {
553     for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
554     event->ptime_prev = t;
555   }
556   PetscFunctionReturn(0);
557 }
558 
559 PetscErrorCode TSAdjointEventHandler(TS ts)
560 {
561   TSEvent   event;
562   PetscReal t;
563   Vec       U;
564   PetscInt  ctr;
565   PetscBool forwardsolve = PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
566 
567   PetscFunctionBegin;
568   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
569   if (!ts->event) PetscFunctionReturn(0);
570   event = ts->event;
571 
572   PetscCall(TSGetTime(ts, &t));
573   PetscCall(TSGetSolution(ts, &U));
574 
575   ctr = event->recorder.ctr - 1;
576   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
577     /* Call the user postevent function */
578     if (event->postevent) {
579       PetscCall((*event->postevent)(ts, event->recorder.nevents[ctr], event->recorder.eventidx[ctr], t, U, forwardsolve, event->ctx));
580       event->recorder.ctr--;
581     }
582   }
583 
584   PetscBarrier((PetscObject)ts);
585   PetscFunctionReturn(0);
586 }
587 
588 /*@
589   TSGetNumEvents - Get the numbers of events set
590 
591   Logically Collective
592 
593   Input Parameter:
594 . ts - the TS context
595 
596   Output Parameter:
597 . nevents - number of events
598 
599   Level: intermediate
600 
601 .seealso: `TSSetEventHandler()`
602 
603 @*/
604 PetscErrorCode TSGetNumEvents(TS ts, PetscInt *nevents)
605 {
606   PetscFunctionBegin;
607   *nevents = ts->event->nevents;
608   PetscFunctionReturn(0);
609 }
610