xref: /petsc/src/ts/event/tsevent.c (revision 0369a392ee39ea4d01e519e42e389892d4e1e2e4)
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   PetscErrorCode ierr;
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(PetscNewLog(ts,&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   ierr = PetscOptionsBegin(((PetscObject)ts)->comm,((PetscObject)ts)->prefix,"TS Event options","TS");PetscCall(ierr);
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   ierr = PetscOptionsEnd();PetscCall(ierr);
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++) {
234     PetscCall(PetscMalloc1(event->nevents,&event->recorder.eventidx[i]));
235   }
236   /* Initialize the event recorder */
237   event->recorder.ctr = 0;
238 
239   for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
240   if (flg) PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor));
241 
242   PetscCall(TSEventDestroy(&ts->event));
243   ts->event = event;
244   ts->event->refct = 1;
245   PetscFunctionReturn(0);
246 }
247 
248 /*
249   TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
250                           is reached.
251 */
252 static PetscErrorCode TSEventRecorderResize(TSEvent event)
253 {
254   PetscReal      *time;
255   PetscInt       *stepnum;
256   PetscInt       *nevents;
257   PetscInt       **eventidx;
258   PetscInt       i,fact=2;
259 
260   PetscFunctionBegin;
261 
262   /* Create large arrays */
263   PetscCall(PetscMalloc1(fact*event->recsize,&time));
264   PetscCall(PetscMalloc1(fact*event->recsize,&stepnum));
265   PetscCall(PetscMalloc1(fact*event->recsize,&nevents));
266   PetscCall(PetscMalloc1(fact*event->recsize,&eventidx));
267   for (i=0; i < fact*event->recsize; i++) {
268     PetscCall(PetscMalloc1(event->nevents,&eventidx[i]));
269   }
270 
271   /* Copy over data */
272   PetscCall(PetscArraycpy(time,event->recorder.time,event->recsize));
273   PetscCall(PetscArraycpy(stepnum,event->recorder.stepnum,event->recsize));
274   PetscCall(PetscArraycpy(nevents,event->recorder.nevents,event->recsize));
275   for (i=0; i < event->recsize; i++) {
276     PetscCall(PetscArraycpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]));
277   }
278 
279   /* Destroy old arrays */
280   for (i=0; i < event->recsize; i++) {
281     PetscCall(PetscFree(event->recorder.eventidx[i]));
282   }
283   PetscCall(PetscFree(event->recorder.eventidx));
284   PetscCall(PetscFree(event->recorder.nevents));
285   PetscCall(PetscFree(event->recorder.stepnum));
286   PetscCall(PetscFree(event->recorder.time));
287 
288   /* Set pointers */
289   event->recorder.time = time;
290   event->recorder.stepnum = stepnum;
291   event->recorder.nevents = nevents;
292   event->recorder.eventidx = eventidx;
293 
294   /* Double size */
295   event->recsize *= fact;
296 
297   PetscFunctionReturn(0);
298 }
299 
300 /*
301    Helper routine to handle user postevents and recording
302 */
303 static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
304 {
305   TSEvent        event = ts->event;
306   PetscBool      terminate = PETSC_FALSE;
307   PetscBool      restart = PETSC_FALSE;
308   PetscInt       i,ctr,stepnum;
309   PetscBool      inflag[2],outflag[2];
310   PetscBool      forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
311 
312   PetscFunctionBegin;
313   if (event->postevent) {
314     PetscObjectState state_prev,state_post;
315     PetscCall(PetscObjectStateGet((PetscObject)U,&state_prev));
316     PetscCall((*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx));
317     PetscCall(PetscObjectStateGet((PetscObject)U,&state_post));
318     if (state_prev != state_post) restart = PETSC_TRUE;
319   }
320 
321   /* Handle termination events and step restart */
322   for (i=0; i<event->nevents_zero; i++) if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
323   inflag[0] = restart; inflag[1] = terminate;
324   PetscCall(MPIU_Allreduce(inflag,outflag,2,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm));
325   restart = outflag[0]; terminate = outflag[1];
326   if (restart) PetscCall(TSRestartStep(ts));
327   if (terminate) PetscCall(TSSetConvergedReason(ts,TS_CONVERGED_EVENT));
328   event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
329 
330   /* Reset event residual functions as states might get changed by the postevent callback */
331   if (event->postevent) {
332     PetscCall(VecLockReadPush(U));
333     PetscCall((*event->eventhandler)(ts,t,U,event->fvalue,event->ctx));
334     PetscCall(VecLockReadPop(U));
335   }
336 
337   /* Cache current time and event residual functions */
338   event->ptime_prev = t;
339   for (i=0; i<event->nevents; i++)
340     event->fvalue_prev[i] = event->fvalue[i];
341 
342   /* Record the event in the event recorder */
343   PetscCall(TSGetStepNumber(ts,&stepnum));
344   ctr = event->recorder.ctr;
345   if (ctr == event->recsize) {
346     PetscCall(TSEventRecorderResize(event));
347   }
348   event->recorder.time[ctr] = t;
349   event->recorder.stepnum[ctr] = stepnum;
350   event->recorder.nevents[ctr] = event->nevents_zero;
351   for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
352   event->recorder.ctr++;
353   PetscFunctionReturn(0);
354 }
355 
356 /* Uses Anderson-Bjorck variant of regula falsi method */
357 static inline PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt)
358 {
359   PetscReal new_dt, scal = 1.0;
360   if (PetscRealPart(fleft)*PetscRealPart(f) < 0) {
361     if (side == 1) {
362       scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
363       if (scal < PETSC_SMALL) scal = 0.5;
364     }
365     new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
366   } else {
367     if (side == -1) {
368       scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
369       if (scal < PETSC_SMALL) scal = 0.5;
370     }
371     new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t;
372   }
373   return PetscMin(dt,new_dt);
374 }
375 
376 static PetscErrorCode TSEventDetection(TS ts)
377 {
378   TSEvent        event = ts->event;
379   PetscReal      t;
380   PetscInt       i;
381   PetscInt       fvalue_sign,fvalueprev_sign;
382   PetscInt       in,out;
383 
384   PetscFunctionBegin;
385   PetscCall(TSGetTime(ts,&t));
386   for (i=0; i < event->nevents; i++) {
387     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
388       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
389       event->status = TSEVENT_LOCATED_INTERVAL;
390       if (event->monitor) {
391         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));
392       }
393       continue;
394     }
395     if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */
396     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
397     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
398     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) {
399       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
400       event->status = TSEVENT_LOCATED_INTERVAL;
401       if (event->monitor) {
402         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));
403       }
404     }
405   }
406   in = (PetscInt)event->status;
407   PetscCall(MPIU_Allreduce(&in,&out,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts)));
408   event->status = (TSEventStatus)out;
409   PetscFunctionReturn(0);
410 }
411 
412 static PetscErrorCode TSEventLocation(TS ts,PetscReal *dt)
413 {
414   TSEvent        event = ts->event;
415   PetscInt       i;
416   PetscReal      t;
417   PetscInt       fvalue_sign,fvalueprev_sign;
418   PetscInt       rollback=0,in[2],out[2];
419 
420   PetscFunctionBegin;
421   PetscCall(TSGetTime(ts,&t));
422   event->nevents_zero = 0;
423   for (i=0; i < event->nevents; i++) {
424     if (event->zerocrossing[i]) {
425       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 */
426         event->status = TSEVENT_ZERO;
427         event->fvalue_right[i] = event->fvalue[i];
428         continue;
429       }
430       /* Compute new time step */
431       *dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],*dt);
432       fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
433       fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
434       switch (event->direction[i]) {
435       case -1:
436         if (fvalue_sign < 0) {
437           rollback = 1;
438           event->fvalue_right[i] = event->fvalue[i];
439           event->side[i] = 1;
440         }
441         break;
442       case 1:
443         if (fvalue_sign > 0) {
444           rollback = 1;
445           event->fvalue_right[i] = event->fvalue[i];
446           event->side[i] = 1;
447         }
448         break;
449       case 0:
450         if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */
451           rollback = 1;
452           event->fvalue_right[i] = event->fvalue[i];
453           event->side[i] = 1;
454         }
455         break;
456       }
457       if (event->status == TSEVENT_PROCESSING) event->side[i] = -1;
458     }
459   }
460   in[0] = (PetscInt)event->status; in[1] = rollback;
461   PetscCall(MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts)));
462   event->status = (TSEventStatus)out[0]; rollback = out[1];
463   /* 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 */
464   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
465   if (event->status == TSEVENT_ZERO) {
466     for (i=0; i < event->nevents; i++) {
467       if (event->zerocrossing[i]) {
468         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 */
469           event->events_zero[event->nevents_zero++] = i;
470           if (event->monitor) {
471             PetscCall(PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D zero crossing located at time %g\n",event->iterctr,i,(double)t));
472           }
473           event->zerocrossing[i] = PETSC_FALSE;
474         }
475       }
476       event->side[i] = 0;
477     }
478   }
479   PetscFunctionReturn(0);
480 }
481 
482 PetscErrorCode TSEventHandler(TS ts)
483 {
484   TSEvent        event;
485   PetscReal      t;
486   Vec            U;
487   PetscInt       i;
488   PetscReal      dt,dt_min,dt_reset = 0.0;
489 
490   PetscFunctionBegin;
491   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
492   if (!ts->event) PetscFunctionReturn(0);
493   event = ts->event;
494 
495   PetscCall(TSGetTime(ts,&t));
496   PetscCall(TSGetTimeStep(ts,&dt));
497   PetscCall(TSGetSolution(ts,&U));
498 
499   if (event->status == TSEVENT_NONE) {
500     event->timestep_prev = dt;
501     event->ptime_end = t;
502   }
503   if (event->status == TSEVENT_RESET_NEXTSTEP) {
504     /* user has specified a PostEventInterval dt */
505     dt = event->timestep_posteventinterval;
506     if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
507       PetscReal maxdt = ts->max_time-t;
508       dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt);
509     }
510     PetscCall(TSSetTimeStep(ts,dt));
511     event->status = TSEVENT_NONE;
512   }
513 
514   PetscCall(VecLockReadPush(U));
515   PetscCall((*event->eventhandler)(ts,t,U,event->fvalue,event->ctx));
516   PetscCall(VecLockReadPop(U));
517 
518   /* Detect the events */
519   PetscCall(TSEventDetection(ts));
520 
521   /* Locate the events */
522   if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) {
523     /* Approach the zero crosing by setting a new step size */
524     PetscCall(TSEventLocation(ts,&dt));
525     /* Roll back when new events are detected */
526     if (event->status == TSEVENT_LOCATED_INTERVAL) {
527       PetscCall(TSRollBack(ts));
528       PetscCall(TSSetConvergedReason(ts,TS_CONVERGED_ITERATING));
529       event->iterctr++;
530     }
531     PetscCall(MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts)));
532     if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset;
533     PetscCall(TSSetTimeStep(ts,dt_min));
534     /* Found the zero crossing */
535     if (event->status == TSEVENT_ZERO) {
536       PetscCall(TSPostEvent(ts,t,U));
537 
538       dt = event->ptime_end - t;
539       if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */
540         dt = event->timestep_prev;
541         event->status = TSEVENT_NONE;
542       }
543       if (event->timestep_postevent) { /* user has specified a PostEvent dt*/
544         dt = event->timestep_postevent;
545       }
546       if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
547         PetscReal maxdt = ts->max_time-t;
548         dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt);
549       }
550       PetscCall(TSSetTimeStep(ts,dt));
551       event->iterctr = 0;
552     }
553     /* Have not found the zero crosing yet */
554     if (event->status == TSEVENT_PROCESSING) {
555       if (event->monitor) {
556         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));
557       }
558       event->iterctr++;
559     }
560   }
561   if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */
562     event->status = TSEVENT_PROCESSING;
563     event->ptime_right = t;
564   } else {
565     for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
566     event->ptime_prev = t;
567   }
568   PetscFunctionReturn(0);
569 }
570 
571 PetscErrorCode TSAdjointEventHandler(TS ts)
572 {
573   TSEvent        event;
574   PetscReal      t;
575   Vec            U;
576   PetscInt       ctr;
577   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
578 
579   PetscFunctionBegin;
580   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
581   if (!ts->event) PetscFunctionReturn(0);
582   event = ts->event;
583 
584   PetscCall(TSGetTime(ts,&t));
585   PetscCall(TSGetSolution(ts,&U));
586 
587   ctr = event->recorder.ctr-1;
588   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
589     /* Call the user postevent function */
590     if (event->postevent) {
591       PetscCall((*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx));
592       event->recorder.ctr--;
593     }
594   }
595 
596   PetscBarrier((PetscObject)ts);
597   PetscFunctionReturn(0);
598 }
599 
600 /*@
601   TSGetNumEvents - Get the numbers of events set
602 
603   Logically Collective
604 
605   Input Parameter:
606 . ts - the TS context
607 
608   Output Parameter:
609 . nevents - number of events
610 
611   Level: intermediate
612 
613 .seealso: TSSetEventHandler()
614 
615 @*/
616 PetscErrorCode TSGetNumEvents(TS ts,PetscInt * nevents)
617 {
618   PetscFunctionBegin;
619   *nevents = ts->event->nevents;
620   PetscFunctionReturn(0);
621 }
622