1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 22dc7a7e3SShri Abhyankar 36fea3669SShri Abhyankar /* 4*ca4445c7SIlya Fursov Fills array sign[] with signs of array f[]. If abs(f[i]) < vtol[i], the zero sign is taken. 5*ca4445c7SIlya Fursov All arrays should have length 'nev' 6*ca4445c7SIlya Fursov */ 7*ca4445c7SIlya Fursov static inline void TSEventCalcSigns(PetscInt nev, const PetscReal *f, const PetscReal *vtol, PetscInt *sign) 8*ca4445c7SIlya Fursov { 9*ca4445c7SIlya Fursov for (PetscInt i = 0; i < nev; i++) { 10*ca4445c7SIlya Fursov if (PetscAbsReal(f[i]) < vtol[i]) sign[i] = 0; 11*ca4445c7SIlya Fursov else sign[i] = PetscSign(f[i]); 12*ca4445c7SIlya Fursov } 13*ca4445c7SIlya Fursov } 14*ca4445c7SIlya Fursov 15*ca4445c7SIlya Fursov /* 166427ac75SLisandro Dalcin TSEventInitialize - Initializes TSEvent for TSSolve 176fea3669SShri Abhyankar */ 18d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventInitialize(TSEvent event, TS ts, PetscReal t, Vec U) 19d71ae5a4SJacob Faibussowitsch { 206fea3669SShri Abhyankar PetscFunctionBegin; 213ba16761SJacob Faibussowitsch if (!event) PetscFunctionReturn(PETSC_SUCCESS); 224f572ea9SToby Isaac PetscAssertPointer(event, 1); 236427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 246427ac75SLisandro Dalcin PetscValidHeaderSpecific(U, VEC_CLASSID, 4); 256fea3669SShri Abhyankar event->ptime_prev = t; 2638bf2713SShri Abhyankar event->iterctr = 0; 27*ca4445c7SIlya Fursov event->processing = PETSC_FALSE; 28*ca4445c7SIlya Fursov event->revisit_right = PETSC_FALSE; 29*ca4445c7SIlya Fursov PetscCallBack("TSEvent indicator", (*event->indicator)(ts, t, U, event->fvalue_prev, event->ctx)); 30*ca4445c7SIlya Fursov TSEventCalcSigns(event->nevents, event->fvalue_prev, event->vtol, event->fsign_prev); // by this time event->vtol should have been defined 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 326fea3669SShri Abhyankar } 336fea3669SShri Abhyankar 34d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventDestroy(TSEvent *event) 35d71ae5a4SJacob Faibussowitsch { 367dbe0728SLisandro Dalcin PetscFunctionBegin; 374f572ea9SToby Isaac PetscAssertPointer(event, 1); 383ba16761SJacob Faibussowitsch if (!*event) PetscFunctionReturn(PETSC_SUCCESS); 399371c9d4SSatish Balay if (--(*event)->refct > 0) { 409371c9d4SSatish Balay *event = NULL; 413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 429371c9d4SSatish Balay } 43e7069c78SShri 449566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->fvalue_prev)); 45*ca4445c7SIlya Fursov PetscCall(PetscFree((*event)->fvalue)); 469566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->fvalue_right)); 47*ca4445c7SIlya Fursov PetscCall(PetscFree((*event)->fsign_prev)); 48*ca4445c7SIlya Fursov PetscCall(PetscFree((*event)->fsign)); 49*ca4445c7SIlya Fursov PetscCall(PetscFree((*event)->fsign_right)); 509566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->side)); 51*ca4445c7SIlya Fursov PetscCall(PetscFree((*event)->side_prev)); 52*ca4445c7SIlya Fursov PetscCall(PetscFree((*event)->justrefined_AB)); 53*ca4445c7SIlya Fursov PetscCall(PetscFree((*event)->gamma_AB)); 549566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->direction)); 559566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->terminate)); 569566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->events_zero)); 579566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->vtol)); 58a4ffd976SShri Abhyankar 59*ca4445c7SIlya Fursov for (PetscInt i = 0; i < (*event)->recsize; i++) PetscCall(PetscFree((*event)->recorder.eventidx[i])); 609566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.eventidx)); 619566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.nevents)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.stepnum)); 639566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.time)); 64a4ffd976SShri Abhyankar 659566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*event)->monitor)); 669566063dSJacob Faibussowitsch PetscCall(PetscFree(*event)); 673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 687dbe0728SLisandro Dalcin } 697dbe0728SLisandro Dalcin 70e3005195SShri Abhyankar /*@ 71*ca4445c7SIlya Fursov TSSetPostEventStep - Set the time step to use immediately following the event 72*ca4445c7SIlya Fursov 73*ca4445c7SIlya Fursov Logically Collective 74*ca4445c7SIlya Fursov 75*ca4445c7SIlya Fursov Input Parameters: 76*ca4445c7SIlya Fursov + ts - time integration context 77*ca4445c7SIlya Fursov - dt - post event step 78*ca4445c7SIlya Fursov 79*ca4445c7SIlya Fursov Options Database Key: 80*ca4445c7SIlya Fursov . -ts_event_post_event_step <dt> - time step after the event; zero value - to keep using previous time steps 81*ca4445c7SIlya Fursov 82*ca4445c7SIlya Fursov Level: advanced 83*ca4445c7SIlya Fursov 84*ca4445c7SIlya Fursov Notes: 85*ca4445c7SIlya Fursov `TSSetPostEventStep()` allows one to set a time step that is used immediately following an event. 86*ca4445c7SIlya Fursov If a positive real number is specified, it will be applied as is. 87*ca4445c7SIlya Fursov However, if `TSAdapt` is allowed to interfere, a large 'dt' may get truncated, resulting in a smaller actual post-event step. 88*ca4445c7SIlya Fursov 89*ca4445c7SIlya Fursov The post-event time step should be selected based on the post-event dynamics. 90*ca4445c7SIlya Fursov If the dynamics are stiff, or a significant jump in the equations or the state vector has taken place at the event, 91*ca4445c7SIlya Fursov a conservative (small) step should be employed. If not, then a larger time step may be appropriate. 92*ca4445c7SIlya Fursov 93*ca4445c7SIlya Fursov In the latter case, instead of explicitly setting the post-event time step, 94*ca4445c7SIlya Fursov the user may also choose a special option of keeping the time steps used before the event, which is a sort of 'petsc-decide' 95*ca4445c7SIlya Fursov strategy. For this, use special value 0. It will signal the TS to directly step to the time point it planned to visit 96*ca4445c7SIlya Fursov prior to the event interval detection. E.g. if a step t0 -> t1 was planned originally, and an event 'te' occurred, t0 < te < t1, 97*ca4445c7SIlya Fursov then after the event the TS will step: te -> t1. 98*ca4445c7SIlya Fursov Moreover, in this situation the originally planned subsequent step t1 -> t2 will also be preserved. 99*ca4445c7SIlya Fursov 100*ca4445c7SIlya Fursov This function can be called not only in the initial setup, but also inside the `postevent()` callback set with `TSSetEventHandler()`, 101*ca4445c7SIlya Fursov affecting the post-event step for the current event, and the subsequent ones. 102*ca4445c7SIlya Fursov So, the strategy of the post-event time step definition can be adjusted on the fly. 103*ca4445c7SIlya Fursov Even if several events have been triggered in the given time point, only a single postevent handler is invoked, 104*ca4445c7SIlya Fursov and the user is to determine what post-event time step is more appropriate in this situation. 105*ca4445c7SIlya Fursov 106*ca4445c7SIlya Fursov By default (on `TSSetEventHandler()` call), the post-event time step is set equal to the (initial) `TS` time step. 107*ca4445c7SIlya Fursov 108*ca4445c7SIlya Fursov .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()` 109*ca4445c7SIlya Fursov @*/ 110*ca4445c7SIlya Fursov PetscErrorCode TSSetPostEventStep(TS ts, PetscReal dt) 111*ca4445c7SIlya Fursov { 112*ca4445c7SIlya Fursov PetscFunctionBegin; 113*ca4445c7SIlya Fursov ts->event->timestep_postevent = dt; 114*ca4445c7SIlya Fursov PetscFunctionReturn(PETSC_SUCCESS); 115*ca4445c7SIlya Fursov } 116*ca4445c7SIlya Fursov 117*ca4445c7SIlya Fursov /*@ 11820f4b53cSBarry Smith TSSetPostEventIntervalStep - Set the time-step used immediately following an event interval 119458122a4SShri Abhyankar 120458122a4SShri Abhyankar Logically Collective 121458122a4SShri Abhyankar 1224165533cSJose E. Roman Input Parameters: 123458122a4SShri Abhyankar + ts - time integration context 1249b3d3759SBarry Smith - dt - post-event interval step 125458122a4SShri Abhyankar 126bcf0153eSBarry Smith Options Database Key: 127b43aa488SJacob Faibussowitsch . -ts_event_post_eventinterval_step <dt> - time-step after event interval 128458122a4SShri Abhyankar 129bcf0153eSBarry Smith Level: advanced 130458122a4SShri Abhyankar 131bcf0153eSBarry Smith Notes: 132*ca4445c7SIlya Fursov This function is deprecated, and its invocation will throw a runtime error. Use `TSSetPostEventStep()`. 133458122a4SShri Abhyankar 1341cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()` 135458122a4SShri Abhyankar @*/ 136d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt) 137d71ae5a4SJacob Faibussowitsch { 138458122a4SShri Abhyankar PetscFunctionBegin; 139*ca4445c7SIlya Fursov //ts->event->timestep_posteventinterval = dt; 140*ca4445c7SIlya Fursov /* 141*ca4445c7SIlya Fursov This deprecated function is set to always throw a runtime error. Attempting to reproduce its original behaviour, 142*ca4445c7SIlya Fursov i.e. setting the second (not the first) step after event, would break the logic of the TSEventHandler new code. 143*ca4445c7SIlya Fursov */ 144*ca4445c7SIlya Fursov PetscCheck(PETSC_FALSE, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSSetPostEventIntervalStep() is deprecated --> TSSetPostEventStep() should be used instead"); 1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 146458122a4SShri Abhyankar } 147458122a4SShri Abhyankar 148458122a4SShri Abhyankar /*@ 149*ca4445c7SIlya Fursov TSSetEventTolerances - Set tolerances for event (indicator function) zero crossings 150e3005195SShri Abhyankar 151e3005195SShri Abhyankar Logically Collective 152e3005195SShri Abhyankar 1534165533cSJose E. Roman Input Parameters: 154e3005195SShri Abhyankar + ts - time integration context 155*ca4445c7SIlya Fursov . tol - tolerance, `PETSC_DECIDE` or `PETSC_DEFAULT` to leave the current value 1569b3d3759SBarry Smith - vtol - array of tolerances or `NULL`, used in preference to `tol` if present 157e3005195SShri Abhyankar 158bcf0153eSBarry Smith Options Database Key: 159*ca4445c7SIlya Fursov . -ts_event_tol <tol> - tolerance for event (indicator function) zero crossing 160e3005195SShri Abhyankar 161e3005195SShri Abhyankar Level: beginner 162e3005195SShri Abhyankar 163bcf0153eSBarry Smith Notes: 164*ca4445c7SIlya Fursov One must call `TSSetEventHandler()` before setting the tolerances. 165bcf0153eSBarry Smith 166*ca4445c7SIlya Fursov The size of `vtol` should be equal to the number of events on the given process. 16720f4b53cSBarry Smith 168*ca4445c7SIlya Fursov This function can be also called from the `postevent()` callback set with `TSSetEventHandler()`, 169*ca4445c7SIlya Fursov to adjust the tolerances on the fly. 170bcf0153eSBarry Smith 1711cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()` 172e3005195SShri Abhyankar @*/ 173d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetEventTolerances(TS ts, PetscReal tol, PetscReal vtol[]) 174d71ae5a4SJacob Faibussowitsch { 1756427ac75SLisandro Dalcin TSEvent event; 176e3005195SShri Abhyankar 177e3005195SShri Abhyankar PetscFunctionBegin; 1786427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1794f572ea9SToby Isaac if (vtol) PetscAssertPointer(vtol, 3); 1803c633725SBarry Smith PetscCheck(ts->event, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set the events first by calling TSSetEventHandler()"); 1816427ac75SLisandro Dalcin 1826427ac75SLisandro Dalcin event = ts->event; 183e3005195SShri Abhyankar if (vtol) { 184*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 185e3005195SShri Abhyankar } else { 186*ca4445c7SIlya Fursov if (tol != (PetscReal)PETSC_DECIDE && tol != (PetscReal)PETSC_DEFAULT) { 187*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) event->vtol[i] = tol; 188e3005195SShri Abhyankar } 189e3005195SShri Abhyankar } 1903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 191e3005195SShri Abhyankar } 192e3005195SShri Abhyankar 1932dc7a7e3SShri Abhyankar /*@C 194*ca4445c7SIlya Fursov TSSetEventHandler - Sets functions and parameters used for indicating events and handling them 1952dc7a7e3SShri Abhyankar 196c3339decSBarry Smith Logically Collective 1972dc7a7e3SShri Abhyankar 1982dc7a7e3SShri Abhyankar Input Parameters: 199bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 200*ca4445c7SIlya Fursov . nevents - number of local events (i.e. managed by the given MPI process) 201*ca4445c7SIlya Fursov . direction - direction of zero crossing to be detected (one for each local event). 202*ca4445c7SIlya Fursov `-1` => zero crossing in negative direction, 203*ca4445c7SIlya Fursov `+1` => zero crossing in positive direction, `0` => both ways 204*ca4445c7SIlya Fursov . terminate - flag to indicate whether time stepping should be terminated after 205*ca4445c7SIlya Fursov an event is detected (one for each local event) 206*ca4445c7SIlya Fursov . indicator - callback defininig the user indicator functions whose sign changes (see `direction`) mark presence of the events 207*ca4445c7SIlya Fursov . postevent - [optional] user post-event callback; it can change the solution, ODE etc at the time of the event 208*ca4445c7SIlya Fursov - ctx - [optional] user-defined context for private data for the 209*ca4445c7SIlya Fursov `indicator()` and `postevent()` routines (use `NULL` if no 210*ca4445c7SIlya Fursov context is desired) 2112dc7a7e3SShri Abhyankar 2129b3d3759SBarry Smith Calling sequence of `indicator`: 2132fe279fdSBarry Smith + ts - the `TS` context 2142dc7a7e3SShri Abhyankar . t - current time 215*ca4445c7SIlya Fursov . U - current solution 216*ca4445c7SIlya Fursov . fvalue - output array with values of local indicator functions (length == `nevents`) for time t and state-vector U 2179b3d3759SBarry Smith - ctx - the context passed as the final argument to `TSSetEventHandler()` 2182dc7a7e3SShri Abhyankar 21920f4b53cSBarry Smith Calling sequence of `postevent`: 2202fe279fdSBarry Smith + ts - the `TS` context 221*ca4445c7SIlya Fursov . nevents_zero - number of triggered local events (whose indicator function is marked as crossing zero, and direction is appropriate) 222*ca4445c7SIlya Fursov . events_zero - indices of the triggered local events 2232dc7a7e3SShri Abhyankar . t - current time 2242dc7a7e3SShri Abhyankar . U - current solution 225*ca4445c7SIlya Fursov . forwardsolve - flag to indicate whether `TS` is doing a forward solve (`PETSC_TRUE`) or adjoint solve (`PETSC_FALSE`) 2269b3d3759SBarry Smith - ctx - the context passed as the final argument to `TSSetEventHandler()` 2272dc7a7e3SShri Abhyankar 228*ca4445c7SIlya Fursov Options Database Keys: 229*ca4445c7SIlya Fursov + -ts_event_tol <tol> - tolerance for zero crossing check of indicator functions 230*ca4445c7SIlya Fursov . -ts_event_monitor - print choices made by event handler 231*ca4445c7SIlya Fursov . -ts_event_recorder_initial_size <recsize> - initial size of event recorder 232*ca4445c7SIlya Fursov . -ts_event_post_event_step <dt> - time step after event 233*ca4445c7SIlya Fursov - -ts_event_dt_min <dt> - minimum time step considered for TSEvent 234*ca4445c7SIlya Fursov 2352dc7a7e3SShri Abhyankar Level: intermediate 2362dc7a7e3SShri Abhyankar 237*ca4445c7SIlya Fursov Notes: 238*ca4445c7SIlya Fursov The indicator functions should be defined in the `indicator` callback using the components of solution `U` and/or time `t`. 239*ca4445c7SIlya Fursov Note that `U` is `PetscScalar`-valued, and the indicator functions are `PetscReal`-valued. It is the user's responsibility to 240*ca4445c7SIlya Fursov properly handle this difference, e.g. by applying `PetscRealPart()` or other appropriate conversion means. 241*ca4445c7SIlya Fursov 242*ca4445c7SIlya Fursov The full set of events is distributed (by the user design) across MPI processes, with each process defining its own local sub-set of events. 243*ca4445c7SIlya Fursov However, the `postevent()` callback invocation is performed synchronously on all processes, including 244*ca4445c7SIlya Fursov those processes which have not currently triggered any events. 245*ca4445c7SIlya Fursov 2461cc06b55SBarry Smith .seealso: [](ch_ts), `TSEvent`, `TSCreate()`, `TSSetTimeStep()`, `TSSetConvergedReason()` 2472dc7a7e3SShri Abhyankar @*/ 248*ca4445c7SIlya Fursov PetscErrorCode TSSetEventHandler(TS ts, PetscInt nevents, PetscInt direction[], PetscBool terminate[], PetscErrorCode (*indicator)(TS ts, PetscReal t, Vec U, PetscReal fvalue[], void *ctx), PetscErrorCode (*postevent)(TS ts, PetscInt nevents_zero, PetscInt events_zero[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx), void *ctx) 249d71ae5a4SJacob Faibussowitsch { 250d294eb03SHong Zhang TSAdapt adapt; 251d294eb03SHong Zhang PetscReal hmin; 2522dc7a7e3SShri Abhyankar TSEvent event; 253006e6a18SShri Abhyankar PetscBool flg; 254a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE 255a6c783d2SShri Abhyankar PetscReal tol = 1e-4; 256a6c783d2SShri Abhyankar #else 257d569cc17SSatish Balay PetscReal tol = 1e-6; 258a6c783d2SShri Abhyankar #endif 2592dc7a7e3SShri Abhyankar 2602dc7a7e3SShri Abhyankar PetscFunctionBegin; 2616427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 2620a82b154SShri if (nevents) { 2634f572ea9SToby Isaac PetscAssertPointer(direction, 3); 2644f572ea9SToby Isaac PetscAssertPointer(terminate, 4); 2650a82b154SShri } 2664dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&event)); 2679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue_prev)); 268*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->fvalue)); 2699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue_right)); 270*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->fsign_prev)); 271*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->fsign)); 272*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->fsign_right)); 2739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->side)); 274*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->side_prev)); 275*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->justrefined_AB)); 276*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->gamma_AB)); 2779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->direction)); 2789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->terminate)); 279*ca4445c7SIlya Fursov PetscCall(PetscMalloc1(nevents, &event->events_zero)); 2809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->vtol)); 281*ca4445c7SIlya Fursov for (PetscInt i = 0; i < nevents; i++) { 282d94325d3SShri Abhyankar event->direction[i] = direction[i]; 283d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 284*ca4445c7SIlya Fursov event->justrefined_AB[i] = PETSC_FALSE; 285*ca4445c7SIlya Fursov event->gamma_AB[i] = 1; 286*ca4445c7SIlya Fursov event->side[i] = 2; 287*ca4445c7SIlya Fursov event->side_prev[i] = 0; 288d94325d3SShri Abhyankar } 289*ca4445c7SIlya Fursov event->iterctr = 0; 290*ca4445c7SIlya Fursov event->processing = PETSC_FALSE; 291*ca4445c7SIlya Fursov event->revisit_right = PETSC_FALSE; 2922589fa24SLisandro Dalcin event->nevents = nevents; 2939b3d3759SBarry Smith event->indicator = indicator; 2942dc7a7e3SShri Abhyankar event->postevent = postevent; 2956427ac75SLisandro Dalcin event->ctx = ctx; 296*ca4445c7SIlya Fursov event->timestep_postevent = ts->time_step; 2979566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 2989566063dSJacob Faibussowitsch PetscCall(TSAdaptGetStepLimits(adapt, &hmin, NULL)); 299d294eb03SHong Zhang event->timestep_min = hmin; 3002dc7a7e3SShri Abhyankar 30102749585SLisandro Dalcin event->recsize = 8; /* Initial size of the recorder */ 302d0609cedSBarry Smith PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "TS Event options", "TS"); 3032dc7a7e3SShri Abhyankar { 304*ca4445c7SIlya Fursov PetscCall(PetscOptionsReal("-ts_event_tol", "Tolerance for zero crossing check of indicator functions", "TSSetEventTolerances", tol, &tol, NULL)); 3059566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_event_monitor", "Print choices made by event handler", "", &flg)); 3069566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size", "Initial size of event recorder", "", event->recsize, &event->recsize, NULL)); 3079566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_post_event_step", "Time step after event", "", event->timestep_postevent, &event->timestep_postevent, NULL)); 308*ca4445c7SIlya Fursov PetscCall(PetscOptionsDeprecated("-ts_event_post_eventinterval_step", NULL, "3.20", "Use -ts_event_post_event_step")); 3099566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_dt_min", "Minimum time step considered for TSEvent", "", event->timestep_min, &event->timestep_min, NULL)); 3102dc7a7e3SShri Abhyankar } 311d0609cedSBarry Smith PetscOptionsEnd(); 312a4ffd976SShri Abhyankar 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.time)); 3149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.stepnum)); 3159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.nevents)); 3169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.eventidx)); 317*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &event->recorder.eventidx[i])); 318e7069c78SShri /* Initialize the event recorder */ 319e7069c78SShri event->recorder.ctr = 0; 320a4ffd976SShri Abhyankar 321*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) event->vtol[i] = tol; 3229566063dSJacob Faibussowitsch if (flg) PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF, "stdout", &event->monitor)); 3239e12be75SShri Abhyankar 3249566063dSJacob Faibussowitsch PetscCall(TSEventDestroy(&ts->event)); 325d94325d3SShri Abhyankar ts->event = event; 326e7069c78SShri ts->event->refct = 1; 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3282dc7a7e3SShri Abhyankar } 3292dc7a7e3SShri Abhyankar 330a4ffd976SShri Abhyankar /* 331a4ffd976SShri Abhyankar TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize) 332a4ffd976SShri Abhyankar is reached. 333a4ffd976SShri Abhyankar */ 334d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventRecorderResize(TSEvent event) 335d71ae5a4SJacob Faibussowitsch { 336a4ffd976SShri Abhyankar PetscReal *time; 337*ca4445c7SIlya Fursov PetscInt *stepnum, *nevents; 338a4ffd976SShri Abhyankar PetscInt **eventidx; 339*ca4445c7SIlya Fursov PetscInt fact = 2; 340a4ffd976SShri Abhyankar 341a4ffd976SShri Abhyankar PetscFunctionBegin; 342*ca4445c7SIlya Fursov /* Create larger arrays */ 3439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &time)); 3449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &stepnum)); 3459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &nevents)); 3469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &eventidx)); 347*ca4445c7SIlya Fursov for (PetscInt i = 0; i < fact * event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &eventidx[i])); 348a4ffd976SShri Abhyankar 349a4ffd976SShri Abhyankar /* Copy over data */ 3509566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(time, event->recorder.time, event->recsize)); 3519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(stepnum, event->recorder.stepnum, event->recsize)); 3529566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nevents, event->recorder.nevents, event->recsize)); 353*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->recsize; i++) PetscCall(PetscArraycpy(eventidx[i], event->recorder.eventidx[i], event->recorder.nevents[i])); 354a4ffd976SShri Abhyankar 355a4ffd976SShri Abhyankar /* Destroy old arrays */ 356*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->recsize; i++) PetscCall(PetscFree(event->recorder.eventidx[i])); 3579566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.eventidx)); 3589566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.nevents)); 3599566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.stepnum)); 3609566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.time)); 361a4ffd976SShri Abhyankar 362a4ffd976SShri Abhyankar /* Set pointers */ 363a4ffd976SShri Abhyankar event->recorder.time = time; 364a4ffd976SShri Abhyankar event->recorder.stepnum = stepnum; 365a4ffd976SShri Abhyankar event->recorder.nevents = nevents; 366a4ffd976SShri Abhyankar event->recorder.eventidx = eventidx; 367a4ffd976SShri Abhyankar 368*ca4445c7SIlya Fursov /* Update the size */ 369a4ffd976SShri Abhyankar event->recsize *= fact; 3703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 371a4ffd976SShri Abhyankar } 372a4ffd976SShri Abhyankar 373031fbad4SShri Abhyankar /* 374*ca4445c7SIlya Fursov Helper routine to handle user postevents and recording 375031fbad4SShri Abhyankar */ 376d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPostEvent(TS ts, PetscReal t, Vec U) 377d71ae5a4SJacob Faibussowitsch { 3782dc7a7e3SShri Abhyankar TSEvent event = ts->event; 37928d5b5d6SLisandro Dalcin PetscBool restart = PETSC_FALSE; 380*ca4445c7SIlya Fursov PetscBool terminate = PETSC_FALSE; 381*ca4445c7SIlya Fursov PetscBool statechanged = PETSC_FALSE; 382*ca4445c7SIlya Fursov PetscInt ctr, stepnum; 383*ca4445c7SIlya Fursov PetscBool inflag[3], outflag[3]; 384*ca4445c7SIlya Fursov PetscBool forwardsolve = PETSC_TRUE; // Flag indicating that TS is doing a forward solve 3852dc7a7e3SShri Abhyankar 3862dc7a7e3SShri Abhyankar PetscFunctionBegin; 3872dc7a7e3SShri Abhyankar if (event->postevent) { 38828d5b5d6SLisandro Dalcin PetscObjectState state_prev, state_post; 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &state_prev)); 390*ca4445c7SIlya Fursov PetscCallBack("TSEvent post-event processing", (*event->postevent)(ts, event->nevents_zero, event->events_zero, t, U, forwardsolve, event->ctx)); // TODO update 'restart' here? 3919566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &state_post)); 392*ca4445c7SIlya Fursov if (state_prev != state_post) { 393*ca4445c7SIlya Fursov restart = PETSC_TRUE; 394*ca4445c7SIlya Fursov statechanged = PETSC_TRUE; 395*ca4445c7SIlya Fursov } 3962dc7a7e3SShri Abhyankar } 3974597913aSLisandro Dalcin 398*ca4445c7SIlya Fursov // Handle termination events and step restart 399*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents_zero; i++) 4009371c9d4SSatish Balay if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE; 4019371c9d4SSatish Balay inflag[0] = restart; 4029371c9d4SSatish Balay inflag[1] = terminate; 403*ca4445c7SIlya Fursov inflag[2] = statechanged; 404*ca4445c7SIlya Fursov PetscCall(MPIU_Allreduce(inflag, outflag, 3, MPIU_BOOL, MPI_LOR, ((PetscObject)ts)->comm)); 4059371c9d4SSatish Balay restart = outflag[0]; 4069371c9d4SSatish Balay terminate = outflag[1]; 407*ca4445c7SIlya Fursov statechanged = outflag[2]; 4089566063dSJacob Faibussowitsch if (restart) PetscCall(TSRestartStep(ts)); 4099566063dSJacob Faibussowitsch if (terminate) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_EVENT)); 410f7aea88cSShri Abhyankar 411*ca4445c7SIlya Fursov /* 412*ca4445c7SIlya Fursov Recalculate the indicator functions and signs if the state has been changed by the user postevent callback. 413*ca4445c7SIlya Fursov Note! If the state HAS NOT changed, the existing event->fsign (equal to zero) is kept, which: 414*ca4445c7SIlya Fursov - might have been defined using the previous (now-possibly-overridden) event->vtol, 415*ca4445c7SIlya Fursov - might have been set to zero on reaching a small time step rather than using the vtol criterion. 416*ca4445c7SIlya Fursov This will enforce keeping event->fsign = 0 where the zero-crossings were actually marked, 417*ca4445c7SIlya Fursov resulting in a more consistent behaviour of fsign's. 418*ca4445c7SIlya Fursov */ 419*ca4445c7SIlya Fursov if (statechanged) { 420*ca4445c7SIlya Fursov if (event->monitor) PetscCall(PetscPrintf(((PetscObject)ts)->comm, "TSEvent: at time %g the vector state has been changed by PostEvent, recalculating fvalues and signs\n", (double)t)); 4219566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(U)); 422*ca4445c7SIlya Fursov PetscCallBack("TSEvent indicator", (*event->indicator)(ts, t, U, event->fvalue, event->ctx)); 4239566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(U)); 424*ca4445c7SIlya Fursov TSEventCalcSigns(event->nevents, event->fvalue, event->vtol, event->fsign); // note, event->vtol might have been changed by the postevent() 425f443add6SLisandro Dalcin } 426f443add6SLisandro Dalcin 427*ca4445c7SIlya Fursov // Record the event in the event recorder 4289566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &stepnum)); 429f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 43048a46eb9SPierre Jolivet if (ctr == event->recsize) PetscCall(TSEventRecorderResize(event)); 431f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 432d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 4334597913aSLisandro Dalcin event->recorder.nevents[ctr] = event->nevents_zero; 434*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 435f7aea88cSShri Abhyankar event->recorder.ctr++; 4363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4372dc7a7e3SShri Abhyankar } 4382dc7a7e3SShri Abhyankar 439*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-sowing-chars 440*ca4445c7SIlya Fursov /* 441*ca4445c7SIlya Fursov (modified) Anderson-Bjorck variant of regula falsi method, refines [tleft, t] or [t, tright] based on 'side' (-1 or +1). 442*ca4445c7SIlya Fursov The scaling parameter is defined based on the 'justrefined' flag, the history of repeats of 'side', and the threshold. 443*ca4445c7SIlya Fursov To escape certain failure modes, the algorithm may drift towards the bisection rule. 444*ca4445c7SIlya Fursov The value pointed to by 'side_prev' gets updated. 445*ca4445c7SIlya Fursov This function returns the new time step. 446*ca4445c7SIlya Fursov 447*ca4445c7SIlya Fursov The underlying pure Anderson-Bjorck algorithm was taken as described in 448*ca4445c7SIlya Fursov J.M. Fernandez-Diaz, C.O. Menendez-Perez "A common framework for modified Regula Falsi methods and new methods of this kind", 2023. 449*ca4445c7SIlya Fursov The modifications subsequently introduced have little effect on the behaviour for simple cases requiring only a few iterations 450*ca4445c7SIlya Fursov (some minor convergence slowdown may take place though), but the effect becomes more pronounced for tough cases requiring many iterations. 451*ca4445c7SIlya Fursov For the latter situation the speed-up may be order(s) of magnitude compared to the classical Anderson-Bjorck. 452*ca4445c7SIlya Fursov The modifications (the threshold trick, and the drift towards bisection) were tested and tweaked 453*ca4445c7SIlya Fursov based on a number of test functions from the mentioned paper. 454*ca4445c7SIlya Fursov */ 455*ca4445c7SIlya Fursov static inline PetscReal RefineAndersonBjorck(PetscReal tleft, PetscReal t, PetscReal tright, PetscReal fleft, PetscReal f, PetscReal fright, PetscInt side, PetscInt *side_prev, PetscBool justrefined, PetscReal *gamma) 456d71ae5a4SJacob Faibussowitsch { 457*ca4445c7SIlya Fursov PetscReal new_dt, scal = 1.0, scalB = 1.0, threshold = 0.0, power; 458*ca4445c7SIlya Fursov PetscInt reps = 0; 459*ca4445c7SIlya Fursov const PetscInt REPS_CAP = 8; // an upper bound to be imposed on 'reps' (set to 8, somewhat arbitrary number, found after some tweaking) 460*ca4445c7SIlya Fursov 461*ca4445c7SIlya Fursov // Preparations 462*ca4445c7SIlya Fursov if (justrefined) { 463*ca4445c7SIlya Fursov if (*side_prev * side > 0) *side_prev += side; // the side keeps repeating -> increment the side counter (-ve or +ve) 464*ca4445c7SIlya Fursov else *side_prev = side; // reset the counter 465*ca4445c7SIlya Fursov reps = PetscMin(*side_prev * side, REPS_CAP); // the length of the recent side-repeat series, including the current 'side' 466*ca4445c7SIlya Fursov threshold = PetscPowReal(0.5, reps) * 0.1; // ad-hoc strategy for threshold calculation (involved some tweaking) 467*ca4445c7SIlya Fursov } else *side_prev = side; // initial reset of the counter 468*ca4445c7SIlya Fursov 469*ca4445c7SIlya Fursov // Time step calculation 470e2cdd850SShri Abhyankar if (side == -1) { 471*ca4445c7SIlya Fursov if (justrefined && fright != 0.0 && fleft != 0.0) { 472*ca4445c7SIlya Fursov scal = (fright - f) / fright; 473*ca4445c7SIlya Fursov scalB = -f / fleft; 474e2cdd850SShri Abhyankar } 475*ca4445c7SIlya Fursov } else { // must be side == +1 476*ca4445c7SIlya Fursov if (justrefined && fleft != 0.0 && fright != 0.0) { 477*ca4445c7SIlya Fursov scal = (fleft - f) / fleft; 478*ca4445c7SIlya Fursov scalB = -f / fright; 479e2cdd850SShri Abhyankar } 48038bf2713SShri Abhyankar } 481e2cdd850SShri Abhyankar 482*ca4445c7SIlya Fursov if (scal < threshold) scal = 0.5; 483*ca4445c7SIlya Fursov if (reps > 1) *gamma *= scal; // side did not switch since the last time, accumulate gamma 484*ca4445c7SIlya Fursov else *gamma = 1.0; // side switched -> reset gamma 485*ca4445c7SIlya Fursov power = PetscMax(0.0, (reps - 2.0) / (REPS_CAP - 2.0)); 486*ca4445c7SIlya Fursov scal = PetscPowReal(scalB / *gamma, power) * (*gamma); // mix the Anderson-Bjorck scaling and Bisection scaling 487*ca4445c7SIlya Fursov 488*ca4445c7SIlya Fursov if (side == -1) new_dt = scal * fleft / (scal * fleft - f) * (t - tleft); 489*ca4445c7SIlya Fursov else new_dt = f / (f - scal * fright) * (tright - t); 490*ca4445c7SIlya Fursov /* 491*ca4445c7SIlya Fursov In tough cases (e.g. a polynomial of high order), there is a failure mode for the standard Anderson-Bjorck, 492*ca4445c7SIlya Fursov when the new proposed point jumps from one end-point of the bracket to the other, however the bracket 493*ca4445c7SIlya Fursov is contracting very slowly. A larger threshold for 'scal' prevents entering this mode. 494*ca4445c7SIlya Fursov On the other hand, if the iteration gets stuck near one end-point of the bracket, and the 'side' does not switch for a while, 495*ca4445c7SIlya Fursov the 'scal' drifts towards the bisection approach (via scalB), ensuring stable convergence. 496*ca4445c7SIlya Fursov */ 497*ca4445c7SIlya Fursov return new_dt; 498*ca4445c7SIlya Fursov } 499*ca4445c7SIlya Fursov 500*ca4445c7SIlya Fursov /* 501*ca4445c7SIlya Fursov Checks if the current point (t) is the zero-crossing location, based on the indicator function signs and direction[]: 502*ca4445c7SIlya Fursov - using the dt_min criterion, 503*ca4445c7SIlya Fursov - using the vtol criterion. 504*ca4445c7SIlya Fursov The situation (fsign_prev, fsign) = (0, 0) is treated as staying in the near-zero-zone of the previous zero-crossing, 505*ca4445c7SIlya Fursov and is not marked as a new zero-crossing. 506*ca4445c7SIlya Fursov This function may update event->side[]. 507*ca4445c7SIlya Fursov */ 508*ca4445c7SIlya Fursov static PetscErrorCode TSEventTestZero(TS ts, PetscReal t) 509d71ae5a4SJacob Faibussowitsch { 510d294eb03SHong Zhang TSEvent event = ts->event; 511d294eb03SHong Zhang 512d294eb03SHong Zhang PetscFunctionBegin; 513*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) { 514*ca4445c7SIlya Fursov const PetscBool bracket_is_left = (event->fsign_prev[i] * event->fsign[i] < 0 && event->fsign[i] * event->direction[i] >= 0) ? PETSC_TRUE : PETSC_FALSE; 515d294eb03SHong Zhang 516*ca4445c7SIlya Fursov if (bracket_is_left && ((t - event->ptime_prev <= event->timestep_min) || event->revisit_right)) event->side[i] = 0; // mark zero-crossing from dt_min; 'bracket_is_left' accounts for direction 517*ca4445c7SIlya Fursov if (event->fsign[i] == 0 && event->fsign_prev[i] != 0 && event->fsign_prev[i] * event->direction[i] <= 0) event->side[i] = 0; // mark zero-crossing from vtol 518d294eb03SHong Zhang } 5193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 520d294eb03SHong Zhang } 52138bf2713SShri Abhyankar 522*ca4445c7SIlya Fursov /* 523*ca4445c7SIlya Fursov Checks if [fleft, f] or [f, fright] are 'brackets', i.e. intervals with the sign change, satisfying the 'direction'. 524*ca4445c7SIlya Fursov The right interval is only checked if iterctr > 0 (i.e. Anderson-Bjorck refinement has started). 525*ca4445c7SIlya Fursov The intervals like [0, x] and [x, 0] are not counted as brackets, i.e. intervals with the sign change. 526*ca4445c7SIlya Fursov The function returns the 'side' value: -1 (left, or both are brackets), +1 (only right one), +2 (neither). 527*ca4445c7SIlya Fursov */ 528*ca4445c7SIlya Fursov static inline PetscInt TSEventTestBracket(PetscInt fsign_left, PetscInt fsign, PetscInt fsign_right, PetscInt direction, PetscInt iterctr) 529*ca4445c7SIlya Fursov { 530*ca4445c7SIlya Fursov PetscInt side = 2; 531*ca4445c7SIlya Fursov if (fsign_left * fsign < 0 && fsign * direction >= 0) side = -1; 532*ca4445c7SIlya Fursov if (side != -1 && iterctr > 0 && fsign * fsign_right < 0 && fsign_right * direction >= 0) side = 1; 533*ca4445c7SIlya Fursov return side; 534*ca4445c7SIlya Fursov } 535*ca4445c7SIlya Fursov 536*ca4445c7SIlya Fursov /* 537*ca4445c7SIlya Fursov Caps the time steps, accounting for time span points. 538*ca4445c7SIlya Fursov It uses 'event->timestep_cache' as a time step to calculate the tolerance for tspan points detection. This 539*ca4445c7SIlya Fursov is done since the event resolution may result in significant time step refinement, and we don't use these small steps for tolerances. 540*ca4445c7SIlya Fursov To enhance the consistency of tspan points detection, tolerance 'tspan->worktol' is reused later in the TSSolve iteration. 541*ca4445c7SIlya Fursov If a user-defined step is cut by this function, the input uncut step is saved to adapt->dt_span_cached. 542*ca4445c7SIlya Fursov Flag 'user_dt' indicates if the step was defined by user. 543*ca4445c7SIlya Fursov */ 544*ca4445c7SIlya Fursov static inline PetscReal TSEvent_dt_cap(TS ts, PetscReal t, PetscReal dt, PetscBool user_dt) 545*ca4445c7SIlya Fursov { 546*ca4445c7SIlya Fursov PetscReal res = dt; 547*ca4445c7SIlya Fursov if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 548*ca4445c7SIlya Fursov PetscReal maxdt = ts->max_time - t; // this may be overriden by tspan 549*ca4445c7SIlya Fursov PetscBool cut_made = PETSC_FALSE; 550*ca4445c7SIlya Fursov PetscReal eps = 10 * PETSC_MACHINE_EPSILON; 551*ca4445c7SIlya Fursov if (ts->tspan) { 552*ca4445c7SIlya Fursov PetscInt ctr = ts->tspan->spanctr; 553*ca4445c7SIlya Fursov PetscInt Ns = ts->tspan->num_span_times; 554*ca4445c7SIlya Fursov PetscReal *st = ts->tspan->span_times; 555*ca4445c7SIlya Fursov 556*ca4445c7SIlya Fursov if (ts->tspan->worktol == 0) ts->tspan->worktol = ts->tspan->reltol * ts->event->timestep_cache + ts->tspan->abstol; // in case TSAdaptChoose() has not defined it 557*ca4445c7SIlya Fursov if (ctr < Ns && PetscIsCloseAtTol(t, st[ctr], ts->tspan->worktol, 0)) { // just hit a time span point 558*ca4445c7SIlya Fursov if (ctr + 1 < Ns) maxdt = st[ctr + 1] - t; // ok to use the next time span point 559*ca4445c7SIlya Fursov else maxdt = ts->max_time - t; // can't use the next time span point: they have finished 560*ca4445c7SIlya Fursov } else if (ctr < Ns) maxdt = st[ctr] - t; // haven't hit a time span point, use the nearest one 561*ca4445c7SIlya Fursov } 562*ca4445c7SIlya Fursov maxdt = PetscMin(maxdt, ts->max_time - t); 563*ca4445c7SIlya Fursov PetscCheck((maxdt > eps) || (PetscAbsReal(maxdt) <= eps && PetscIsCloseAtTol(t, ts->max_time, eps, 0)), PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Unexpected state: bad maxdt in TSEvent_dt_cap()"); 564*ca4445c7SIlya Fursov 565*ca4445c7SIlya Fursov if (PetscIsCloseAtTol(dt, maxdt, eps, 0)) res = maxdt; // no cut 566*ca4445c7SIlya Fursov else { 567*ca4445c7SIlya Fursov if (dt > maxdt) { 568*ca4445c7SIlya Fursov res = maxdt; // yes cut 569*ca4445c7SIlya Fursov cut_made = PETSC_TRUE; 570*ca4445c7SIlya Fursov } else res = dt; // no cut 571*ca4445c7SIlya Fursov } 572*ca4445c7SIlya Fursov if (ts->adapt && user_dt) { // only update dt_span_cached for the user-defined step 573*ca4445c7SIlya Fursov if (cut_made) ts->adapt->dt_span_cached = dt; 574*ca4445c7SIlya Fursov else ts->adapt->dt_span_cached = 0; 575*ca4445c7SIlya Fursov } 576*ca4445c7SIlya Fursov } 577*ca4445c7SIlya Fursov return res; 578*ca4445c7SIlya Fursov } 579*ca4445c7SIlya Fursov 580*ca4445c7SIlya Fursov /* 581*ca4445c7SIlya Fursov Updates the left-end values 582*ca4445c7SIlya Fursov */ 583*ca4445c7SIlya Fursov static inline void TSEvent_update_left(TSEvent event, PetscReal t) 584*ca4445c7SIlya Fursov { 585*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) { 586*ca4445c7SIlya Fursov event->fvalue_prev[i] = event->fvalue[i]; 587*ca4445c7SIlya Fursov event->fsign_prev[i] = event->fsign[i]; 588*ca4445c7SIlya Fursov } 589*ca4445c7SIlya Fursov event->ptime_prev = t; 590*ca4445c7SIlya Fursov } 591*ca4445c7SIlya Fursov 592*ca4445c7SIlya Fursov /* 593*ca4445c7SIlya Fursov Updates the right-end values 594*ca4445c7SIlya Fursov */ 595*ca4445c7SIlya Fursov static inline void TSEvent_update_right(TSEvent event, PetscReal t) 596*ca4445c7SIlya Fursov { 597*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) { 598*ca4445c7SIlya Fursov event->fvalue_right[i] = event->fvalue[i]; 599*ca4445c7SIlya Fursov event->fsign_right[i] = event->fsign[i]; 600*ca4445c7SIlya Fursov } 601*ca4445c7SIlya Fursov event->ptime_right = t; 602*ca4445c7SIlya Fursov } 603*ca4445c7SIlya Fursov 604*ca4445c7SIlya Fursov /* 605*ca4445c7SIlya Fursov Updates the current values from the right-end values 606*ca4445c7SIlya Fursov */ 607*ca4445c7SIlya Fursov static inline PetscReal TSEvent_update_from_right(TSEvent event) 608*ca4445c7SIlya Fursov { 609*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) { 610*ca4445c7SIlya Fursov event->fvalue[i] = event->fvalue_right[i]; 611*ca4445c7SIlya Fursov event->fsign[i] = event->fsign_right[i]; 612*ca4445c7SIlya Fursov } 613*ca4445c7SIlya Fursov return event->ptime_right; 614*ca4445c7SIlya Fursov } 615*ca4445c7SIlya Fursov 616*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-spacing 617*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-unknown 618*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-spelling 619*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-missing 620*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-section-header-fishy-header 621*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation 622*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-synopsis-missing-description 623*ca4445c7SIlya Fursov // PetscClangLinter pragma disable: -fdoc-sowing-chars 624*ca4445c7SIlya Fursov /* 625*ca4445c7SIlya Fursov TSEventHandler - the main function to perform a single iteration of event detection. 626*ca4445c7SIlya Fursov 627*ca4445c7SIlya Fursov Developer notes: 628*ca4445c7SIlya Fursov a) The 'event->iterctr > 0' is used as an indicator that Anderson-Bjorck refinement has started. 629*ca4445c7SIlya Fursov b) If event->iterctr == 0, then justrefined_AB[i] is always false. 630*ca4445c7SIlya Fursov c) The right-end quantities: ptime_right, fvalue_right[i] and fsign_right[i] are only guaranteed to be valid 631*ca4445c7SIlya Fursov for event->iterctr > 0. 632*ca4445c7SIlya Fursov d) If event->iterctr > 0, then event->processing is PETSC_TRUE; the opposite may not hold. 633*ca4445c7SIlya Fursov e) event->side[i] may take values: 0 <=> point t is a zero-crossing for indicator function i (via vtol/dt_min criterion); 634*ca4445c7SIlya Fursov -1/+1 <=> detected a bracket to the left/right of t for indicator function i; +2 <=> no brackets/zero-crossings. 635*ca4445c7SIlya Fursov f) The signs event->fsign[i] (with values 0/-1/+1) are calculated for each new point. Zero sign is set if the function value is 636*ca4445c7SIlya Fursov smaller than the tolerance. Besides, zero sign is enforced after marking a zero-crossing due to small bracket size criterion. 637*ca4445c7SIlya Fursov 638*ca4445c7SIlya Fursov The intervals with the indicator function sign change (i.e. containing the potential zero-crossings) are called 'brackets'. 639*ca4445c7SIlya Fursov To find a zero-crossing, the algorithm first locates a bracket, and then sequentially subdivides it, generating a sequence 640*ca4445c7SIlya Fursov of brackets whose length tends to zero. The bracket subdivision involves the (modified) Anderson-Bjorck method. 641*ca4445c7SIlya Fursov 642*ca4445c7SIlya Fursov Apart from the comments scattered throughout the code to clarify different lines and blocks, 643*ca4445c7SIlya Fursov a few tricky aspects of the algorithm (and the underlying reasoning) are discussed in detail below: 644*ca4445c7SIlya Fursov 645*ca4445c7SIlya Fursov =Sign tracking= 646*ca4445c7SIlya Fursov When a zero-crossing is found, the sign variable (event->fsign[i]) is set to zero for the current point t. 647*ca4445c7SIlya Fursov This happens both for zero-crossings triggered via the vtol criterion, and those triggered via the dt_min 648*ca4445c7SIlya Fursov criterion. After the event, as the TS steps forward, the current sign values are handed over to event->fsign_prev[i]. 649*ca4445c7SIlya Fursov The recalculation of signs is avoided if possible: e.g. if a 'vtol' criterion resulted in a zero-crossing at point t, 650*ca4445c7SIlya Fursov but the subsequent call to postevent() handler decreased 'vtol', making the indicator function no longer "close to zero" 651*ca4445c7SIlya Fursov at point t, the fsign[i] will still consistently keep the zero value. This allows avoiding the erroneous duplication 652*ca4445c7SIlya Fursov of events: 653*ca4445c7SIlya Fursov 654*ca4445c7SIlya Fursov E.g. consider a bracket [t0, t2], where f0 < 0, f2 > 0, which resulted in a zero-crossing t1 with f1 < 0, abs(f1) < vtol. 655*ca4445c7SIlya Fursov Suppose the postevent() handler changes vtol to vtol*, such that abs(f1) > vtol*. The TS makes a step t1 -> t3, where 656*ca4445c7SIlya Fursov again f1 < 0, f3 > 0, and the event handler will find a new event near t1, which is actually a duplication of the 657*ca4445c7SIlya Fursov original event at t1. The duplications are avoided by NOT counting the sign progressions 0 -> +1, or 0 -> -1 658*ca4445c7SIlya Fursov as brackets. Tracking (instead of recalculating) the sign values makes this procedure work more consistently. 659*ca4445c7SIlya Fursov 660*ca4445c7SIlya Fursov The sign values are however recalculated if the postevent() callback has changed the current solution vector U 661*ca4445c7SIlya Fursov (such a change resets everything). 662*ca4445c7SIlya Fursov The sign value is also set to zero if the dt_min criterion has triggered the event. This allows the algorithm to 663*ca4445c7SIlya Fursov work consistently, irrespective of the type of criterion involved (vtol/dt_min). 664*ca4445c7SIlya Fursov 665*ca4445c7SIlya Fursov =Event from min bracket= 666*ca4445c7SIlya Fursov When the event handler ends up with a bracket [t0, t1] with size <= dt_min, a zero crossing is reported at t1, 667*ca4445c7SIlya Fursov and never at t0. If such a bracket is discovered when TS is staying at t0, one more step forward (to t1) is necessary 668*ca4445c7SIlya Fursov to mark the found event. This is the situation of revisiting t1, which is described below (see Revisiting). 669*ca4445c7SIlya Fursov 670*ca4445c7SIlya Fursov Why t0 is not reported as event location? Suppose it is, and let f0 < 0, f1 > 0. Also suppose that the 671*ca4445c7SIlya Fursov postevent() handler has slightly changed the solution U, so the sign at t0 is recalculated: it equals -1. As the TS steps 672*ca4445c7SIlya Fursov further: t0 -> t2, with sign0 == -1, and sign2 == +1, the event handler will locate the bracket [t0, t2], eventually 673*ca4445c7SIlya Fursov resolving a new event near t1, i.e. finding a duplicate event. 674*ca4445c7SIlya Fursov This situation is avoided by reporting the event at t1 in the first place. 675*ca4445c7SIlya Fursov 676*ca4445c7SIlya Fursov =Revisiting= 677*ca4445c7SIlya Fursov When handling the situation with small bracket size, the TS solver may happen to visit the same point twice, 678*ca4445c7SIlya Fursov but with different results. 679*ca4445c7SIlya Fursov 680*ca4445c7SIlya Fursov E.g. originally it discovered a bracket with sign change [t0, t10], and started resolving the zero-crossing, 681*ca4445c7SIlya Fursov visiting the points t1,...,t9 : t0 < t1 < ... < t9 < t10. Suppose that at t9 the algorithm discovers 682*ca4445c7SIlya Fursov that [t9, t10] is a bracket with the sign change it was looking for, and that |t10 - t9| is too small. 683*ca4445c7SIlya Fursov So point t10 should be revisited and marked as the zero crossing (by the minimum bracket size criterion). 684*ca4445c7SIlya Fursov On re-visiting t10, via the refined sequence of steps t0,...,t10, the TS solver may arrive at a solution U* 685*ca4445c7SIlya Fursov different from the solution U it found at t10 originally. Hence, the indicator functions at t10 may become different, 686*ca4445c7SIlya Fursov and the condition of the sign change, which existed originally, may disappear, breaking the logic of the algorithm. 687*ca4445c7SIlya Fursov 688*ca4445c7SIlya Fursov To handle such (-=unlikely=-, but possible) situations, two strategies can be considered: 689*ca4445c7SIlya Fursov 1) [not used here] Allow the brackets with sign change to disappear during iterations. The algorithm should be able 690*ca4445c7SIlya Fursov to cleanly exit the iteration and leave all the objects/variables/caches involved in a valid state. 691*ca4445c7SIlya Fursov 2) [ADOPTED HERE!] On revisiting t10, the event handler reuses the indicator functions previously calculated for the 692*ca4445c7SIlya Fursov original solution U. This U may be less precise than U*, but this trick does not allow the algorithm logic to break down. 693*ca4445c7SIlya Fursov HOWEVER, the original U is not stored anywhere, it is essentially lost since the TS performed the rollback from it. 694*ca4445c7SIlya Fursov On revisiting t10, the updated solution U* will inevitably be found and used everywhere EXCEPT the current 695*ca4445c7SIlya Fursov indicator functions calculation, e.g. U* will be used in the postevent() handler call. Since t10 is the event location, 696*ca4445c7SIlya Fursov the appropriate indicator-function-signs will be enforced to be 0 (regardless if the solution was U or U*). 697*ca4445c7SIlya Fursov If the solution is then changed by the postevent(), the indicator-function-signs will be recalculated. 698*ca4445c7SIlya Fursov 699*ca4445c7SIlya Fursov Whether the algorithm is revisiting a point in the current TSEventHandler() call is flagged by 'event->revisit_right'. 700*ca4445c7SIlya Fursov */ 701d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventHandler(TS ts) 702d71ae5a4SJacob Faibussowitsch { 7036427ac75SLisandro Dalcin TSEvent event; 704*ca4445c7SIlya Fursov PetscReal t, dt_next = 0.0; 7052dc7a7e3SShri Abhyankar Vec U; 706*ca4445c7SIlya Fursov PetscInt minsidein = 2, minsideout = 2; // minsideout is sync on all ranks 707*ca4445c7SIlya Fursov PetscBool finished = PETSC_FALSE; // should stay sync on all ranks 708*ca4445c7SIlya Fursov PetscBool revisit_right_cache; // [sync] flag for inner consistency checks 7092dc7a7e3SShri Abhyankar 7102dc7a7e3SShri Abhyankar PetscFunctionBegin; 7116427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 712*ca4445c7SIlya Fursov 7133ba16761SJacob Faibussowitsch if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS); 7146427ac75SLisandro Dalcin event = ts->event; 715*ca4445c7SIlya Fursov event->nevents_zero = 0; 716*ca4445c7SIlya Fursov revisit_right_cache = event->revisit_right; 717*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) event->side[i] = 2; // side's are reset on each new iteration 718*ca4445c7SIlya Fursov if (event->iterctr == 0) 719*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) event->justrefined_AB[i] = PETSC_FALSE; 7202dc7a7e3SShri Abhyankar 7219566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 722*ca4445c7SIlya Fursov if (!event->processing) { // update the caches 723*ca4445c7SIlya Fursov PetscReal dt; 7249566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 725*ca4445c7SIlya Fursov event->ptime_cache = t; 726*ca4445c7SIlya Fursov event->timestep_cache = dt; // the next TS move is planned to be: t -> t+dt 7272dc7a7e3SShri Abhyankar } 7282dc7a7e3SShri Abhyankar 729*ca4445c7SIlya Fursov PetscCall(TSGetSolution(ts, &U)); // if revisiting, this will be the updated U* (see discussion on "Revisiting" in the Developer notes above) 730*ca4445c7SIlya Fursov if (event->revisit_right) { 731*ca4445c7SIlya Fursov PetscReal tr = TSEvent_update_from_right(event); 732*ca4445c7SIlya Fursov PetscCheck(PetscAbsReal(tr - t) < PETSC_SMALL, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Inconsistent time value when performing 'revisiting' in TSEventHandler()"); 7332dc7a7e3SShri Abhyankar } else { 734*ca4445c7SIlya Fursov PetscCall(VecLockReadPush(U)); 735*ca4445c7SIlya Fursov PetscCallBack("TSEvent indicator", (*event->indicator)(ts, t, U, event->fvalue, event->ctx)); // fill fvalue's at point 't' 736*ca4445c7SIlya Fursov PetscCall(VecLockReadPop(U)); 737*ca4445c7SIlya Fursov TSEventCalcSigns(event->nevents, event->fvalue, event->vtol, event->fsign); // fill fvalue signs 738*ca4445c7SIlya Fursov } 739*ca4445c7SIlya Fursov PetscCall(TSEventTestZero(ts, t)); // check if the current point 't' is the event location; event->side[] may get updated 740*ca4445c7SIlya Fursov 741*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) { // check for brackets on the left/right of 't' 742*ca4445c7SIlya Fursov if (event->side[i] != 0) event->side[i] = TSEventTestBracket(event->fsign_prev[i], event->fsign[i], event->fsign_right[i], event->direction[i], event->iterctr); 743*ca4445c7SIlya Fursov minsidein = PetscMin(minsidein, event->side[i]); 744*ca4445c7SIlya Fursov } 745*ca4445c7SIlya Fursov PetscCall(MPIU_Allreduce(&minsidein, &minsideout, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)ts))); 746*ca4445c7SIlya Fursov /* 747*ca4445c7SIlya Fursov minsideout (sync on all ranks) indicates the minimum of the following states: 748*ca4445c7SIlya Fursov -1 : [ptime_prev, t] is a bracket for some indicator-function-i 749*ca4445c7SIlya Fursov +1 : [t, ptime_right] is a bracket for some indicator-function-i 750*ca4445c7SIlya Fursov 0 : t is a zero-crossing for some indicator-function-i 751*ca4445c7SIlya Fursov 2 : none of the above 752*ca4445c7SIlya Fursov */ 753*ca4445c7SIlya Fursov PetscCheck(!event->revisit_right || minsideout == 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "minsideout != 0 when performing 'revisiting' in TSEventHandler()"); 754*ca4445c7SIlya Fursov 755*ca4445c7SIlya Fursov if (minsideout == -1 || minsideout == +1) { // this if-branch will refine the left/right bracket 756*ca4445c7SIlya Fursov const PetscReal bracket_size = (minsideout == -1) ? t - event->ptime_prev : event->ptime_right - t; // sync on all ranks 757*ca4445c7SIlya Fursov 758*ca4445c7SIlya Fursov if (minsideout == +1 && bracket_size <= event->timestep_min) { // check if the bracket (right) is small 759*ca4445c7SIlya Fursov // [--------------------|-] 760*ca4445c7SIlya Fursov dt_next = bracket_size; // need one more step to get to event->ptime_right 761*ca4445c7SIlya Fursov event->revisit_right = PETSC_TRUE; 762*ca4445c7SIlya Fursov TSEvent_update_left(event, t); 763*ca4445c7SIlya Fursov if (event->monitor) 764*ca4445c7SIlya Fursov PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - reached too small bracket [%g - %g], next stepping to its right end %g (revisiting)\n", PetscGlobalRank, event->iterctr, (double)event->ptime_prev, 765*ca4445c7SIlya Fursov (double)event->ptime_right, (double)(event->ptime_prev + dt_next))); 766*ca4445c7SIlya Fursov } else { // the bracket is not very small -> refine it 767*ca4445c7SIlya Fursov // [--------|-------------] 768*ca4445c7SIlya Fursov if (bracket_size <= 2 * event->timestep_min) dt_next = bracket_size / 2; // the bracket is almost small -> bisect it 769*ca4445c7SIlya Fursov else { // the bracket is not small -> use Anderson-Bjorck 770*ca4445c7SIlya Fursov PetscReal dti_min = PETSC_MAX_REAL; 771*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) { 772*ca4445c7SIlya Fursov if (event->side[i] == minsideout) { // only refine the appropriate brackets 773*ca4445c7SIlya Fursov PetscReal dti = RefineAndersonBjorck(event->ptime_prev, t, event->ptime_right, event->fvalue_prev[i], event->fvalue[i], event->fvalue_right[i], event->side[i], &event->side_prev[i], event->justrefined_AB[i], &event->gamma_AB[i]); 774*ca4445c7SIlya Fursov dti_min = PetscMin(dti_min, dti); 775*ca4445c7SIlya Fursov } 776*ca4445c7SIlya Fursov } 777*ca4445c7SIlya Fursov PetscCall(MPIU_Allreduce(&dti_min, &dt_next, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 778*ca4445c7SIlya Fursov if (dt_next < event->timestep_min) dt_next = event->timestep_min; 779*ca4445c7SIlya Fursov if (bracket_size - dt_next < event->timestep_min) dt_next = bracket_size - event->timestep_min; 780*ca4445c7SIlya Fursov } 781*ca4445c7SIlya Fursov 782*ca4445c7SIlya Fursov if (minsideout == -1) { // minsideout == -1, update the right-end values, retain the left-end values 783*ca4445c7SIlya Fursov TSEvent_update_right(event, t); 784*ca4445c7SIlya Fursov PetscCall(TSRollBack(ts)); 785*ca4445c7SIlya Fursov PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); // e.g. to override TS_CONVERGED_TIME on reaching ts->max_time 786*ca4445c7SIlya Fursov } else TSEvent_update_left(event, t); // minsideout == +1, update the left-end values, retain the right-end values 787*ca4445c7SIlya Fursov 788*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) { // update the "Anderson-Bjorck" flags 789*ca4445c7SIlya Fursov if (event->side[i] == minsideout) { 790*ca4445c7SIlya Fursov event->justrefined_AB[i] = PETSC_TRUE; // only for these i's Anderson-Bjorck was invoked 791*ca4445c7SIlya Fursov if (event->monitor) 792*ca4445c7SIlya Fursov PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " refining the bracket with sign change [%g - %g], next stepping to %g\n", PetscGlobalRank, event->iterctr, i, (double)event->ptime_prev, 793*ca4445c7SIlya Fursov (double)event->ptime_right, (double)(event->ptime_prev + dt_next))); 794*ca4445c7SIlya Fursov } else event->justrefined_AB[i] = PETSC_FALSE; // for these i's Anderson-Bjorck was not invoked 795*ca4445c7SIlya Fursov } 796*ca4445c7SIlya Fursov } 797*ca4445c7SIlya Fursov event->iterctr++; 798*ca4445c7SIlya Fursov event->processing = PETSC_TRUE; 799*ca4445c7SIlya Fursov } else if (minsideout == 0) { // found the appropriate zero-crossing (and no brackets to the left), finishing! 800*ca4445c7SIlya Fursov // [--------0-------------] 801*ca4445c7SIlya Fursov finished = PETSC_TRUE; 802*ca4445c7SIlya Fursov event->revisit_right = PETSC_FALSE; 803*ca4445c7SIlya Fursov for (PetscInt i = 0; i < event->nevents; i++) 804*ca4445c7SIlya Fursov if (event->side[i] == minsideout) { 805*ca4445c7SIlya Fursov event->events_zero[event->nevents_zero++] = i; 806*ca4445c7SIlya Fursov if (event->fsign[i] == 0) { // vtol was engaged 807*ca4445c7SIlya Fursov if (event->monitor) 808*ca4445c7SIlya Fursov PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " zero crossing located at time %g (tol=%g)\n", PetscGlobalRank, event->iterctr, i, (double)t, (double)event->vtol[i])); 809*ca4445c7SIlya Fursov } else { // dt_min was engaged 810*ca4445c7SIlya Fursov event->fsign[i] = 0; // sign = 0 is enforced further 811*ca4445c7SIlya Fursov if (event->monitor) 812*ca4445c7SIlya Fursov PetscCall(PetscViewerASCIIPrintf(event->monitor, "[%d] TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " accepting time %g as event location, due to reaching too small bracket [%g - %g]\n", PetscGlobalRank, event->iterctr, i, (double)t, 813*ca4445c7SIlya Fursov (double)event->ptime_prev, (double)t)); 814*ca4445c7SIlya Fursov } 815*ca4445c7SIlya Fursov } 816*ca4445c7SIlya Fursov event->iterctr++; 817*ca4445c7SIlya Fursov event->processing = PETSC_TRUE; 818*ca4445c7SIlya Fursov } else { // minsideout == 2: no brackets, no zero-crossings 819*ca4445c7SIlya Fursov // [----------------------] 820*ca4445c7SIlya Fursov PetscCheck(event->iterctr == 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Unexpected state (event->iterctr != 0) in TSEventHandler()"); 821*ca4445c7SIlya Fursov if (event->processing) PetscCall(TSSetTimeStep(ts, TSEvent_dt_cap(ts, t, event->timestep_cache, PETSC_FALSE))); 822*ca4445c7SIlya Fursov event->processing = PETSC_FALSE; 823*ca4445c7SIlya Fursov } 824*ca4445c7SIlya Fursov 825*ca4445c7SIlya Fursov // if 'revisit_right' was flagged before the current iteration started, the iteration is expected to finish 826*ca4445c7SIlya Fursov PetscCheck(!revisit_right_cache || finished, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Unexpected state of 'revisit_right_cache' in TSEventHandler()"); 827*ca4445c7SIlya Fursov 828*ca4445c7SIlya Fursov if (finished) { // finished handling the current event 829*ca4445c7SIlya Fursov PetscCall(TSPostEvent(ts, t, U)); 830*ca4445c7SIlya Fursov 831*ca4445c7SIlya Fursov PetscReal dt; 832*ca4445c7SIlya Fursov PetscBool user_dt = PETSC_FALSE; 833*ca4445c7SIlya Fursov if (event->timestep_postevent > 0) { 834*ca4445c7SIlya Fursov dt = event->timestep_postevent; // user-provided post-event dt 835*ca4445c7SIlya Fursov event->processing = PETSC_FALSE; 836*ca4445c7SIlya Fursov user_dt = PETSC_TRUE; 837*ca4445c7SIlya Fursov } else { 838*ca4445c7SIlya Fursov dt = event->ptime_cache - t; // 'petsc-decide' the post-event dt 839*ca4445c7SIlya Fursov event->processing = PETSC_TRUE; 840*ca4445c7SIlya Fursov if (PetscAbsReal(dt) < PETSC_SMALL) { 841*ca4445c7SIlya Fursov dt = event->timestep_cache; // we hit the event, continue with the cached time step 842*ca4445c7SIlya Fursov event->processing = PETSC_FALSE; 843*ca4445c7SIlya Fursov } 844*ca4445c7SIlya Fursov } 845*ca4445c7SIlya Fursov PetscCall(TSSetTimeStep(ts, TSEvent_dt_cap(ts, t, dt, user_dt))); 846*ca4445c7SIlya Fursov event->iterctr = 0; 847*ca4445c7SIlya Fursov } // if-finished 848*ca4445c7SIlya Fursov 849*ca4445c7SIlya Fursov if (event->iterctr == 0) TSEvent_update_left(event, t); // not found an event, or finished the event 850*ca4445c7SIlya Fursov else { 851*ca4445c7SIlya Fursov PetscCall(TSGetTime(ts, &t)); // update 't' to account for potential rollback 852*ca4445c7SIlya Fursov PetscCall(TSSetTimeStep(ts, TSEvent_dt_cap(ts, t, dt_next, PETSC_FALSE))); // continue resolving the event 85338bf2713SShri Abhyankar } 8543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8552dc7a7e3SShri Abhyankar } 8562dc7a7e3SShri Abhyankar 857d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointEventHandler(TS ts) 858d71ae5a4SJacob Faibussowitsch { 8596427ac75SLisandro Dalcin TSEvent event; 860d0578d90SShri Abhyankar PetscReal t; 861d0578d90SShri Abhyankar Vec U; 862d0578d90SShri Abhyankar PetscInt ctr; 863*ca4445c7SIlya Fursov PetscBool forwardsolve = PETSC_FALSE; // Flag indicating that TS is doing an adjoint solve 864d0578d90SShri Abhyankar 865d0578d90SShri Abhyankar PetscFunctionBegin; 8666427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 8673ba16761SJacob Faibussowitsch if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS); 8686427ac75SLisandro Dalcin event = ts->event; 869d0578d90SShri Abhyankar 8709566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 8719566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 872d0578d90SShri Abhyankar 873d0578d90SShri Abhyankar ctr = event->recorder.ctr - 1; 874bcbf8bb3SShri Abhyankar if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 875*ca4445c7SIlya Fursov // Call the user post-event function 876d0578d90SShri Abhyankar if (event->postevent) { 877*ca4445c7SIlya Fursov PetscCallBack("TSEvent post-event processing", (*event->postevent)(ts, event->recorder.nevents[ctr], event->recorder.eventidx[ctr], t, U, forwardsolve, event->ctx)); 878d0578d90SShri Abhyankar event->recorder.ctr--; 879d0578d90SShri Abhyankar } 880d0578d90SShri Abhyankar } 8813ba16761SJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)ts)); 8823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 883d0578d90SShri Abhyankar } 8841ea83e56SMiguel 8851ea83e56SMiguel /*@ 886*ca4445c7SIlya Fursov TSGetNumEvents - Get the number of events defined on the given MPI process 8871ea83e56SMiguel 8881ea83e56SMiguel Logically Collective 8891ea83e56SMiguel 8904165533cSJose E. Roman Input Parameter: 891bcf0153eSBarry Smith . ts - the `TS` context 8921ea83e56SMiguel 8934165533cSJose E. Roman Output Parameter: 894*ca4445c7SIlya Fursov . nevents - the number of local events on each MPI process 8951ea83e56SMiguel 8961ea83e56SMiguel Level: intermediate 8971ea83e56SMiguel 8981cc06b55SBarry Smith .seealso: [](ch_ts), `TSEvent`, `TSSetEventHandler()` 8991ea83e56SMiguel @*/ 900d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetNumEvents(TS ts, PetscInt *nevents) 901d71ae5a4SJacob Faibussowitsch { 9021ea83e56SMiguel PetscFunctionBegin; 9031ea83e56SMiguel *nevents = ts->event->nevents; 9043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9051ea83e56SMiguel } 906