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