xref: /petsc/src/ts/event/tsevent.c (revision 20f4b53cbb5e9bd9ef12b76a8697d60d197cda17)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
22dc7a7e3SShri Abhyankar 
36fea3669SShri Abhyankar /*
46427ac75SLisandro Dalcin   TSEventInitialize - Initializes TSEvent for TSSolve
56fea3669SShri Abhyankar */
6d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventInitialize(TSEvent event, TS ts, PetscReal t, Vec U)
7d71ae5a4SJacob Faibussowitsch {
86fea3669SShri Abhyankar   PetscFunctionBegin;
93ba16761SJacob Faibussowitsch   if (!event) PetscFunctionReturn(PETSC_SUCCESS);
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));
163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
176fea3669SShri Abhyankar }
186fea3669SShri Abhyankar 
19d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventDestroy(TSEvent *event)
20d71ae5a4SJacob Faibussowitsch {
217dbe0728SLisandro Dalcin   PetscInt i;
227dbe0728SLisandro Dalcin 
237dbe0728SLisandro Dalcin   PetscFunctionBegin;
247dbe0728SLisandro Dalcin   PetscValidPointer(event, 1);
253ba16761SJacob Faibussowitsch   if (!*event) PetscFunctionReturn(PETSC_SUCCESS);
269371c9d4SSatish Balay   if (--(*event)->refct > 0) {
279371c9d4SSatish Balay     *event = NULL;
283ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
299371c9d4SSatish Balay   }
30e7069c78SShri 
319566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue));
329566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue_prev));
339566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->fvalue_right));
349566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->zerocrossing));
359566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->side));
369566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->direction));
379566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->terminate));
389566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->events_zero));
399566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->vtol));
40a4ffd976SShri Abhyankar 
4148a46eb9SPierre Jolivet   for (i = 0; i < (*event)->recsize; i++) PetscCall(PetscFree((*event)->recorder.eventidx[i]));
429566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.eventidx));
439566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.nevents));
449566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.stepnum));
459566063dSJacob Faibussowitsch   PetscCall(PetscFree((*event)->recorder.time));
46a4ffd976SShri Abhyankar 
479566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*event)->monitor));
489566063dSJacob Faibussowitsch   PetscCall(PetscFree(*event));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
507dbe0728SLisandro Dalcin }
517dbe0728SLisandro Dalcin 
52e3005195SShri Abhyankar /*@
53*20f4b53cSBarry Smith   TSSetPostEventIntervalStep - Set the time-step used immediately following an event interval
54458122a4SShri Abhyankar 
55458122a4SShri Abhyankar   Logically Collective
56458122a4SShri Abhyankar 
574165533cSJose E. Roman   Input Parameters:
58458122a4SShri Abhyankar + ts - time integration context
59458122a4SShri Abhyankar - dt - post event interval step
60458122a4SShri Abhyankar 
61bcf0153eSBarry Smith   Options Database Key:
62458122a4SShri Abhyankar . -ts_event_post_eventinterval_step <dt> time-step after event interval
63458122a4SShri Abhyankar 
64bcf0153eSBarry Smith   Level: advanced
65458122a4SShri Abhyankar 
66bcf0153eSBarry Smith   Notes:
67bcf0153eSBarry Smith  `TSSetPostEventIntervalStep()` allows one to set a time-step that is used immediately following an event interval.
68bcf0153eSBarry Smith 
69bcf0153eSBarry Smith   This function should be called from the postevent function set with `TSSetEventHandler()`.
70458122a4SShri Abhyankar 
71458122a4SShri Abhyankar   The post event interval time-step should be selected based on the dynamics following the event.
72458122a4SShri Abhyankar   If the dynamics are stiff, a conservative (small) step should be used.
73458122a4SShri Abhyankar   If not, then a larger time-step can be used.
74458122a4SShri Abhyankar 
75bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSEvent`, `TSSetEventHandler()`
76458122a4SShri Abhyankar @*/
77d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt)
78d71ae5a4SJacob Faibussowitsch {
79458122a4SShri Abhyankar   PetscFunctionBegin;
80458122a4SShri Abhyankar   ts->event->timestep_posteventinterval = dt;
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82458122a4SShri Abhyankar }
83458122a4SShri Abhyankar 
84458122a4SShri Abhyankar /*@
85*20f4b53cSBarry Smith    TSSetEventTolerances - Set tolerances for event zero crossings
86e3005195SShri Abhyankar 
87e3005195SShri Abhyankar    Logically Collective
88e3005195SShri Abhyankar 
894165533cSJose E. Roman    Input Parameters:
90e3005195SShri Abhyankar +  ts - time integration context
91bcf0153eSBarry Smith .  tol - scalar tolerance, `PETSC_DECIDE` to leave current value
92*20f4b53cSBarry Smith -  vtol - array of tolerances or `NULL`, used in preference to tol if present
93e3005195SShri Abhyankar 
94bcf0153eSBarry Smith    Options Database Key:
95147403d9SBarry Smith .  -ts_event_tol <tol> - tolerance for event zero crossing
96e3005195SShri Abhyankar 
97e3005195SShri Abhyankar    Level: beginner
98e3005195SShri Abhyankar 
99bcf0153eSBarry Smith    Notes:
100bcf0153eSBarry Smith    Must call `TSSetEventHandler(`) before setting the tolerances.
101bcf0153eSBarry Smith 
102*20f4b53cSBarry Smith    The size of `vtol` is equal to the number of events.
103*20f4b53cSBarry Smith 
104*20f4b53cSBarry Smith    The tolerance is some measure of how close the event function is to zero for the event detector to stop
105*20f4b53cSBarry Smith    and declare the time of the event has been detected.
106bcf0153eSBarry Smith 
107bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSEvent`, `TSSetEventHandler()`
108e3005195SShri Abhyankar @*/
109d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetEventTolerances(TS ts, PetscReal tol, PetscReal vtol[])
110d71ae5a4SJacob Faibussowitsch {
1116427ac75SLisandro Dalcin   TSEvent  event;
112e3005195SShri Abhyankar   PetscInt i;
113e3005195SShri Abhyankar 
114e3005195SShri Abhyankar   PetscFunctionBegin;
1156427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1166427ac75SLisandro Dalcin   if (vtol) PetscValidRealPointer(vtol, 3);
1173c633725SBarry Smith   PetscCheck(ts->event, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set the events first by calling TSSetEventHandler()");
1186427ac75SLisandro Dalcin 
1196427ac75SLisandro Dalcin   event = ts->event;
120e3005195SShri Abhyankar   if (vtol) {
121e3005195SShri Abhyankar     for (i = 0; i < event->nevents; i++) event->vtol[i] = vtol[i];
122e3005195SShri Abhyankar   } else {
123f3fa974cSJacob Faibussowitsch     if ((tol != (PetscReal)PETSC_DECIDE) || (tol != (PetscReal)PETSC_DEFAULT)) {
124e3005195SShri Abhyankar       for (i = 0; i < event->nevents; i++) event->vtol[i] = tol;
125e3005195SShri Abhyankar     }
126e3005195SShri Abhyankar   }
1273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128e3005195SShri Abhyankar }
129e3005195SShri Abhyankar 
1302dc7a7e3SShri Abhyankar /*@C
131ac6a796dSBarry Smith    TSSetEventHandler - Sets a function used for detecting events
1322dc7a7e3SShri Abhyankar 
133c3339decSBarry Smith    Logically Collective
1342dc7a7e3SShri Abhyankar 
1352dc7a7e3SShri Abhyankar    Input Parameters:
136bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
1372dc7a7e3SShri Abhyankar .  nevents - number of local events
138d94325d3SShri Abhyankar .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
139d94325d3SShri Abhyankar                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
140d94325d3SShri Abhyankar .  terminate - flag to indicate whether time stepping should be terminated after
141d94325d3SShri Abhyankar                event is detected (one for each event)
142*20f4b53cSBarry Smith .  eventhandler - a change in sign of this function (see `direction`) is used to determine an even has occurred
143*20f4b53cSBarry Smith .  postevent - [optional] post-event function, this function can change properties of the solution, ODE etc at the time of the event
1442589fa24SLisandro Dalcin -  ctx       - [optional] user-defined context for private data for the
145*20f4b53cSBarry Smith                event detector and post event routine (use `NULL` if no
1462dc7a7e3SShri Abhyankar                context is desired)
1472dc7a7e3SShri Abhyankar 
148*20f4b53cSBarry Smith    Calling sequence of `eventhandler`:
149*20f4b53cSBarry Smith $   PetscErrorCode eventhandler(TS ts, PetscReal t, Vec U, PetscScalar fvalue[], void* ctx)
1502dc7a7e3SShri Abhyankar +  ts  - the TS context
1512dc7a7e3SShri Abhyankar .  t   - current time
1522dc7a7e3SShri Abhyankar .  U   - current iterate
153*20f4b53cSBarry Smith .  ctx - [optional] context passed with eventhandler
154*20f4b53cSBarry Smith -  fvalue    - function value of events at time t
1552dc7a7e3SShri Abhyankar 
156*20f4b53cSBarry Smith    Calling sequence of `postevent`:
157*20f4b53cSBarry Smith $   PetscErrorCode postevent(TS ts, PetscInt nevents_zero, PetscInt events_zero[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
1582dc7a7e3SShri Abhyankar +  ts - the TS context
1592dc7a7e3SShri Abhyankar .  nevents_zero - number of local events whose event function is zero
1602dc7a7e3SShri Abhyankar .  events_zero  - indices of local events which have reached zero
1612dc7a7e3SShri Abhyankar .  t            - current time
1622dc7a7e3SShri Abhyankar .  U            - current solution
163031fbad4SShri Abhyankar .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
1646427ac75SLisandro Dalcin -  ctx          - the context passed with eventhandler
1652dc7a7e3SShri Abhyankar 
1662dc7a7e3SShri Abhyankar    Level: intermediate
1672dc7a7e3SShri Abhyankar 
168*20f4b53cSBarry Smith    Note:
169*20f4b53cSBarry Smith    The `eventhandler` is actually the event detector function and the `postevent` function actually handles the desired changes that
170*20f4b53cSBarry Smith    should take place at the time of the event
171*20f4b53cSBarry Smith 
172bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSEvent`, `TSCreate()`, `TSSetTimeStep()`, `TSSetConvergedReason()`
1732dc7a7e3SShri Abhyankar @*/
174d71ae5a4SJacob Faibussowitsch 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)
175d71ae5a4SJacob Faibussowitsch {
176d294eb03SHong Zhang   TSAdapt   adapt;
177d294eb03SHong Zhang   PetscReal hmin;
1782dc7a7e3SShri Abhyankar   TSEvent   event;
179d94325d3SShri Abhyankar   PetscInt  i;
180006e6a18SShri Abhyankar   PetscBool flg;
181a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE
182a6c783d2SShri Abhyankar   PetscReal tol = 1e-4;
183a6c783d2SShri Abhyankar #else
184d569cc17SSatish Balay   PetscReal tol = 1e-6;
185a6c783d2SShri Abhyankar #endif
1862dc7a7e3SShri Abhyankar 
1872dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1886427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1890a82b154SShri   if (nevents) {
190064a246eSJacob Faibussowitsch     PetscValidIntPointer(direction, 3);
191064a246eSJacob Faibussowitsch     PetscValidBoolPointer(terminate, 4);
1920a82b154SShri   }
1936427ac75SLisandro Dalcin 
1944dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&event));
1959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->fvalue));
1969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->fvalue_prev));
1979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->fvalue_right));
1989566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->zerocrossing));
1999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->side));
2009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->direction));
2019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->terminate));
2029566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->vtol));
203d94325d3SShri Abhyankar   for (i = 0; i < nevents; i++) {
204d94325d3SShri Abhyankar     event->direction[i]    = direction[i];
205d94325d3SShri Abhyankar     event->terminate[i]    = terminate[i];
206e2cdd850SShri Abhyankar     event->zerocrossing[i] = PETSC_FALSE;
207e2cdd850SShri Abhyankar     event->side[i]         = 0;
208d94325d3SShri Abhyankar   }
2099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nevents, &event->events_zero));
2102589fa24SLisandro Dalcin   event->nevents                    = nevents;
2116427ac75SLisandro Dalcin   event->eventhandler               = eventhandler;
2122dc7a7e3SShri Abhyankar   event->postevent                  = postevent;
2136427ac75SLisandro Dalcin   event->ctx                        = ctx;
214458122a4SShri Abhyankar   event->timestep_posteventinterval = ts->time_step;
2159566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &adapt));
2169566063dSJacob Faibussowitsch   PetscCall(TSAdaptGetStepLimits(adapt, &hmin, NULL));
217d294eb03SHong Zhang   event->timestep_min = hmin;
2182dc7a7e3SShri Abhyankar 
21902749585SLisandro Dalcin   event->recsize = 8; /* Initial size of the recorder */
220d0609cedSBarry Smith   PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "TS Event options", "TS");
2212dc7a7e3SShri Abhyankar   {
2229566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_tol", "Scalar event tolerance for zero crossing check", "TSSetEventTolerances", tol, &tol, NULL));
2239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsName("-ts_event_monitor", "Print choices made by event handler", "", &flg));
2249566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size", "Initial size of event recorder", "", event->recsize, &event->recsize, NULL));
2259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_post_eventinterval_step", "Time step after event interval", "", event->timestep_posteventinterval, &event->timestep_posteventinterval, NULL));
2269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_post_event_step", "Time step after event", "", event->timestep_postevent, &event->timestep_postevent, NULL));
2279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_event_dt_min", "Minimum time step considered for TSEvent", "", event->timestep_min, &event->timestep_min, NULL));
2282dc7a7e3SShri Abhyankar   }
229d0609cedSBarry Smith   PetscOptionsEnd();
230a4ffd976SShri Abhyankar 
2319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.time));
2329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.stepnum));
2339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.nevents));
2349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(event->recsize, &event->recorder.eventidx));
23548a46eb9SPierre Jolivet   for (i = 0; i < event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &event->recorder.eventidx[i]));
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;
2453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 */
252d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventRecorderResize(TSEvent event)
253d71ae5a4SJacob Faibussowitsch {
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));
26748a46eb9SPierre Jolivet   for (i = 0; i < fact * event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &eventidx[i]));
268a4ffd976SShri Abhyankar 
269a4ffd976SShri Abhyankar   /* Copy over data */
2709566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(time, event->recorder.time, event->recsize));
2719566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(stepnum, event->recorder.stepnum, event->recsize));
2729566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(nevents, event->recorder.nevents, event->recsize));
27348a46eb9SPierre Jolivet   for (i = 0; i < event->recsize; i++) PetscCall(PetscArraycpy(eventidx[i], event->recorder.eventidx[i], event->recorder.nevents[i]));
274a4ffd976SShri Abhyankar 
275a4ffd976SShri Abhyankar   /* Destroy old arrays */
27648a46eb9SPierre Jolivet   for (i = 0; i < event->recsize; i++) PetscCall(PetscFree(event->recorder.eventidx[i]));
2779566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.eventidx));
2789566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.nevents));
2799566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.stepnum));
2809566063dSJacob Faibussowitsch   PetscCall(PetscFree(event->recorder.time));
281a4ffd976SShri Abhyankar 
282a4ffd976SShri Abhyankar   /* Set pointers */
283a4ffd976SShri Abhyankar   event->recorder.time     = time;
284a4ffd976SShri Abhyankar   event->recorder.stepnum  = stepnum;
285a4ffd976SShri Abhyankar   event->recorder.nevents  = nevents;
286a4ffd976SShri Abhyankar   event->recorder.eventidx = eventidx;
287a4ffd976SShri Abhyankar 
288a4ffd976SShri Abhyankar   /* Double size */
289a4ffd976SShri Abhyankar   event->recsize *= fact;
290a4ffd976SShri Abhyankar 
2913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
292a4ffd976SShri Abhyankar }
293a4ffd976SShri Abhyankar 
294031fbad4SShri Abhyankar /*
295ac6a796dSBarry Smith    Helper routine to handle user postevents and recording
296031fbad4SShri Abhyankar */
297d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPostEvent(TS ts, PetscReal t, Vec U)
298d71ae5a4SJacob Faibussowitsch {
2992dc7a7e3SShri Abhyankar   TSEvent   event     = ts->event;
3002dc7a7e3SShri Abhyankar   PetscBool terminate = PETSC_FALSE;
30128d5b5d6SLisandro Dalcin   PetscBool restart   = PETSC_FALSE;
302d0578d90SShri Abhyankar   PetscInt  i, ctr, stepnum;
3037324a0ffSLisandro Dalcin   PetscBool inflag[2], outflag[2];
3044597913aSLisandro Dalcin   PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
3052dc7a7e3SShri Abhyankar 
3062dc7a7e3SShri Abhyankar   PetscFunctionBegin;
3072dc7a7e3SShri Abhyankar   if (event->postevent) {
30828d5b5d6SLisandro Dalcin     PetscObjectState state_prev, state_post;
3099566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateGet((PetscObject)U, &state_prev));
3109566063dSJacob Faibussowitsch     PetscCall((*event->postevent)(ts, event->nevents_zero, event->events_zero, t, U, forwardsolve, event->ctx));
3119566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateGet((PetscObject)U, &state_post));
31228d5b5d6SLisandro Dalcin     if (state_prev != state_post) restart = PETSC_TRUE;
3132dc7a7e3SShri Abhyankar   }
3144597913aSLisandro Dalcin 
31528d5b5d6SLisandro Dalcin   /* Handle termination events and step restart */
3169371c9d4SSatish Balay   for (i = 0; i < event->nevents_zero; i++)
3179371c9d4SSatish Balay     if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
3189371c9d4SSatish Balay   inflag[0] = restart;
3199371c9d4SSatish Balay   inflag[1] = terminate;
3201c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(inflag, outflag, 2, MPIU_BOOL, MPI_LOR, ((PetscObject)ts)->comm));
3219371c9d4SSatish Balay   restart   = outflag[0];
3229371c9d4SSatish Balay   terminate = outflag[1];
3239566063dSJacob Faibussowitsch   if (restart) PetscCall(TSRestartStep(ts));
3249566063dSJacob Faibussowitsch   if (terminate) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_EVENT));
3257324a0ffSLisandro Dalcin   event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
326f7aea88cSShri Abhyankar 
3274597913aSLisandro Dalcin   /* Reset event residual functions as states might get changed by the postevent callback */
328f443add6SLisandro Dalcin   if (event->postevent) {
3299566063dSJacob Faibussowitsch     PetscCall(VecLockReadPush(U));
3309566063dSJacob Faibussowitsch     PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx));
3319566063dSJacob Faibussowitsch     PetscCall(VecLockReadPop(U));
332f443add6SLisandro Dalcin   }
333f443add6SLisandro Dalcin 
334f443add6SLisandro Dalcin   /* Cache current time and event residual functions */
335f443add6SLisandro Dalcin   event->ptime_prev = t;
3369371c9d4SSatish Balay   for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
3374597913aSLisandro Dalcin 
338d0578d90SShri Abhyankar   /* Record the event in the event recorder */
3399566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &stepnum));
340f7aea88cSShri Abhyankar   ctr = event->recorder.ctr;
34148a46eb9SPierre Jolivet   if (ctr == event->recsize) PetscCall(TSEventRecorderResize(event));
342f7aea88cSShri Abhyankar   event->recorder.time[ctr]    = t;
343d0578d90SShri Abhyankar   event->recorder.stepnum[ctr] = stepnum;
3444597913aSLisandro Dalcin   event->recorder.nevents[ctr] = event->nevents_zero;
3454597913aSLisandro Dalcin   for (i = 0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
346f7aea88cSShri Abhyankar   event->recorder.ctr++;
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3482dc7a7e3SShri Abhyankar }
3492dc7a7e3SShri Abhyankar 
35002749585SLisandro Dalcin /* Uses Anderson-Bjorck variant of regula falsi method */
351d71ae5a4SJacob Faibussowitsch static inline PetscReal TSEventComputeStepSize(PetscReal tleft, PetscReal t, PetscReal tright, PetscScalar fleft, PetscScalar f, PetscScalar fright, PetscInt side, PetscReal dt)
352d71ae5a4SJacob Faibussowitsch {
35302749585SLisandro Dalcin   PetscReal new_dt, scal = 1.0;
354e2cdd850SShri Abhyankar   if (PetscRealPart(fleft) * PetscRealPart(f) < 0) {
355e2cdd850SShri Abhyankar     if (side == 1) {
356a6c783d2SShri Abhyankar       scal = (PetscRealPart(fright) - PetscRealPart(f)) / PetscRealPart(fright);
357a6c783d2SShri Abhyankar       if (scal < PETSC_SMALL) scal = 0.5;
358e2cdd850SShri Abhyankar     }
35902749585SLisandro Dalcin     new_dt = (scal * PetscRealPart(fleft) * t - PetscRealPart(f) * tleft) / (scal * PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
360e2cdd850SShri Abhyankar   } else {
361e2cdd850SShri Abhyankar     if (side == -1) {
362a6c783d2SShri Abhyankar       scal = (PetscRealPart(fleft) - PetscRealPart(f)) / PetscRealPart(fleft);
363a6c783d2SShri Abhyankar       if (scal < PETSC_SMALL) scal = 0.5;
364e2cdd850SShri Abhyankar     }
36502749585SLisandro Dalcin     new_dt = (PetscRealPart(f) * tright - scal * PetscRealPart(fright) * t) / (PetscRealPart(f) - scal * PetscRealPart(fright)) - t;
366e2cdd850SShri Abhyankar   }
36702749585SLisandro Dalcin   return PetscMin(dt, new_dt);
36838bf2713SShri Abhyankar }
369e2cdd850SShri Abhyankar 
370d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventDetection(TS ts)
371d71ae5a4SJacob Faibussowitsch {
372d294eb03SHong Zhang   TSEvent   event = ts->event;
373d294eb03SHong Zhang   PetscReal t;
374d294eb03SHong Zhang   PetscInt  i;
375d294eb03SHong Zhang   PetscInt  fvalue_sign, fvalueprev_sign;
376ea3dac1cSHong Zhang   PetscInt  in, out;
377d294eb03SHong Zhang 
378d294eb03SHong Zhang   PetscFunctionBegin;
3799566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
380d294eb03SHong Zhang   for (i = 0; i < event->nevents; i++) {
381d294eb03SHong Zhang     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
382d294eb03SHong Zhang       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
383d294eb03SHong Zhang       event->status = TSEVENT_LOCATED_INTERVAL;
384d294eb03SHong Zhang       if (event->monitor) {
38563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " interval detected due to zero value (tol=%g) [%g - %g]\n", event->iterctr, i, (double)event->vtol[i], (double)event->ptime_prev, (double)t));
386d294eb03SHong Zhang       }
387d294eb03SHong Zhang       continue;
388d294eb03SHong Zhang     }
389ea3dac1cSHong Zhang     if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */
390d294eb03SHong Zhang     fvalue_sign     = PetscSign(PetscRealPart(event->fvalue[i]));
391d294eb03SHong Zhang     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
392d294eb03SHong Zhang     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) {
393d294eb03SHong Zhang       if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
394d294eb03SHong Zhang       event->status = TSEVENT_LOCATED_INTERVAL;
39548a46eb9SPierre Jolivet       if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " interval detected due to sign change [%g - %g]\n", event->iterctr, i, (double)event->ptime_prev, (double)t));
396d294eb03SHong Zhang     }
397d294eb03SHong Zhang   }
398d5f990dbSBarry Smith   in = (PetscInt)event->status;
3991c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&in, &out, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts)));
400ea3dac1cSHong Zhang   event->status = (TSEventStatus)out;
4013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
402d294eb03SHong Zhang }
403d294eb03SHong Zhang 
404d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventLocation(TS ts, PetscReal *dt)
405d71ae5a4SJacob Faibussowitsch {
406d294eb03SHong Zhang   TSEvent   event = ts->event;
407d294eb03SHong Zhang   PetscInt  i;
408d294eb03SHong Zhang   PetscReal t;
409ea3dac1cSHong Zhang   PetscInt  fvalue_sign, fvalueprev_sign;
410ea3dac1cSHong Zhang   PetscInt  rollback = 0, in[2], out[2];
411d294eb03SHong Zhang 
412d294eb03SHong Zhang   PetscFunctionBegin;
4139566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
414d294eb03SHong Zhang   event->nevents_zero = 0;
415d294eb03SHong Zhang   for (i = 0; i < event->nevents; i++) {
416d294eb03SHong Zhang     if (event->zerocrossing[i]) {
417d294eb03SHong 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 */
418d294eb03SHong Zhang         event->status          = TSEVENT_ZERO;
419d294eb03SHong Zhang         event->fvalue_right[i] = event->fvalue[i];
420d294eb03SHong Zhang         continue;
421d294eb03SHong Zhang       }
422d294eb03SHong Zhang       /* Compute new time step */
423d294eb03SHong 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);
424d294eb03SHong Zhang       fvalue_sign     = PetscSign(PetscRealPart(event->fvalue[i]));
425ea3dac1cSHong Zhang       fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
426d294eb03SHong Zhang       switch (event->direction[i]) {
427d294eb03SHong Zhang       case -1:
428d294eb03SHong Zhang         if (fvalue_sign < 0) {
429ea3dac1cSHong Zhang           rollback               = 1;
430d294eb03SHong Zhang           event->fvalue_right[i] = event->fvalue[i];
431d294eb03SHong Zhang           event->side[i]         = 1;
432d294eb03SHong Zhang         }
433d294eb03SHong Zhang         break;
434d294eb03SHong Zhang       case 1:
435d294eb03SHong Zhang         if (fvalue_sign > 0) {
436ea3dac1cSHong Zhang           rollback               = 1;
437d294eb03SHong Zhang           event->fvalue_right[i] = event->fvalue[i];
438d294eb03SHong Zhang           event->side[i]         = 1;
439d294eb03SHong Zhang         }
440d294eb03SHong Zhang         break;
441d294eb03SHong Zhang       case 0:
442ea3dac1cSHong Zhang         if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */
443ea3dac1cSHong Zhang           rollback               = 1;
444d294eb03SHong Zhang           event->fvalue_right[i] = event->fvalue[i];
445d294eb03SHong Zhang           event->side[i]         = 1;
446ea3dac1cSHong Zhang         }
447d294eb03SHong Zhang         break;
448d294eb03SHong Zhang       }
449ea3dac1cSHong Zhang       if (event->status == TSEVENT_PROCESSING) event->side[i] = -1;
450d294eb03SHong Zhang     }
451d294eb03SHong Zhang   }
4529371c9d4SSatish Balay   in[0] = (PetscInt)event->status;
4539371c9d4SSatish Balay   in[1] = rollback;
4541c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts)));
4559371c9d4SSatish Balay   event->status = (TSEventStatus)out[0];
4569371c9d4SSatish Balay   rollback      = out[1];
457da81f932SPierre Jolivet   /* 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 correct order */
458ea3dac1cSHong Zhang   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
459d294eb03SHong Zhang   if (event->status == TSEVENT_ZERO) {
460ea3dac1cSHong Zhang     for (i = 0; i < event->nevents; i++) {
461ea3dac1cSHong Zhang       if (event->zerocrossing[i]) {
462ea3dac1cSHong 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 */
463ea3dac1cSHong Zhang           event->events_zero[event->nevents_zero++] = i;
46448a46eb9SPierre Jolivet           if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " zero crossing located at time %g\n", event->iterctr, i, (double)t));
465ea3dac1cSHong Zhang           event->zerocrossing[i] = PETSC_FALSE;
466ea3dac1cSHong Zhang         }
467ea3dac1cSHong Zhang       }
468ea3dac1cSHong Zhang       event->side[i] = 0;
469ea3dac1cSHong Zhang     }
470d294eb03SHong Zhang   }
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
472d294eb03SHong Zhang }
47338bf2713SShri Abhyankar 
474d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventHandler(TS ts)
475d71ae5a4SJacob Faibussowitsch {
4766427ac75SLisandro Dalcin   TSEvent   event;
4772dc7a7e3SShri Abhyankar   PetscReal t;
4782dc7a7e3SShri Abhyankar   Vec       U;
4792dc7a7e3SShri Abhyankar   PetscInt  i;
480ea3dac1cSHong Zhang   PetscReal dt, dt_min, dt_reset = 0.0;
4812dc7a7e3SShri Abhyankar 
4822dc7a7e3SShri Abhyankar   PetscFunctionBegin;
4836427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4843ba16761SJacob Faibussowitsch   if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS);
4856427ac75SLisandro Dalcin   event = ts->event;
4862dc7a7e3SShri Abhyankar 
4879566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
4889566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
4899566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
4902dc7a7e3SShri Abhyankar 
4917dbe0728SLisandro Dalcin   if (event->status == TSEVENT_NONE) {
4927dbe0728SLisandro Dalcin     event->timestep_prev = dt;
493d294eb03SHong Zhang     event->ptime_end     = t;
4947dbe0728SLisandro Dalcin   }
4952dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_RESET_NEXTSTEP) {
496d294eb03SHong Zhang     /* user has specified a PostEventInterval dt */
497458122a4SShri Abhyankar     dt = event->timestep_posteventinterval;
498e97c63d7SStefano Zampini     if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
499e97c63d7SStefano Zampini       PetscReal maxdt = ts->max_time - t;
500e97c63d7SStefano Zampini       dt              = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt);
501e97c63d7SStefano Zampini     }
5029566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
5032dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
5042dc7a7e3SShri Abhyankar   }
5052dc7a7e3SShri Abhyankar 
5069566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(U));
5079566063dSJacob Faibussowitsch   PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx));
5089566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(U));
5099e12be75SShri Abhyankar 
510d294eb03SHong Zhang   /* Detect the events */
5119566063dSJacob Faibussowitsch   PetscCall(TSEventDetection(ts));
512d294eb03SHong Zhang 
513d294eb03SHong Zhang   /* Locate the events */
514d294eb03SHong Zhang   if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) {
515d294eb03SHong Zhang     /* Approach the zero crosing by setting a new step size */
5169566063dSJacob Faibussowitsch     PetscCall(TSEventLocation(ts, &dt));
517d294eb03SHong Zhang     /* Roll back when new events are detected */
518d294eb03SHong Zhang     if (event->status == TSEVENT_LOCATED_INTERVAL) {
5199566063dSJacob Faibussowitsch       PetscCall(TSRollBack(ts));
5209566063dSJacob Faibussowitsch       PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING));
521d294eb03SHong Zhang       event->iterctr++;
522006e6a18SShri Abhyankar     }
5231c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&dt, &dt_min, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts)));
524ea3dac1cSHong Zhang     if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset;
5259566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt_min));
526d294eb03SHong Zhang     /* Found the zero crossing */
5279e12be75SShri Abhyankar     if (event->status == TSEVENT_ZERO) {
5289566063dSJacob Faibussowitsch       PetscCall(TSPostEvent(ts, t, U));
5299e12be75SShri Abhyankar 
530e2cdd850SShri Abhyankar       dt = event->ptime_end - t;
531ad508104SStefano Zampini       if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */
532ad508104SStefano Zampini         dt            = event->timestep_prev;
533ad508104SStefano Zampini         event->status = TSEVENT_NONE;
534ad508104SStefano Zampini       }
535d294eb03SHong Zhang       if (event->timestep_postevent) { /* user has specified a PostEvent dt*/
536d294eb03SHong Zhang         dt = event->timestep_postevent;
537d294eb03SHong Zhang       }
538e97c63d7SStefano Zampini       if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
539e97c63d7SStefano Zampini         PetscReal maxdt = ts->max_time - t;
540e97c63d7SStefano Zampini         dt              = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt);
541e97c63d7SStefano Zampini       }
5429566063dSJacob Faibussowitsch       PetscCall(TSSetTimeStep(ts, dt));
54338bf2713SShri Abhyankar       event->iterctr = 0;
5449e12be75SShri Abhyankar     }
545d294eb03SHong Zhang     /* Have not found the zero crosing yet */
546d294eb03SHong Zhang     if (event->status == TSEVENT_PROCESSING) {
54748a46eb9SPierre Jolivet       if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Stepping forward as no event detected in interval [%g - %g]\n", event->iterctr, (double)event->ptime_prev, (double)t));
548d294eb03SHong Zhang       event->iterctr++;
549d294eb03SHong Zhang     }
550d294eb03SHong Zhang   }
551d294eb03SHong Zhang   if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */
5522dc7a7e3SShri Abhyankar     event->status      = TSEVENT_PROCESSING;
553e2cdd850SShri Abhyankar     event->ptime_right = t;
5542dc7a7e3SShri Abhyankar   } else {
555d294eb03SHong Zhang     for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
55638bf2713SShri Abhyankar     event->ptime_prev = t;
55738bf2713SShri Abhyankar   }
5583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5592dc7a7e3SShri Abhyankar }
5602dc7a7e3SShri Abhyankar 
561d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointEventHandler(TS ts)
562d71ae5a4SJacob Faibussowitsch {
5636427ac75SLisandro Dalcin   TSEvent   event;
564d0578d90SShri Abhyankar   PetscReal t;
565d0578d90SShri Abhyankar   Vec       U;
566d0578d90SShri Abhyankar   PetscInt  ctr;
567d0578d90SShri Abhyankar   PetscBool forwardsolve = PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
568d0578d90SShri Abhyankar 
569d0578d90SShri Abhyankar   PetscFunctionBegin;
5706427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5713ba16761SJacob Faibussowitsch   if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS);
5726427ac75SLisandro Dalcin   event = ts->event;
573d0578d90SShri Abhyankar 
5749566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
5759566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
576d0578d90SShri Abhyankar 
577d0578d90SShri Abhyankar   ctr = event->recorder.ctr - 1;
578bcbf8bb3SShri Abhyankar   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
579d0578d90SShri Abhyankar     /* Call the user postevent function */
580d0578d90SShri Abhyankar     if (event->postevent) {
5819566063dSJacob Faibussowitsch       PetscCall((*event->postevent)(ts, event->recorder.nevents[ctr], event->recorder.eventidx[ctr], t, U, forwardsolve, event->ctx));
582d0578d90SShri Abhyankar       event->recorder.ctr--;
583d0578d90SShri Abhyankar     }
584d0578d90SShri Abhyankar   }
585d0578d90SShri Abhyankar 
5863ba16761SJacob Faibussowitsch   PetscCall(PetscBarrier((PetscObject)ts));
5873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
588d0578d90SShri Abhyankar }
5891ea83e56SMiguel 
5901ea83e56SMiguel /*@
591*20f4b53cSBarry Smith   TSGetNumEvents - Get the numbers of events currently set to be detected
5921ea83e56SMiguel 
5931ea83e56SMiguel   Logically Collective
5941ea83e56SMiguel 
5954165533cSJose E. Roman   Input Parameter:
596bcf0153eSBarry Smith . ts - the `TS` context
5971ea83e56SMiguel 
5984165533cSJose E. Roman   Output Parameter:
599*20f4b53cSBarry Smith . nevents - the number of events
6001ea83e56SMiguel 
6011ea83e56SMiguel   Level: intermediate
6021ea83e56SMiguel 
603bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSEvent`, `TSSetEventHandler()`
6041ea83e56SMiguel @*/
605d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetNumEvents(TS ts, PetscInt *nevents)
606d71ae5a4SJacob Faibussowitsch {
6071ea83e56SMiguel   PetscFunctionBegin;
6081ea83e56SMiguel   *nevents = ts->event->nevents;
6093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6101ea83e56SMiguel }
611