xref: /petsc/src/ts/event/tsevent.c (revision 9c334d8fdb557fc53fd345d68cbb3545b09ccab8)
1 #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "TSEventInitialize"
5 /*
6   TSEventInitialize - Initializes TSEvent for TSSolve
7 */
8 PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U)
9 {
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   if (!event) PetscFunctionReturn(0);
14   PetscValidPointer(event,1);
15   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
16   PetscValidHeaderSpecific(U,VEC_CLASSID,4);
17   event->ptime_prev = t;
18   event->iterctr = 0;
19   ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr);
20   /* Initialize the event recorder */
21   event->recorder.ctr = 0;
22   PetscFunctionReturn(0);
23 }
24 
25 #undef __FUNCT__
26 #define __FUNCT__ "TSEventDestroy"
27 PetscErrorCode TSEventDestroy(TSEvent *event)
28 {
29   PetscErrorCode ierr;
30   PetscInt       i;
31 
32   PetscFunctionBegin;
33   PetscValidPointer(event,1);
34   if (!*event) PetscFunctionReturn(0);
35   ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr);
36   ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr);
37   ierr = PetscFree((*event)->fvalue_right);CHKERRQ(ierr);
38   ierr = PetscFree((*event)->zerocrossing);CHKERRQ(ierr);
39   ierr = PetscFree((*event)->side);CHKERRQ(ierr);
40   ierr = PetscFree((*event)->direction);CHKERRQ(ierr);
41   ierr = PetscFree((*event)->terminate);CHKERRQ(ierr);
42   ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr);
43   ierr = PetscFree((*event)->vtol);CHKERRQ(ierr);
44 
45   for (i=0; i < (*event)->recsize; i++) {
46     ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr);
47   }
48   ierr = PetscFree((*event)->recorder.eventidx);CHKERRQ(ierr);
49   ierr = PetscFree((*event)->recorder.nevents);CHKERRQ(ierr);
50   ierr = PetscFree((*event)->recorder.stepnum);CHKERRQ(ierr);
51   ierr = PetscFree((*event)->recorder.time);CHKERRQ(ierr);
52 
53   ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr);
54   ierr = PetscFree(*event);CHKERRQ(ierr);
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "TSSetEventTolerances"
60 /*@
61    TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
62 
63    Logically Collective
64 
65    Input Arguments:
66 +  ts - time integration context
67 .  tol - scalar tolerance, PETSC_DECIDE to leave current value
68 -  vtol - array of tolerances or NULL, used in preference to tol if present
69 
70    Options Database Keys:
71 .  -ts_event_tol <tol> tolerance for event zero crossing
72 
73    Notes:
74    Must call TSSetEventHandler() before setting the tolerances.
75 
76    The size of vtol is equal to the number of events.
77 
78    Level: beginner
79 
80 .seealso: TS, TSEvent, TSSetEventHandler()
81 @*/
82 PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
83 {
84   TSEvent        event;
85   PetscInt       i;
86 
87   PetscFunctionBegin;
88   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
89   if (vtol) PetscValidRealPointer(vtol,3);
90   if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
91 
92   event = ts->event;
93   if (vtol) {
94     for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
95   } else {
96     if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
97       for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
98     }
99   }
100   PetscFunctionReturn(0);
101 }
102 
103 #undef __FUNCT__
104 #define __FUNCT__ "TSSetEventHandler"
105 /*@C
106    TSSetEventHandler - Sets a monitoring function used for detecting events
107 
108    Logically Collective on TS
109 
110    Input Parameters:
111 +  ts - the TS context obtained from TSCreate()
112 .  nevents - number of local events
113 .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
114                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
115 .  terminate - flag to indicate whether time stepping should be terminated after
116                event is detected (one for each event)
117 .  eventhandler - event monitoring routine
118 .  postevent - [optional] post-event function
119 -  ctx       - [optional] user-defined context for private data for the
120                event monitor and post event routine (use NULL if no
121                context is desired)
122 
123    Calling sequence of eventhandler:
124    PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
125 
126    Input Parameters:
127 +  ts  - the TS context
128 .  t   - current time
129 .  U   - current iterate
130 -  ctx - [optional] context passed with eventhandler
131 
132    Output parameters:
133 .  fvalue    - function value of events at time t
134 
135    Calling sequence of postevent:
136    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
137 
138    Input Parameters:
139 +  ts - the TS context
140 .  nevents_zero - number of local events whose event function is zero
141 .  events_zero  - indices of local events which have reached zero
142 .  t            - current time
143 .  U            - current solution
144 .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
145 -  ctx          - the context passed with eventhandler
146 
147    Level: intermediate
148 
149 .keywords: TS, event, set
150 
151 .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
152 @*/
153 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)
154 {
155   PetscErrorCode ierr;
156   TSEvent        event;
157   PetscInt       i;
158   PetscBool      flg;
159 #if defined PETSC_USE_REAL_SINGLE
160   PetscReal      tol=1e-4;
161 #else
162   PetscReal      tol=1e-6;
163 #endif
164 
165   PetscFunctionBegin;
166   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
167   PetscValidIntPointer(direction,2);
168   PetscValidIntPointer(terminate,3);
169 
170   ierr = PetscNewLog(ts,&event);CHKERRQ(ierr);
171   ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr);
172   ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr);
173   ierr = PetscMalloc1(nevents,&event->fvalue_right);CHKERRQ(ierr);
174   ierr = PetscMalloc1(nevents,&event->zerocrossing);CHKERRQ(ierr);
175   ierr = PetscMalloc1(nevents,&event->side);CHKERRQ(ierr);
176   ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr);
177   ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr);
178   ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr);
179   for (i=0; i < nevents; i++) {
180     event->direction[i] = direction[i];
181     event->terminate[i] = terminate[i];
182     event->zerocrossing[i] = PETSC_FALSE;
183     event->side[i] = 0;
184   }
185   ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr);
186   event->nevents = nevents;
187   event->eventhandler = eventhandler;
188   event->postevent = postevent;
189   event->ctx = ctx;
190 
191   event->recsize = 8;  /* Initial size of the recorder */
192   ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr);
193   {
194     ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL);CHKERRQ(ierr);
195     ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr);
196     ierr = PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);CHKERRQ(ierr);
197   }
198   PetscOptionsEnd();
199 
200   ierr = PetscMalloc1(event->recsize,&event->recorder.time);CHKERRQ(ierr);
201   ierr = PetscMalloc1(event->recsize,&event->recorder.stepnum);CHKERRQ(ierr);
202   ierr = PetscMalloc1(event->recsize,&event->recorder.nevents);CHKERRQ(ierr);
203   ierr = PetscMalloc1(event->recsize,&event->recorder.eventidx);CHKERRQ(ierr);
204   for (i=0; i < event->recsize; i++) {
205     ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr);
206   }
207 
208   for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
209   if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);}
210 
211   ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr);
212   ts->event = event;
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "TSEventRecorderResize"
218 /*
219   TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
220                           is reached.
221 */
222 static PetscErrorCode TSEventRecorderResize(TSEvent event)
223 {
224   PetscErrorCode ierr;
225   PetscReal      *time;
226   PetscInt       *stepnum;
227   PetscInt       *nevents;
228   PetscInt       **eventidx;
229   PetscInt       i,fact=2;
230 
231   PetscFunctionBegin;
232 
233   /* Create large arrays */
234   ierr = PetscMalloc1(fact*event->recsize,&time);CHKERRQ(ierr);
235   ierr = PetscMalloc1(fact*event->recsize,&stepnum);CHKERRQ(ierr);
236   ierr = PetscMalloc1(fact*event->recsize,&nevents);CHKERRQ(ierr);
237   ierr = PetscMalloc1(fact*event->recsize,&eventidx);CHKERRQ(ierr);
238   for (i=0; i < fact*event->recsize; i++) {
239     ierr = PetscMalloc1(event->nevents,&eventidx[i]);CHKERRQ(ierr);
240   }
241 
242   /* Copy over data */
243   ierr = PetscMemcpy(time,event->recorder.time,event->recsize*sizeof(PetscReal));CHKERRQ(ierr);
244   ierr = PetscMemcpy(stepnum,event->recorder.stepnum,event->recsize*sizeof(PetscInt));CHKERRQ(ierr);
245   ierr = PetscMemcpy(nevents,event->recorder.nevents,event->recsize*sizeof(PetscInt));CHKERRQ(ierr);
246   for (i=0; i < event->recsize; i++) {
247     ierr = PetscMemcpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]*sizeof(PetscInt));CHKERRQ(ierr);
248   }
249 
250   /* Destroy old arrays */
251   for (i=0; i < event->recsize; i++) {
252     ierr = PetscFree(event->recorder.eventidx[i]);CHKERRQ(ierr);
253   }
254   ierr = PetscFree(event->recorder.eventidx);CHKERRQ(ierr);
255   ierr = PetscFree(event->recorder.nevents);CHKERRQ(ierr);
256   ierr = PetscFree(event->recorder.stepnum);CHKERRQ(ierr);
257   ierr = PetscFree(event->recorder.time);CHKERRQ(ierr);
258 
259   /* Set pointers */
260   event->recorder.time = time;
261   event->recorder.stepnum = stepnum;
262   event->recorder.nevents = nevents;
263   event->recorder.eventidx = eventidx;
264 
265   /* Double size */
266   event->recsize *= fact;
267 
268   PetscFunctionReturn(0);
269 }
270 
271 /*
272    Helper rutine to handle user postenvents and recording
273 */
274 #undef __FUNCT__
275 #define __FUNCT__ "TSPostEvent"
276 static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
277 {
278   PetscErrorCode ierr;
279   TSEvent        event = ts->event;
280   PetscBool      terminate = PETSC_FALSE;
281   PetscInt       i,ctr,stepnum;
282   PetscBool      ts_terminate;
283   PetscBool      forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
284 
285   PetscFunctionBegin;
286   if (event->postevent) {
287     ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
288   }
289 
290   /* Handle termination events */
291   for (i=0; i < event->nevents_zero; i++) {
292     terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]);
293   }
294   ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr);
295   if (ts_terminate) {
296     ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);
297     event->status = TSEVENT_NONE;
298   } else {
299     event->status = TSEVENT_RESET_NEXTSTEP;
300   }
301 
302   event->ptime_prev = t;
303   /* Reset event residual functions as states might get changed by the postevent callback */
304   if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);}
305   /* Cache current event residual functions */
306   for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
307 
308   /* Record the event in the event recorder */
309   ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr);
310   ctr = event->recorder.ctr;
311   if (ctr == event->recsize) {
312     ierr = TSEventRecorderResize(event);CHKERRQ(ierr);
313   }
314   event->recorder.time[ctr] = t;
315   event->recorder.stepnum[ctr] = stepnum;
316   event->recorder.nevents[ctr] = event->nevents_zero;
317   for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
318   event->recorder.ctr++;
319   PetscFunctionReturn(0);
320 }
321 
322 /* Uses Anderson-Bjorck variant of regula falsi method */
323 PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt)
324 {
325   PetscReal new_dt, scal = 1.0;
326   if (PetscRealPart(fleft)*PetscRealPart(f) < 0) {
327     if (side == 1) {
328       scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
329       if (scal < PETSC_SMALL) scal = 0.5;
330     }
331     new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
332   } else {
333     if (side == -1) {
334       scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
335       if (scal < PETSC_SMALL) scal = 0.5;
336     }
337     new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t;
338   }
339   return PetscMin(dt,new_dt);
340 }
341 
342 
343 #undef __FUNCT__
344 #define __FUNCT__ "TSEventHandler"
345 PetscErrorCode TSEventHandler(TS ts)
346 {
347   PetscErrorCode ierr;
348   TSEvent        event;
349   PetscReal      t;
350   Vec            U;
351   PetscInt       i;
352   PetscReal      dt,dt_min;
353   PetscInt       rollback=0,in[2],out[2];
354   PetscInt       fvalue_sign,fvalueprev_sign;
355 
356   PetscFunctionBegin;
357   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
358   if (!ts->event) PetscFunctionReturn(0);
359   event = ts->event;
360 
361   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
362   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
363   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
364 
365   if (event->status == TSEVENT_NONE) {
366     event->timestep_prev = dt;
367   }
368 
369   if (event->status == TSEVENT_RESET_NEXTSTEP) {
370     /* Restore previous time step */
371     dt = event->timestep_prev;
372     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
373     event->status = TSEVENT_NONE;
374   }
375 
376   if (event->status == TSEVENT_NONE) {
377     event->ptime_end = t;
378   }
379 
380   ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);
381 
382   for (i=0; i < event->nevents; i++) {
383     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
384       event->status = TSEVENT_ZERO;
385       event->fvalue_right[i] = event->fvalue[i];
386       continue;
387     }
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) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
391       switch (event->direction[i]) {
392       case -1:
393         if (fvalue_sign < 0) {
394           rollback = 1;
395 
396           /* Compute new time step */
397           dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
398 
399           if (event->monitor) {
400             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
401           }
402           event->fvalue_right[i] = event->fvalue[i];
403           event->side[i] = 1;
404 
405           if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
406           event->status = TSEVENT_LOCATED_INTERVAL;
407         }
408         break;
409       case 1:
410         if (fvalue_sign > 0) {
411           rollback = 1;
412 
413           /* Compute new time step */
414           dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
415 
416           if (event->monitor) {
417             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
418           }
419           event->fvalue_right[i] = event->fvalue[i];
420           event->side[i] = 1;
421 
422           if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
423           event->status = TSEVENT_LOCATED_INTERVAL;
424         }
425         break;
426       case 0:
427         rollback = 1;
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 
432         if (event->monitor) {
433           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
434         }
435         event->fvalue_right[i] = event->fvalue[i];
436         event->side[i] = 1;
437 
438         if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
439         event->status = TSEVENT_LOCATED_INTERVAL;
440 
441         break;
442       }
443     }
444   }
445 
446   in[0] = event->status; in[1] = rollback;
447   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
448   event->status = (TSEventStatus)out[0]; rollback = out[1];
449   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
450 
451   event->nevents_zero = 0;
452   if (event->status == TSEVENT_ZERO) {
453     for (i=0; i < event->nevents; i++) {
454       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
455         event->events_zero[event->nevents_zero++] = i;
456         if (event->monitor) {
457           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);CHKERRQ(ierr);
458         }
459         event->zerocrossing[i] = PETSC_FALSE;
460       }
461       event->side[i] = 0;
462     }
463     ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr);
464 
465     dt = event->ptime_end - t;
466     if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */
467     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
468     event->iterctr = 0;
469     PetscFunctionReturn(0);
470   }
471 
472   if (event->status == TSEVENT_LOCATED_INTERVAL) {
473     ierr = TSRollBack(ts);CHKERRQ(ierr);
474     ts->steps--;
475     ts->total_steps--;
476     ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr);
477     event->status = TSEVENT_PROCESSING;
478     event->ptime_right = t;
479   } else {
480     for (i=0; i < event->nevents; i++) {
481       if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) {
482         /* Compute new time step */
483         dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
484         event->side[i] = -1;
485       }
486       event->fvalue_prev[i] = event->fvalue[i];
487     }
488     if (event->monitor && event->status == TSEVENT_PROCESSING) {
489       ierr = 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);CHKERRQ(ierr);
490     }
491     event->ptime_prev = t;
492   }
493 
494   if (event->status == TSEVENT_PROCESSING) event->iterctr++;
495 
496   ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
497   ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 
501 #undef __FUNCT__
502 #define __FUNCT__ "TSAdjointEventHandler"
503 PetscErrorCode TSAdjointEventHandler(TS ts)
504 {
505   PetscErrorCode ierr;
506   TSEvent        event;
507   PetscReal      t;
508   Vec            U;
509   PetscInt       ctr;
510   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
511 
512   PetscFunctionBegin;
513   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
514   if (!ts->event) PetscFunctionReturn(0);
515   event = ts->event;
516 
517   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
518   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
519 
520   ctr = event->recorder.ctr-1;
521   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
522     /* Call the user postevent function */
523     if (event->postevent) {
524       ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
525       event->recorder.ctr--;
526     }
527   }
528 
529   PetscBarrier((PetscObject)ts);
530   PetscFunctionReturn(0);
531 }
532