xref: /petsc/src/ts/event/tsevent.c (revision 1c2dc1cbabe5212f80cf309c4ca5a26f0cadc660)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
22dc7a7e3SShri Abhyankar 
36fea3669SShri Abhyankar /*
46427ac75SLisandro Dalcin   TSEventInitialize - Initializes TSEvent for TSSolve
56fea3669SShri Abhyankar */
66427ac75SLisandro Dalcin PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U)
76fea3669SShri Abhyankar {
86fea3669SShri Abhyankar   PetscFunctionBegin;
96427ac75SLisandro Dalcin   if (!event) PetscFunctionReturn(0);
106427ac75SLisandro Dalcin   PetscValidPointer(event,1);
116427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
126427ac75SLisandro Dalcin   PetscValidHeaderSpecific(U,VEC_CLASSID,4);
136fea3669SShri Abhyankar   event->ptime_prev = t;
1438bf2713SShri Abhyankar   event->iterctr = 0;
159566063dSJacob Faibussowitsch   PetscCall((*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx));
166fea3669SShri Abhyankar   PetscFunctionReturn(0);
176fea3669SShri Abhyankar }
186fea3669SShri Abhyankar 
197dbe0728SLisandro Dalcin PetscErrorCode TSEventDestroy(TSEvent *event)
207dbe0728SLisandro Dalcin {
217dbe0728SLisandro Dalcin   PetscInt       i;
227dbe0728SLisandro Dalcin 
237dbe0728SLisandro Dalcin   PetscFunctionBegin;
247dbe0728SLisandro Dalcin   PetscValidPointer(event,1);
257dbe0728SLisandro Dalcin   if (!*event) PetscFunctionReturn(0);
26c793f718SLisandro Dalcin   if (--(*event)->refct > 0) {*event = NULL; PetscFunctionReturn(0);}
27e7069c78SShri 
289566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue));
299566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue_prev));
309566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue_right));
319566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->zerocrossing));
329566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->side));
339566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->direction));
349566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->terminate));
359566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->events_zero));
369566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->vtol));
37a4ffd976SShri Abhyankar 
38a4ffd976SShri Abhyankar   for (i=0; i < (*event)->recsize; i++) {
399566063dSJacob Faibussowitsch     PetscCall(PetscFree((*event)->recorder.eventidx[i]));
407dbe0728SLisandro Dalcin   }
419566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.eventidx));
429566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.nevents));
439566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.stepnum));
449566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.time));
45a4ffd976SShri Abhyankar 
469566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*event)->monitor));
479566063dSJacob Faibussowitsch   PetscCall(PetscFree(*event));
487dbe0728SLisandro Dalcin   PetscFunctionReturn(0);
497dbe0728SLisandro Dalcin }
507dbe0728SLisandro Dalcin 
51e3005195SShri Abhyankar /*@
52ac6a796dSBarry Smith   TSSetPostEventIntervalStep - Set the time-step used immediately following the event interval
53458122a4SShri Abhyankar 
54458122a4SShri Abhyankar   Logically Collective
55458122a4SShri Abhyankar 
564165533cSJose E. Roman   Input Parameters:
57458122a4SShri Abhyankar + ts - time integration context
58458122a4SShri Abhyankar - dt - post event interval step
59458122a4SShri Abhyankar 
60458122a4SShri Abhyankar   Options Database Keys:
61458122a4SShri Abhyankar . -ts_event_post_eventinterval_step <dt> time-step after event interval
62458122a4SShri Abhyankar 
63458122a4SShri Abhyankar   Notes:
64ac6a796dSBarry Smith   TSSetPostEventIntervalStep allows one to set a time-step that is used immediately following an event interval.
65458122a4SShri Abhyankar 
66458122a4SShri Abhyankar   This function should be called from the postevent function set with TSSetEventHandler().
67458122a4SShri Abhyankar 
68458122a4SShri Abhyankar   The post event interval time-step should be selected based on the dynamics following the event.
69458122a4SShri Abhyankar   If the dynamics are stiff, a conservative (small) step should be used.
70458122a4SShri Abhyankar   If not, then a larger time-step can be used.
71458122a4SShri Abhyankar 
72458122a4SShri Abhyankar   Level: Advanced
73458122a4SShri Abhyankar   .seealso: TS, TSEvent, TSSetEventHandler()
74458122a4SShri Abhyankar @*/
75458122a4SShri Abhyankar PetscErrorCode TSSetPostEventIntervalStep(TS ts,PetscReal dt)
76458122a4SShri Abhyankar {
77458122a4SShri Abhyankar   PetscFunctionBegin;
78458122a4SShri Abhyankar   ts->event->timestep_posteventinterval = dt;
79458122a4SShri Abhyankar   PetscFunctionReturn(0);
80458122a4SShri Abhyankar }
81458122a4SShri Abhyankar 
82458122a4SShri Abhyankar /*@
83e3005195SShri Abhyankar    TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
84e3005195SShri Abhyankar 
85e3005195SShri Abhyankar    Logically Collective
86e3005195SShri Abhyankar 
874165533cSJose E. Roman    Input Parameters:
88e3005195SShri Abhyankar +  ts - time integration context
89e3005195SShri Abhyankar .  tol - scalar tolerance, PETSC_DECIDE to leave current value
90e3005195SShri Abhyankar -  vtol - array of tolerances or NULL, used in preference to tol if present
91e3005195SShri Abhyankar 
922589fa24SLisandro Dalcin    Options Database Keys:
93147403d9SBarry Smith .  -ts_event_tol <tol> - tolerance for event zero crossing
94e3005195SShri Abhyankar 
95e3005195SShri Abhyankar    Notes:
96f25fe674SLisandro Dalcin    Must call TSSetEventHandler() before setting the tolerances.
97e3005195SShri Abhyankar 
98e3005195SShri Abhyankar    The size of vtol is equal to the number of events.
99e3005195SShri Abhyankar 
100e3005195SShri Abhyankar    Level: beginner
101e3005195SShri Abhyankar 
102f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler()
103e3005195SShri Abhyankar @*/
1046427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
105e3005195SShri Abhyankar {
1066427ac75SLisandro Dalcin   TSEvent        event;
107e3005195SShri Abhyankar   PetscInt       i;
108e3005195SShri Abhyankar 
109e3005195SShri Abhyankar   PetscFunctionBegin;
1106427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1116427ac75SLisandro Dalcin   if (vtol) PetscValidRealPointer(vtol,3);
1123c633725SBarry Smith   PetscCheck(ts->event,PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
1136427ac75SLisandro Dalcin 
1146427ac75SLisandro Dalcin   event = ts->event;
115e3005195SShri Abhyankar   if (vtol) {
116e3005195SShri Abhyankar     for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
117e3005195SShri Abhyankar   } else {
118e3005195SShri Abhyankar     if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
119e3005195SShri Abhyankar       for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
120e3005195SShri Abhyankar     }
121e3005195SShri Abhyankar   }
122e3005195SShri Abhyankar   PetscFunctionReturn(0);
123e3005195SShri Abhyankar }
124e3005195SShri Abhyankar 
1252dc7a7e3SShri Abhyankar /*@C
126ac6a796dSBarry Smith    TSSetEventHandler - Sets a function used for detecting events
1272dc7a7e3SShri Abhyankar 
1282dc7a7e3SShri Abhyankar    Logically Collective on TS
1292dc7a7e3SShri Abhyankar 
1302dc7a7e3SShri Abhyankar    Input Parameters:
1312dc7a7e3SShri Abhyankar +  ts - the TS context obtained from TSCreate()
1322dc7a7e3SShri Abhyankar .  nevents - number of local events
133d94325d3SShri Abhyankar .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
134d94325d3SShri Abhyankar                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
135d94325d3SShri Abhyankar .  terminate - flag to indicate whether time stepping should be terminated after
136d94325d3SShri Abhyankar                event is detected (one for each event)
1376427ac75SLisandro Dalcin .  eventhandler - event monitoring routine
1382dc7a7e3SShri Abhyankar .  postevent - [optional] post-event function
1392589fa24SLisandro Dalcin -  ctx       - [optional] user-defined context for private data for the
1402dc7a7e3SShri Abhyankar                event monitor and post event routine (use NULL if no
1412dc7a7e3SShri Abhyankar                context is desired)
1422dc7a7e3SShri Abhyankar 
1436427ac75SLisandro Dalcin    Calling sequence of eventhandler:
1443a88037aSBarry Smith    PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
1452dc7a7e3SShri Abhyankar 
1462dc7a7e3SShri Abhyankar    Input Parameters:
1472dc7a7e3SShri Abhyankar +  ts  - the TS context
1482dc7a7e3SShri Abhyankar .  t   - current time
1492dc7a7e3SShri Abhyankar .  U   - current iterate
1506427ac75SLisandro Dalcin -  ctx - [optional] context passed with eventhandler
1512dc7a7e3SShri Abhyankar 
1522dc7a7e3SShri Abhyankar    Output parameters:
153d94325d3SShri Abhyankar .  fvalue    - function value of events at time t
1542dc7a7e3SShri Abhyankar 
1552dc7a7e3SShri Abhyankar    Calling sequence of postevent:
1562589fa24SLisandro Dalcin    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
1572dc7a7e3SShri Abhyankar 
1582dc7a7e3SShri Abhyankar    Input Parameters:
1592dc7a7e3SShri Abhyankar +  ts - the TS context
1602dc7a7e3SShri Abhyankar .  nevents_zero - number of local events whose event function is zero
1612dc7a7e3SShri Abhyankar .  events_zero  - indices of local events which have reached zero
1622dc7a7e3SShri Abhyankar .  t            - current time
1632dc7a7e3SShri Abhyankar .  U            - current solution
164031fbad4SShri Abhyankar .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
1656427ac75SLisandro Dalcin -  ctx          - the context passed with eventhandler
1662dc7a7e3SShri Abhyankar 
1672dc7a7e3SShri Abhyankar    Level: intermediate
1682dc7a7e3SShri Abhyankar 
1692dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
1702dc7a7e3SShri Abhyankar @*/
1712589fa24SLisandro Dalcin 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)
1722dc7a7e3SShri Abhyankar {
1732dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
174d294eb03SHong Zhang   TSAdapt        adapt;
175d294eb03SHong Zhang   PetscReal      hmin;
1762dc7a7e3SShri Abhyankar   TSEvent        event;
177d94325d3SShri Abhyankar   PetscInt       i;
178006e6a18SShri Abhyankar   PetscBool      flg;
179a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE
180a6c783d2SShri Abhyankar   PetscReal      tol=1e-4;
181a6c783d2SShri Abhyankar #else
182d569cc17SSatish Balay   PetscReal      tol=1e-6;
183a6c783d2SShri Abhyankar #endif
1842dc7a7e3SShri Abhyankar 
1852dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1866427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1870a82b154SShri   if (nevents) {
188064a246eSJacob Faibussowitsch     PetscValidIntPointer(direction,3);
189064a246eSJacob Faibussowitsch     PetscValidBoolPointer(terminate,4);
1900a82b154SShri   }
1916427ac75SLisandro Dalcin 
1929566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&event));
1939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->fvalue));
1949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->fvalue_prev));
1959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->fvalue_right));
1969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->zerocrossing));
1979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->side));
1989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->direction));
1999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->terminate));
2009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->vtol));
201d94325d3SShri Abhyankar   for (i=0; i < nevents; i++) {
202d94325d3SShri Abhyankar     event->direction[i] = direction[i];
203d94325d3SShri Abhyankar     event->terminate[i] = terminate[i];
204e2cdd850SShri Abhyankar     event->zerocrossing[i] = PETSC_FALSE;
205e2cdd850SShri Abhyankar     event->side[i] = 0;
206d94325d3SShri Abhyankar   }
2079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents,&event->events_zero));
2082589fa24SLisandro Dalcin   event->nevents = nevents;
2096427ac75SLisandro Dalcin   event->eventhandler = eventhandler;
2102dc7a7e3SShri Abhyankar   event->postevent = postevent;
2116427ac75SLisandro Dalcin   event->ctx = ctx;
212458122a4SShri Abhyankar   event->timestep_posteventinterval = ts->time_step;
2139566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&adapt));
2149566063dSJacob Faibussowitsch   PetscCall(TSAdaptGetStepLimits(adapt,&hmin,NULL));
215d294eb03SHong Zhang   event->timestep_min = hmin;
2162dc7a7e3SShri Abhyankar 
21702749585SLisandro Dalcin   event->recsize = 8;  /* Initial size of the recorder */
2189566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(((PetscObject)ts)->comm,((PetscObject)ts)->prefix,"TS Event options","TS");PetscCall(ierr);
2192dc7a7e3SShri Abhyankar   {
2209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL));
2219566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg));
2229566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL));
2239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_post_eventinterval_step","Time step after event interval","",event->timestep_posteventinterval,&event->timestep_posteventinterval,NULL));
2249566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_post_event_step","Time step after event","",event->timestep_postevent,&event->timestep_postevent,NULL));
2259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_dt_min","Minimum time step considered for TSEvent","",event->timestep_min,&event->timestep_min,NULL));
2262dc7a7e3SShri Abhyankar   }
2279566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
228a4ffd976SShri Abhyankar 
2299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize,&event->recorder.time));
2309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize,&event->recorder.stepnum));
2319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize,&event->recorder.nevents));
2329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize,&event->recorder.eventidx));
233a4ffd976SShri Abhyankar   for (i=0; i < event->recsize; i++) {
2349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(event->nevents,&event->recorder.eventidx[i]));
235a4ffd976SShri Abhyankar   }
236e7069c78SShri   /* Initialize the event recorder */
237e7069c78SShri   event->recorder.ctr = 0;
238a4ffd976SShri Abhyankar 
239e3005195SShri Abhyankar   for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
2409566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor));
2419e12be75SShri Abhyankar 
2429566063dSJacob Faibussowitsch   PetscCall(TSEventDestroy(&ts->event));
243d94325d3SShri Abhyankar   ts->event = event;
244e7069c78SShri   ts->event->refct = 1;
2452dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
2462dc7a7e3SShri Abhyankar }
2472dc7a7e3SShri Abhyankar 
248a4ffd976SShri Abhyankar /*
249a4ffd976SShri Abhyankar   TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
250a4ffd976SShri Abhyankar                           is reached.
251a4ffd976SShri Abhyankar */
25202749585SLisandro Dalcin static PetscErrorCode TSEventRecorderResize(TSEvent event)
253a4ffd976SShri Abhyankar {
254a4ffd976SShri Abhyankar   PetscReal      *time;
255a4ffd976SShri Abhyankar   PetscInt       *stepnum;
256a4ffd976SShri Abhyankar   PetscInt       *nevents;
257a4ffd976SShri Abhyankar   PetscInt       **eventidx;
258a4ffd976SShri Abhyankar   PetscInt       i,fact=2;
259a4ffd976SShri Abhyankar 
260a4ffd976SShri Abhyankar   PetscFunctionBegin;
261a4ffd976SShri Abhyankar 
262a4ffd976SShri Abhyankar   /* Create large arrays */
2639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact*event->recsize,&time));
2649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact*event->recsize,&stepnum));
2659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact*event->recsize,&nevents));
2669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fact*event->recsize,&eventidx));
267a4ffd976SShri Abhyankar   for (i=0; i < fact*event->recsize; i++) {
2689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(event->nevents,&eventidx[i]));
269a4ffd976SShri Abhyankar   }
270a4ffd976SShri Abhyankar 
271a4ffd976SShri Abhyankar   /* Copy over data */
2729566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(time,event->recorder.time,event->recsize));
2739566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(stepnum,event->recorder.stepnum,event->recsize));
2749566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(nevents,event->recorder.nevents,event->recsize));
275a4ffd976SShri Abhyankar   for (i=0; i < event->recsize; i++) {
2769566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]));
277a4ffd976SShri Abhyankar   }
278a4ffd976SShri Abhyankar 
279a4ffd976SShri Abhyankar   /* Destroy old arrays */
280a4ffd976SShri Abhyankar   for (i=0; i < event->recsize; i++) {
2819566063dSJacob Faibussowitsch     PetscCall(PetscFree(event->recorder.eventidx[i]));
282a4ffd976SShri Abhyankar   }
2839566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.eventidx));
2849566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.nevents));
2859566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.stepnum));
2869566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.time));
287a4ffd976SShri Abhyankar 
288a4ffd976SShri Abhyankar   /* Set pointers */
289a4ffd976SShri Abhyankar   event->recorder.time = time;
290a4ffd976SShri Abhyankar   event->recorder.stepnum = stepnum;
291a4ffd976SShri Abhyankar   event->recorder.nevents = nevents;
292a4ffd976SShri Abhyankar   event->recorder.eventidx = eventidx;
293a4ffd976SShri Abhyankar 
294a4ffd976SShri Abhyankar   /* Double size */
295a4ffd976SShri Abhyankar   event->recsize *= fact;
296a4ffd976SShri Abhyankar 
297a4ffd976SShri Abhyankar   PetscFunctionReturn(0);
298a4ffd976SShri Abhyankar }
299a4ffd976SShri Abhyankar 
300031fbad4SShri Abhyankar /*
301ac6a796dSBarry Smith    Helper routine to handle user postevents and recording
302031fbad4SShri Abhyankar */
3034597913aSLisandro Dalcin static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
3042dc7a7e3SShri Abhyankar {
3052dc7a7e3SShri Abhyankar   TSEvent        event = ts->event;
3062dc7a7e3SShri Abhyankar   PetscBool      terminate = PETSC_FALSE;
30728d5b5d6SLisandro Dalcin   PetscBool      restart = PETSC_FALSE;
308d0578d90SShri Abhyankar   PetscInt       i,ctr,stepnum;
3097324a0ffSLisandro Dalcin   PetscBool      inflag[2],outflag[2];
3104597913aSLisandro Dalcin   PetscBool      forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
3112dc7a7e3SShri Abhyankar 
3122dc7a7e3SShri Abhyankar   PetscFunctionBegin;
3132dc7a7e3SShri Abhyankar   if (event->postevent) {
31428d5b5d6SLisandro Dalcin     PetscObjectState state_prev,state_post;
3159566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateGet((PetscObject)U,&state_prev));
3169566063dSJacob Faibussowitsch     PetscCall((*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx));
3179566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateGet((PetscObject)U,&state_post));
31828d5b5d6SLisandro Dalcin     if (state_prev != state_post) restart = PETSC_TRUE;
3192dc7a7e3SShri Abhyankar   }
3204597913aSLisandro Dalcin 
32128d5b5d6SLisandro Dalcin   /* Handle termination events and step restart */
32228d5b5d6SLisandro Dalcin   for (i=0; i<event->nevents_zero; i++) if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
3237324a0ffSLisandro Dalcin   inflag[0] = restart; inflag[1] = terminate;
324*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(inflag,outflag,2,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm));
3257324a0ffSLisandro Dalcin   restart = outflag[0]; terminate = outflag[1];
3269566063dSJacob Faibussowitsch   if (restart) PetscCall(TSRestartStep(ts));
3279566063dSJacob Faibussowitsch   if (terminate) PetscCall(TSSetConvergedReason(ts,TS_CONVERGED_EVENT));
3287324a0ffSLisandro Dalcin   event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
329f7aea88cSShri Abhyankar 
3304597913aSLisandro Dalcin   /* Reset event residual functions as states might get changed by the postevent callback */
331f443add6SLisandro Dalcin   if (event->postevent) {
3329566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(U));
3339566063dSJacob Faibussowitsch     PetscCall((*event->eventhandler)(ts,t,U,event->fvalue,event->ctx));
3349566063dSJacob Faibussowitsch     PetscCall(VecLockReadPop(U));
335f443add6SLisandro Dalcin   }
336f443add6SLisandro Dalcin 
337f443add6SLisandro Dalcin   /* Cache current time and event residual functions */
338f443add6SLisandro Dalcin   event->ptime_prev = t;
339f443add6SLisandro Dalcin   for (i=0; i<event->nevents; i++)
340f443add6SLisandro Dalcin     event->fvalue_prev[i] = event->fvalue[i];
3414597913aSLisandro Dalcin 
342d0578d90SShri Abhyankar   /* Record the event in the event recorder */
3439566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&stepnum));
344f7aea88cSShri Abhyankar   ctr = event->recorder.ctr;
345a4ffd976SShri Abhyankar   if (ctr == event->recsize) {
3469566063dSJacob Faibussowitsch     PetscCall(TSEventRecorderResize(event));
347a4ffd976SShri Abhyankar   }
348f7aea88cSShri Abhyankar   event->recorder.time[ctr] = t;
349d0578d90SShri Abhyankar   event->recorder.stepnum[ctr] = stepnum;
3504597913aSLisandro Dalcin   event->recorder.nevents[ctr] = event->nevents_zero;
3514597913aSLisandro Dalcin   for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
352f7aea88cSShri Abhyankar   event->recorder.ctr++;
3532dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
3542dc7a7e3SShri Abhyankar }
3552dc7a7e3SShri Abhyankar 
35602749585SLisandro Dalcin /* Uses Anderson-Bjorck variant of regula falsi method */
3579fbee547SJacob Faibussowitsch static inline PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt)
35838bf2713SShri Abhyankar {
35902749585SLisandro Dalcin   PetscReal new_dt, scal = 1.0;
360e2cdd850SShri Abhyankar   if (PetscRealPart(fleft)*PetscRealPart(f) < 0) {
361e2cdd850SShri Abhyankar     if (side == 1) {
362a6c783d2SShri Abhyankar       scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
363a6c783d2SShri Abhyankar       if (scal < PETSC_SMALL) scal = 0.5;
364e2cdd850SShri Abhyankar     }
36502749585SLisandro Dalcin     new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
366e2cdd850SShri Abhyankar   } else {
367e2cdd850SShri Abhyankar     if (side == -1) {
368a6c783d2SShri Abhyankar       scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
369a6c783d2SShri Abhyankar       if (scal < PETSC_SMALL) scal = 0.5;
370e2cdd850SShri Abhyankar     }
37102749585SLisandro Dalcin     new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t;
372e2cdd850SShri Abhyankar   }
37302749585SLisandro Dalcin   return PetscMin(dt,new_dt);
37438bf2713SShri Abhyankar }
375e2cdd850SShri Abhyankar 
376d294eb03SHong Zhang static PetscErrorCode TSEventDetection(TS ts)
377d294eb03SHong Zhang {
378d294eb03SHong Zhang   TSEvent        event = ts->event;
379d294eb03SHong Zhang   PetscReal      t;
380d294eb03SHong Zhang   PetscInt       i;
381d294eb03SHong Zhang   PetscInt       fvalue_sign,fvalueprev_sign;
382ea3dac1cSHong Zhang   PetscInt       in,out;
383d294eb03SHong Zhang 
384d294eb03SHong Zhang   PetscFunctionBegin;
3859566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&t));
386d294eb03SHong Zhang   for (i=0; i < event->nevents; i++) {
387d294eb03SHong Zhang     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
388d294eb03SHong Zhang       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
389d294eb03SHong Zhang       event->status = TSEVENT_LOCATED_INTERVAL;
390d294eb03SHong Zhang       if (event->monitor) {
3919566063dSJacob Faibussowitsch         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));
392d294eb03SHong Zhang       }
393d294eb03SHong Zhang       continue;
394d294eb03SHong Zhang     }
395ea3dac1cSHong Zhang     if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */
396d294eb03SHong Zhang     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
397d294eb03SHong Zhang     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
398d294eb03SHong Zhang     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) {
399d294eb03SHong Zhang       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
400d294eb03SHong Zhang       event->status = TSEVENT_LOCATED_INTERVAL;
401d294eb03SHong Zhang       if (event->monitor) {
4029566063dSJacob Faibussowitsch         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));
403d294eb03SHong Zhang       }
404d294eb03SHong Zhang     }
405d294eb03SHong Zhang   }
406d5f990dbSBarry Smith   in = (PetscInt)event->status;
407*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&in,&out,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts)));
408ea3dac1cSHong Zhang   event->status = (TSEventStatus)out;
409d294eb03SHong Zhang   PetscFunctionReturn(0);
410d294eb03SHong Zhang }
411d294eb03SHong Zhang 
412d294eb03SHong Zhang static PetscErrorCode TSEventLocation(TS ts,PetscReal *dt)
413d294eb03SHong Zhang {
414d294eb03SHong Zhang   TSEvent        event = ts->event;
415d294eb03SHong Zhang   PetscInt       i;
416d294eb03SHong Zhang   PetscReal      t;
417ea3dac1cSHong Zhang   PetscInt       fvalue_sign,fvalueprev_sign;
418ea3dac1cSHong Zhang   PetscInt       rollback=0,in[2],out[2];
419d294eb03SHong Zhang 
420d294eb03SHong Zhang   PetscFunctionBegin;
4219566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&t));
422d294eb03SHong Zhang   event->nevents_zero = 0;
423d294eb03SHong Zhang   for (i=0; i < event->nevents; i++) {
424d294eb03SHong Zhang     if (event->zerocrossing[i]) {
425d294eb03SHong Zhang       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 */
426d294eb03SHong Zhang         event->status = TSEVENT_ZERO;
427d294eb03SHong Zhang         event->fvalue_right[i] = event->fvalue[i];
428d294eb03SHong Zhang         continue;
429d294eb03SHong Zhang       }
430d294eb03SHong Zhang       /* Compute new time step */
431d294eb03SHong Zhang       *dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],*dt);
432d294eb03SHong Zhang       fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
433ea3dac1cSHong Zhang       fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
434d294eb03SHong Zhang       switch (event->direction[i]) {
435d294eb03SHong Zhang       case -1:
436d294eb03SHong Zhang         if (fvalue_sign < 0) {
437ea3dac1cSHong Zhang           rollback = 1;
438d294eb03SHong Zhang           event->fvalue_right[i] = event->fvalue[i];
439d294eb03SHong Zhang           event->side[i] = 1;
440d294eb03SHong Zhang         }
441d294eb03SHong Zhang         break;
442d294eb03SHong Zhang       case 1:
443d294eb03SHong Zhang         if (fvalue_sign > 0) {
444ea3dac1cSHong Zhang           rollback = 1;
445d294eb03SHong Zhang           event->fvalue_right[i] = event->fvalue[i];
446d294eb03SHong Zhang           event->side[i] = 1;
447d294eb03SHong Zhang         }
448d294eb03SHong Zhang         break;
449d294eb03SHong Zhang       case 0:
450ea3dac1cSHong Zhang         if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */
451ea3dac1cSHong Zhang           rollback = 1;
452d294eb03SHong Zhang           event->fvalue_right[i] = event->fvalue[i];
453d294eb03SHong Zhang           event->side[i] = 1;
454ea3dac1cSHong Zhang         }
455d294eb03SHong Zhang         break;
456d294eb03SHong Zhang       }
457ea3dac1cSHong Zhang       if (event->status == TSEVENT_PROCESSING) event->side[i] = -1;
458d294eb03SHong Zhang     }
459d294eb03SHong Zhang   }
460d5f990dbSBarry Smith   in[0] = (PetscInt)event->status; in[1] = rollback;
461*1c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts)));
462ea3dac1cSHong Zhang   event->status = (TSEventStatus)out[0]; rollback = out[1];
463ea3dac1cSHong Zhang   /* 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 */
464ea3dac1cSHong Zhang   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
465d294eb03SHong Zhang   if (event->status == TSEVENT_ZERO) {
466ea3dac1cSHong Zhang     for (i=0; i < event->nevents; i++) {
467ea3dac1cSHong Zhang       if (event->zerocrossing[i]) {
468ea3dac1cSHong Zhang         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 */
469ea3dac1cSHong Zhang           event->events_zero[event->nevents_zero++] = i;
470ea3dac1cSHong Zhang           if (event->monitor) {
4719566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D zero crossing located at time %g\n",event->iterctr,i,(double)t));
472ea3dac1cSHong Zhang           }
473ea3dac1cSHong Zhang           event->zerocrossing[i] = PETSC_FALSE;
474ea3dac1cSHong Zhang         }
475ea3dac1cSHong Zhang       }
476ea3dac1cSHong Zhang       event->side[i] = 0;
477ea3dac1cSHong Zhang     }
478d294eb03SHong Zhang   }
479d294eb03SHong Zhang   PetscFunctionReturn(0);
480d294eb03SHong Zhang }
48138bf2713SShri Abhyankar 
4826427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts)
4832dc7a7e3SShri Abhyankar {
4846427ac75SLisandro Dalcin   TSEvent        event;
4852dc7a7e3SShri Abhyankar   PetscReal      t;
4862dc7a7e3SShri Abhyankar   Vec            U;
4872dc7a7e3SShri Abhyankar   PetscInt       i;
488ea3dac1cSHong Zhang   PetscReal      dt,dt_min,dt_reset = 0.0;
4892dc7a7e3SShri Abhyankar 
4902dc7a7e3SShri Abhyankar   PetscFunctionBegin;
4916427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4926427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
4936427ac75SLisandro Dalcin   event = ts->event;
4942dc7a7e3SShri Abhyankar 
4959566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&t));
4969566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts,&dt));
4979566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&U));
4982dc7a7e3SShri Abhyankar 
4997dbe0728SLisandro Dalcin   if (event->status == TSEVENT_NONE) {
5007dbe0728SLisandro Dalcin     event->timestep_prev = dt;
501d294eb03SHong Zhang     event->ptime_end = t;
5027dbe0728SLisandro Dalcin   }
5032dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_RESET_NEXTSTEP) {
504d294eb03SHong Zhang     /* user has specified a PostEventInterval dt */
505458122a4SShri Abhyankar     dt = event->timestep_posteventinterval;
506e97c63d7SStefano Zampini     if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
507e97c63d7SStefano Zampini       PetscReal maxdt = ts->max_time-t;
508e97c63d7SStefano Zampini       dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt);
509e97c63d7SStefano Zampini     }
5109566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts,dt));
5112dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
5122dc7a7e3SShri Abhyankar   }
5132dc7a7e3SShri Abhyankar 
5149566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(U));
5159566063dSJacob Faibussowitsch   PetscCall((*event->eventhandler)(ts,t,U,event->fvalue,event->ctx));
5169566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(U));
5179e12be75SShri Abhyankar 
518d294eb03SHong Zhang   /* Detect the events */
5199566063dSJacob Faibussowitsch   PetscCall(TSEventDetection(ts));
520d294eb03SHong Zhang 
521d294eb03SHong Zhang   /* Locate the events */
522d294eb03SHong Zhang   if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) {
523d294eb03SHong Zhang     /* Approach the zero crosing by setting a new step size */
5249566063dSJacob Faibussowitsch     PetscCall(TSEventLocation(ts,&dt));
525d294eb03SHong Zhang     /* Roll back when new events are detected */
526d294eb03SHong Zhang     if (event->status == TSEVENT_LOCATED_INTERVAL) {
5279566063dSJacob Faibussowitsch       PetscCall(TSRollBack(ts));
5289566063dSJacob Faibussowitsch       PetscCall(TSSetConvergedReason(ts,TS_CONVERGED_ITERATING));
529d294eb03SHong Zhang       event->iterctr++;
530006e6a18SShri Abhyankar     }
531*1c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts)));
532ea3dac1cSHong Zhang     if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset;
5339566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts,dt_min));
534d294eb03SHong Zhang     /* Found the zero crossing */
5359e12be75SShri Abhyankar     if (event->status == TSEVENT_ZERO) {
5369566063dSJacob Faibussowitsch       PetscCall(TSPostEvent(ts,t,U));
5379e12be75SShri Abhyankar 
538e2cdd850SShri Abhyankar       dt = event->ptime_end - t;
539ad508104SStefano Zampini       if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */
540ad508104SStefano Zampini         dt = event->timestep_prev;
541ad508104SStefano Zampini         event->status = TSEVENT_NONE;
542ad508104SStefano Zampini       }
543d294eb03SHong Zhang       if (event->timestep_postevent) { /* user has specified a PostEvent dt*/
544d294eb03SHong Zhang         dt = event->timestep_postevent;
545d294eb03SHong Zhang       }
546e97c63d7SStefano Zampini       if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
547e97c63d7SStefano Zampini         PetscReal maxdt = ts->max_time-t;
548e97c63d7SStefano Zampini         dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt);
549e97c63d7SStefano Zampini       }
5509566063dSJacob Faibussowitsch       PetscCall(TSSetTimeStep(ts,dt));
55138bf2713SShri Abhyankar       event->iterctr = 0;
5529e12be75SShri Abhyankar     }
553d294eb03SHong Zhang     /* Have not found the zero crosing yet */
554d294eb03SHong Zhang     if (event->status == TSEVENT_PROCESSING) {
555d294eb03SHong Zhang       if (event->monitor) {
5569566063dSJacob Faibussowitsch         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));
557d294eb03SHong Zhang       }
558d294eb03SHong Zhang       event->iterctr++;
559d294eb03SHong Zhang     }
560d294eb03SHong Zhang   }
561d294eb03SHong Zhang   if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */
5622dc7a7e3SShri Abhyankar     event->status = TSEVENT_PROCESSING;
563e2cdd850SShri Abhyankar     event->ptime_right = t;
5642dc7a7e3SShri Abhyankar   } else {
565d294eb03SHong Zhang     for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
56638bf2713SShri Abhyankar     event->ptime_prev = t;
56738bf2713SShri Abhyankar   }
5682dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
5692dc7a7e3SShri Abhyankar }
5702dc7a7e3SShri Abhyankar 
5716427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts)
572d0578d90SShri Abhyankar {
5736427ac75SLisandro Dalcin   TSEvent        event;
574d0578d90SShri Abhyankar   PetscReal      t;
575d0578d90SShri Abhyankar   Vec            U;
576d0578d90SShri Abhyankar   PetscInt       ctr;
577d0578d90SShri Abhyankar   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
578d0578d90SShri Abhyankar 
579d0578d90SShri Abhyankar   PetscFunctionBegin;
5806427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5816427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
5826427ac75SLisandro Dalcin   event = ts->event;
583d0578d90SShri Abhyankar 
5849566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&t));
5859566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&U));
586d0578d90SShri Abhyankar 
587d0578d90SShri Abhyankar   ctr = event->recorder.ctr-1;
588bcbf8bb3SShri Abhyankar   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
589d0578d90SShri Abhyankar     /* Call the user postevent function */
590d0578d90SShri Abhyankar     if (event->postevent) {
5919566063dSJacob Faibussowitsch       PetscCall((*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx));
592d0578d90SShri Abhyankar       event->recorder.ctr--;
593d0578d90SShri Abhyankar     }
594d0578d90SShri Abhyankar   }
595d0578d90SShri Abhyankar 
596d0578d90SShri Abhyankar   PetscBarrier((PetscObject)ts);
597d0578d90SShri Abhyankar   PetscFunctionReturn(0);
598d0578d90SShri Abhyankar }
5991ea83e56SMiguel 
6001ea83e56SMiguel /*@
6011ea83e56SMiguel   TSGetNumEvents - Get the numbers of events set
6021ea83e56SMiguel 
6031ea83e56SMiguel   Logically Collective
6041ea83e56SMiguel 
6054165533cSJose E. Roman   Input Parameter:
60699c90e12SSatish Balay . ts - the TS context
6071ea83e56SMiguel 
6084165533cSJose E. Roman   Output Parameter:
6091ea83e56SMiguel . nevents - number of events
6101ea83e56SMiguel 
6111ea83e56SMiguel   Level: intermediate
6121ea83e56SMiguel 
6131ea83e56SMiguel .seealso: TSSetEventHandler()
6141ea83e56SMiguel 
6151ea83e56SMiguel @*/
6161ea83e56SMiguel PetscErrorCode TSGetNumEvents(TS ts,PetscInt * nevents)
6171ea83e56SMiguel {
6181ea83e56SMiguel   PetscFunctionBegin;
6191ea83e56SMiguel   *nevents = ts->event->nevents;
6201ea83e56SMiguel   PetscFunctionReturn(0);
6211ea83e56SMiguel }
622