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