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; 96427ac75SLisandro Dalcin if (!event) PetscFunctionReturn(0); 106427ac75SLisandro Dalcin PetscValidPointer(event, 1); 116427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 126427ac75SLisandro Dalcin PetscValidHeaderSpecific(U, VEC_CLASSID, 4); 136fea3669SShri Abhyankar event->ptime_prev = t; 1438bf2713SShri Abhyankar event->iterctr = 0; 159566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts, t, U, event->fvalue_prev, event->ctx)); 166fea3669SShri Abhyankar PetscFunctionReturn(0); 176fea3669SShri Abhyankar } 186fea3669SShri Abhyankar 19d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventDestroy(TSEvent *event) 20d71ae5a4SJacob Faibussowitsch { 217dbe0728SLisandro Dalcin PetscInt i; 227dbe0728SLisandro Dalcin 237dbe0728SLisandro Dalcin PetscFunctionBegin; 247dbe0728SLisandro Dalcin PetscValidPointer(event, 1); 257dbe0728SLisandro Dalcin if (!*event) PetscFunctionReturn(0); 269371c9d4SSatish Balay if (--(*event)->refct > 0) { 279371c9d4SSatish Balay *event = NULL; 289371c9d4SSatish Balay PetscFunctionReturn(0); 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)); 497dbe0728SLisandro Dalcin PetscFunctionReturn(0); 507dbe0728SLisandro Dalcin } 517dbe0728SLisandro Dalcin 52e3005195SShri Abhyankar /*@ 53ac6a796dSBarry Smith TSSetPostEventIntervalStep - Set the time-step used immediately following the 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 61*bcf0153eSBarry Smith Options Database Key: 62458122a4SShri Abhyankar . -ts_event_post_eventinterval_step <dt> time-step after event interval 63458122a4SShri Abhyankar 64*bcf0153eSBarry Smith Level: advanced 65458122a4SShri Abhyankar 66*bcf0153eSBarry Smith Notes: 67*bcf0153eSBarry Smith `TSSetPostEventIntervalStep()` allows one to set a time-step that is used immediately following an event interval. 68*bcf0153eSBarry Smith 69*bcf0153eSBarry 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 75*bcf0153eSBarry 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; 81458122a4SShri Abhyankar PetscFunctionReturn(0); 82458122a4SShri Abhyankar } 83458122a4SShri Abhyankar 84458122a4SShri Abhyankar /*@ 85e3005195SShri Abhyankar TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 86e3005195SShri Abhyankar 87e3005195SShri Abhyankar Logically Collective 88e3005195SShri Abhyankar 894165533cSJose E. Roman Input Parameters: 90e3005195SShri Abhyankar + ts - time integration context 91*bcf0153eSBarry Smith . tol - scalar tolerance, `PETSC_DECIDE` to leave current value 92e3005195SShri Abhyankar - vtol - array of tolerances or NULL, used in preference to tol if present 93e3005195SShri Abhyankar 94*bcf0153eSBarry Smith Options Database Key: 95147403d9SBarry Smith . -ts_event_tol <tol> - tolerance for event zero crossing 96e3005195SShri Abhyankar 97e3005195SShri Abhyankar Level: beginner 98e3005195SShri Abhyankar 99*bcf0153eSBarry Smith Notes: 100*bcf0153eSBarry Smith Must call `TSSetEventHandler(`) before setting the tolerances. 101*bcf0153eSBarry Smith 102*bcf0153eSBarry Smith The size of vtol is equal to the number of events. 103*bcf0153eSBarry Smith 104*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSEvent`, `TSSetEventHandler()` 105e3005195SShri Abhyankar @*/ 106d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetEventTolerances(TS ts, PetscReal tol, PetscReal vtol[]) 107d71ae5a4SJacob Faibussowitsch { 1086427ac75SLisandro Dalcin TSEvent event; 109e3005195SShri Abhyankar PetscInt i; 110e3005195SShri Abhyankar 111e3005195SShri Abhyankar PetscFunctionBegin; 1126427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1136427ac75SLisandro Dalcin if (vtol) PetscValidRealPointer(vtol, 3); 1143c633725SBarry Smith PetscCheck(ts->event, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set the events first by calling TSSetEventHandler()"); 1156427ac75SLisandro Dalcin 1166427ac75SLisandro Dalcin event = ts->event; 117e3005195SShri Abhyankar if (vtol) { 118e3005195SShri Abhyankar for (i = 0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 119e3005195SShri Abhyankar } else { 120e3005195SShri Abhyankar if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 121e3005195SShri Abhyankar for (i = 0; i < event->nevents; i++) event->vtol[i] = tol; 122e3005195SShri Abhyankar } 123e3005195SShri Abhyankar } 124e3005195SShri Abhyankar PetscFunctionReturn(0); 125e3005195SShri Abhyankar } 126e3005195SShri Abhyankar 1272dc7a7e3SShri Abhyankar /*@C 128ac6a796dSBarry Smith TSSetEventHandler - Sets a function used for detecting events 1292dc7a7e3SShri Abhyankar 130*bcf0153eSBarry Smith Logically Collective on ts 1312dc7a7e3SShri Abhyankar 1322dc7a7e3SShri Abhyankar Input Parameters: 133*bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1342dc7a7e3SShri Abhyankar . nevents - number of local events 135d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 136d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 137d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 138d94325d3SShri Abhyankar event is detected (one for each event) 1396427ac75SLisandro Dalcin . eventhandler - event monitoring routine 1402dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 1412589fa24SLisandro Dalcin - ctx - [optional] user-defined context for private data for the 1422dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 1432dc7a7e3SShri Abhyankar context is desired) 1442dc7a7e3SShri Abhyankar 1456427ac75SLisandro Dalcin Calling sequence of eventhandler: 1463a88037aSBarry Smith PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 1472dc7a7e3SShri Abhyankar 1482dc7a7e3SShri Abhyankar Input Parameters: 1492dc7a7e3SShri Abhyankar + ts - the TS context 1502dc7a7e3SShri Abhyankar . t - current time 1512dc7a7e3SShri Abhyankar . U - current iterate 1526427ac75SLisandro Dalcin - ctx - [optional] context passed with eventhandler 1532dc7a7e3SShri Abhyankar 1542dc7a7e3SShri Abhyankar Output parameters: 155d94325d3SShri Abhyankar . fvalue - function value of events at time t 1562dc7a7e3SShri Abhyankar 1572dc7a7e3SShri Abhyankar Calling sequence of postevent: 1582589fa24SLisandro Dalcin PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 1592dc7a7e3SShri Abhyankar 1602dc7a7e3SShri Abhyankar Input Parameters: 1612dc7a7e3SShri Abhyankar + ts - the TS context 1622dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 1632dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 1642dc7a7e3SShri Abhyankar . t - current time 1652dc7a7e3SShri Abhyankar . U - current solution 166031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 1676427ac75SLisandro Dalcin - ctx - the context passed with eventhandler 1682dc7a7e3SShri Abhyankar 1692dc7a7e3SShri Abhyankar Level: intermediate 1702dc7a7e3SShri Abhyankar 171*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSEvent`, `TSCreate()`, `TSSetTimeStep()`, `TSSetConvergedReason()` 1722dc7a7e3SShri Abhyankar @*/ 173d71ae5a4SJacob 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) 174d71ae5a4SJacob Faibussowitsch { 175d294eb03SHong Zhang TSAdapt adapt; 176d294eb03SHong Zhang PetscReal hmin; 1772dc7a7e3SShri Abhyankar TSEvent event; 178d94325d3SShri Abhyankar PetscInt i; 179006e6a18SShri Abhyankar PetscBool flg; 180a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE 181a6c783d2SShri Abhyankar PetscReal tol = 1e-4; 182a6c783d2SShri Abhyankar #else 183d569cc17SSatish Balay PetscReal tol = 1e-6; 184a6c783d2SShri Abhyankar #endif 1852dc7a7e3SShri Abhyankar 1862dc7a7e3SShri Abhyankar PetscFunctionBegin; 1876427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1880a82b154SShri if (nevents) { 189064a246eSJacob Faibussowitsch PetscValidIntPointer(direction, 3); 190064a246eSJacob Faibussowitsch PetscValidBoolPointer(terminate, 4); 1910a82b154SShri } 1926427ac75SLisandro Dalcin 1934dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&event)); 1949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue)); 1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue_prev)); 1969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue_right)); 1979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->zerocrossing)); 1989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->side)); 1999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->direction)); 2009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->terminate)); 2019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->vtol)); 202d94325d3SShri Abhyankar for (i = 0; i < nevents; i++) { 203d94325d3SShri Abhyankar event->direction[i] = direction[i]; 204d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 205e2cdd850SShri Abhyankar event->zerocrossing[i] = PETSC_FALSE; 206e2cdd850SShri Abhyankar event->side[i] = 0; 207d94325d3SShri Abhyankar } 2089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->events_zero)); 2092589fa24SLisandro Dalcin event->nevents = nevents; 2106427ac75SLisandro Dalcin event->eventhandler = eventhandler; 2112dc7a7e3SShri Abhyankar event->postevent = postevent; 2126427ac75SLisandro Dalcin event->ctx = ctx; 213458122a4SShri Abhyankar event->timestep_posteventinterval = ts->time_step; 2149566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 2159566063dSJacob Faibussowitsch PetscCall(TSAdaptGetStepLimits(adapt, &hmin, NULL)); 216d294eb03SHong Zhang event->timestep_min = hmin; 2172dc7a7e3SShri Abhyankar 21802749585SLisandro Dalcin event->recsize = 8; /* Initial size of the recorder */ 219d0609cedSBarry Smith PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "TS Event options", "TS"); 2202dc7a7e3SShri Abhyankar { 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_tol", "Scalar event tolerance for zero crossing check", "TSSetEventTolerances", tol, &tol, NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_event_monitor", "Print choices made by event handler", "", &flg)); 2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size", "Initial size of event recorder", "", event->recsize, &event->recsize, NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_post_eventinterval_step", "Time step after event interval", "", event->timestep_posteventinterval, &event->timestep_posteventinterval, NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_post_event_step", "Time step after event", "", event->timestep_postevent, &event->timestep_postevent, NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_dt_min", "Minimum time step considered for TSEvent", "", event->timestep_min, &event->timestep_min, NULL)); 2272dc7a7e3SShri Abhyankar } 228d0609cedSBarry Smith PetscOptionsEnd(); 229a4ffd976SShri Abhyankar 2309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.time)); 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.stepnum)); 2329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.nevents)); 2339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.eventidx)); 23448a46eb9SPierre Jolivet for (i = 0; i < event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &event->recorder.eventidx[i])); 235e7069c78SShri /* Initialize the event recorder */ 236e7069c78SShri event->recorder.ctr = 0; 237a4ffd976SShri Abhyankar 238e3005195SShri Abhyankar for (i = 0; i < event->nevents; i++) event->vtol[i] = tol; 2399566063dSJacob Faibussowitsch if (flg) PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF, "stdout", &event->monitor)); 2409e12be75SShri Abhyankar 2419566063dSJacob Faibussowitsch PetscCall(TSEventDestroy(&ts->event)); 242d94325d3SShri Abhyankar ts->event = event; 243e7069c78SShri ts->event->refct = 1; 2442dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2452dc7a7e3SShri Abhyankar } 2462dc7a7e3SShri Abhyankar 247a4ffd976SShri Abhyankar /* 248a4ffd976SShri Abhyankar TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize) 249a4ffd976SShri Abhyankar is reached. 250a4ffd976SShri Abhyankar */ 251d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventRecorderResize(TSEvent event) 252d71ae5a4SJacob Faibussowitsch { 253a4ffd976SShri Abhyankar PetscReal *time; 254a4ffd976SShri Abhyankar PetscInt *stepnum; 255a4ffd976SShri Abhyankar PetscInt *nevents; 256a4ffd976SShri Abhyankar PetscInt **eventidx; 257a4ffd976SShri Abhyankar PetscInt i, fact = 2; 258a4ffd976SShri Abhyankar 259a4ffd976SShri Abhyankar PetscFunctionBegin; 260a4ffd976SShri Abhyankar 261a4ffd976SShri Abhyankar /* Create large arrays */ 2629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &time)); 2639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &stepnum)); 2649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &nevents)); 2659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &eventidx)); 26648a46eb9SPierre Jolivet for (i = 0; i < fact * event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &eventidx[i])); 267a4ffd976SShri Abhyankar 268a4ffd976SShri Abhyankar /* Copy over data */ 2699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(time, event->recorder.time, event->recsize)); 2709566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(stepnum, event->recorder.stepnum, event->recsize)); 2719566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nevents, event->recorder.nevents, event->recsize)); 27248a46eb9SPierre Jolivet for (i = 0; i < event->recsize; i++) PetscCall(PetscArraycpy(eventidx[i], event->recorder.eventidx[i], event->recorder.nevents[i])); 273a4ffd976SShri Abhyankar 274a4ffd976SShri Abhyankar /* Destroy old arrays */ 27548a46eb9SPierre Jolivet for (i = 0; i < event->recsize; i++) PetscCall(PetscFree(event->recorder.eventidx[i])); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.eventidx)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.nevents)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.stepnum)); 2799566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.time)); 280a4ffd976SShri Abhyankar 281a4ffd976SShri Abhyankar /* Set pointers */ 282a4ffd976SShri Abhyankar event->recorder.time = time; 283a4ffd976SShri Abhyankar event->recorder.stepnum = stepnum; 284a4ffd976SShri Abhyankar event->recorder.nevents = nevents; 285a4ffd976SShri Abhyankar event->recorder.eventidx = eventidx; 286a4ffd976SShri Abhyankar 287a4ffd976SShri Abhyankar /* Double size */ 288a4ffd976SShri Abhyankar event->recsize *= fact; 289a4ffd976SShri Abhyankar 290a4ffd976SShri Abhyankar PetscFunctionReturn(0); 291a4ffd976SShri Abhyankar } 292a4ffd976SShri Abhyankar 293031fbad4SShri Abhyankar /* 294ac6a796dSBarry Smith Helper routine to handle user postevents and recording 295031fbad4SShri Abhyankar */ 296d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPostEvent(TS ts, PetscReal t, Vec U) 297d71ae5a4SJacob Faibussowitsch { 2982dc7a7e3SShri Abhyankar TSEvent event = ts->event; 2992dc7a7e3SShri Abhyankar PetscBool terminate = PETSC_FALSE; 30028d5b5d6SLisandro Dalcin PetscBool restart = PETSC_FALSE; 301d0578d90SShri Abhyankar PetscInt i, ctr, stepnum; 3027324a0ffSLisandro Dalcin PetscBool inflag[2], outflag[2]; 3034597913aSLisandro Dalcin PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 3042dc7a7e3SShri Abhyankar 3052dc7a7e3SShri Abhyankar PetscFunctionBegin; 3062dc7a7e3SShri Abhyankar if (event->postevent) { 30728d5b5d6SLisandro Dalcin PetscObjectState state_prev, state_post; 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &state_prev)); 3099566063dSJacob Faibussowitsch PetscCall((*event->postevent)(ts, event->nevents_zero, event->events_zero, t, U, forwardsolve, event->ctx)); 3109566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &state_post)); 31128d5b5d6SLisandro Dalcin if (state_prev != state_post) restart = PETSC_TRUE; 3122dc7a7e3SShri Abhyankar } 3134597913aSLisandro Dalcin 31428d5b5d6SLisandro Dalcin /* Handle termination events and step restart */ 3159371c9d4SSatish Balay for (i = 0; i < event->nevents_zero; i++) 3169371c9d4SSatish Balay if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE; 3179371c9d4SSatish Balay inflag[0] = restart; 3189371c9d4SSatish Balay inflag[1] = terminate; 3191c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(inflag, outflag, 2, MPIU_BOOL, MPI_LOR, ((PetscObject)ts)->comm)); 3209371c9d4SSatish Balay restart = outflag[0]; 3219371c9d4SSatish Balay terminate = outflag[1]; 3229566063dSJacob Faibussowitsch if (restart) PetscCall(TSRestartStep(ts)); 3239566063dSJacob Faibussowitsch if (terminate) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_EVENT)); 3247324a0ffSLisandro Dalcin event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP; 325f7aea88cSShri Abhyankar 3264597913aSLisandro Dalcin /* Reset event residual functions as states might get changed by the postevent callback */ 327f443add6SLisandro Dalcin if (event->postevent) { 3289566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(U)); 3299566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx)); 3309566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(U)); 331f443add6SLisandro Dalcin } 332f443add6SLisandro Dalcin 333f443add6SLisandro Dalcin /* Cache current time and event residual functions */ 334f443add6SLisandro Dalcin event->ptime_prev = t; 3359371c9d4SSatish Balay for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 3364597913aSLisandro Dalcin 337d0578d90SShri Abhyankar /* Record the event in the event recorder */ 3389566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &stepnum)); 339f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 34048a46eb9SPierre Jolivet if (ctr == event->recsize) PetscCall(TSEventRecorderResize(event)); 341f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 342d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 3434597913aSLisandro Dalcin event->recorder.nevents[ctr] = event->nevents_zero; 3444597913aSLisandro Dalcin for (i = 0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 345f7aea88cSShri Abhyankar event->recorder.ctr++; 3462dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3472dc7a7e3SShri Abhyankar } 3482dc7a7e3SShri Abhyankar 34902749585SLisandro Dalcin /* Uses Anderson-Bjorck variant of regula falsi method */ 350d71ae5a4SJacob Faibussowitsch static inline PetscReal TSEventComputeStepSize(PetscReal tleft, PetscReal t, PetscReal tright, PetscScalar fleft, PetscScalar f, PetscScalar fright, PetscInt side, PetscReal dt) 351d71ae5a4SJacob Faibussowitsch { 35202749585SLisandro Dalcin PetscReal new_dt, scal = 1.0; 353e2cdd850SShri Abhyankar if (PetscRealPart(fleft) * PetscRealPart(f) < 0) { 354e2cdd850SShri Abhyankar if (side == 1) { 355a6c783d2SShri Abhyankar scal = (PetscRealPart(fright) - PetscRealPart(f)) / PetscRealPart(fright); 356a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 357e2cdd850SShri Abhyankar } 35802749585SLisandro Dalcin new_dt = (scal * PetscRealPart(fleft) * t - PetscRealPart(f) * tleft) / (scal * PetscRealPart(fleft) - PetscRealPart(f)) - tleft; 359e2cdd850SShri Abhyankar } else { 360e2cdd850SShri Abhyankar if (side == -1) { 361a6c783d2SShri Abhyankar scal = (PetscRealPart(fleft) - PetscRealPart(f)) / PetscRealPart(fleft); 362a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 363e2cdd850SShri Abhyankar } 36402749585SLisandro Dalcin new_dt = (PetscRealPart(f) * tright - scal * PetscRealPart(fright) * t) / (PetscRealPart(f) - scal * PetscRealPart(fright)) - t; 365e2cdd850SShri Abhyankar } 36602749585SLisandro Dalcin return PetscMin(dt, new_dt); 36738bf2713SShri Abhyankar } 368e2cdd850SShri Abhyankar 369d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventDetection(TS ts) 370d71ae5a4SJacob Faibussowitsch { 371d294eb03SHong Zhang TSEvent event = ts->event; 372d294eb03SHong Zhang PetscReal t; 373d294eb03SHong Zhang PetscInt i; 374d294eb03SHong Zhang PetscInt fvalue_sign, fvalueprev_sign; 375ea3dac1cSHong Zhang PetscInt in, out; 376d294eb03SHong Zhang 377d294eb03SHong Zhang PetscFunctionBegin; 3789566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 379d294eb03SHong Zhang for (i = 0; i < event->nevents; i++) { 380d294eb03SHong Zhang if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 381d294eb03SHong Zhang if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 382d294eb03SHong Zhang event->status = TSEVENT_LOCATED_INTERVAL; 383d294eb03SHong Zhang if (event->monitor) { 38463a3b9bcSJacob 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)); 385d294eb03SHong Zhang } 386d294eb03SHong Zhang continue; 387d294eb03SHong Zhang } 388ea3dac1cSHong Zhang if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */ 389d294eb03SHong Zhang fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 390d294eb03SHong Zhang fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 391d294eb03SHong Zhang if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) { 392d294eb03SHong Zhang if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 393d294eb03SHong Zhang event->status = TSEVENT_LOCATED_INTERVAL; 39448a46eb9SPierre 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)); 395d294eb03SHong Zhang } 396d294eb03SHong Zhang } 397d5f990dbSBarry Smith in = (PetscInt)event->status; 3981c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&in, &out, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts))); 399ea3dac1cSHong Zhang event->status = (TSEventStatus)out; 400d294eb03SHong Zhang PetscFunctionReturn(0); 401d294eb03SHong Zhang } 402d294eb03SHong Zhang 403d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventLocation(TS ts, PetscReal *dt) 404d71ae5a4SJacob Faibussowitsch { 405d294eb03SHong Zhang TSEvent event = ts->event; 406d294eb03SHong Zhang PetscInt i; 407d294eb03SHong Zhang PetscReal t; 408ea3dac1cSHong Zhang PetscInt fvalue_sign, fvalueprev_sign; 409ea3dac1cSHong Zhang PetscInt rollback = 0, in[2], out[2]; 410d294eb03SHong Zhang 411d294eb03SHong Zhang PetscFunctionBegin; 4129566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 413d294eb03SHong Zhang event->nevents_zero = 0; 414d294eb03SHong Zhang for (i = 0; i < event->nevents; i++) { 415d294eb03SHong Zhang if (event->zerocrossing[i]) { 416d294eb03SHong 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 */ 417d294eb03SHong Zhang event->status = TSEVENT_ZERO; 418d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 419d294eb03SHong Zhang continue; 420d294eb03SHong Zhang } 421d294eb03SHong Zhang /* Compute new time step */ 422d294eb03SHong 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); 423d294eb03SHong Zhang fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 424ea3dac1cSHong Zhang fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 425d294eb03SHong Zhang switch (event->direction[i]) { 426d294eb03SHong Zhang case -1: 427d294eb03SHong Zhang if (fvalue_sign < 0) { 428ea3dac1cSHong Zhang rollback = 1; 429d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 430d294eb03SHong Zhang event->side[i] = 1; 431d294eb03SHong Zhang } 432d294eb03SHong Zhang break; 433d294eb03SHong Zhang case 1: 434d294eb03SHong Zhang if (fvalue_sign > 0) { 435ea3dac1cSHong Zhang rollback = 1; 436d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 437d294eb03SHong Zhang event->side[i] = 1; 438d294eb03SHong Zhang } 439d294eb03SHong Zhang break; 440d294eb03SHong Zhang case 0: 441ea3dac1cSHong Zhang if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */ 442ea3dac1cSHong Zhang rollback = 1; 443d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 444d294eb03SHong Zhang event->side[i] = 1; 445ea3dac1cSHong Zhang } 446d294eb03SHong Zhang break; 447d294eb03SHong Zhang } 448ea3dac1cSHong Zhang if (event->status == TSEVENT_PROCESSING) event->side[i] = -1; 449d294eb03SHong Zhang } 450d294eb03SHong Zhang } 4519371c9d4SSatish Balay in[0] = (PetscInt)event->status; 4529371c9d4SSatish Balay in[1] = rollback; 4531c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts))); 4549371c9d4SSatish Balay event->status = (TSEventStatus)out[0]; 4559371c9d4SSatish Balay rollback = out[1]; 456ea3dac1cSHong Zhang /* If rollback is true, the status will be overwritten so that an event at the endtime of current time step will be postponed to guarantee corret order */ 457ea3dac1cSHong Zhang if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 458d294eb03SHong Zhang if (event->status == TSEVENT_ZERO) { 459ea3dac1cSHong Zhang for (i = 0; i < event->nevents; i++) { 460ea3dac1cSHong Zhang if (event->zerocrossing[i]) { 461ea3dac1cSHong 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 */ 462ea3dac1cSHong Zhang event->events_zero[event->nevents_zero++] = i; 46348a46eb9SPierre 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)); 464ea3dac1cSHong Zhang event->zerocrossing[i] = PETSC_FALSE; 465ea3dac1cSHong Zhang } 466ea3dac1cSHong Zhang } 467ea3dac1cSHong Zhang event->side[i] = 0; 468ea3dac1cSHong Zhang } 469d294eb03SHong Zhang } 470d294eb03SHong Zhang PetscFunctionReturn(0); 471d294eb03SHong Zhang } 47238bf2713SShri Abhyankar 473d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventHandler(TS ts) 474d71ae5a4SJacob Faibussowitsch { 4756427ac75SLisandro Dalcin TSEvent event; 4762dc7a7e3SShri Abhyankar PetscReal t; 4772dc7a7e3SShri Abhyankar Vec U; 4782dc7a7e3SShri Abhyankar PetscInt i; 479ea3dac1cSHong Zhang PetscReal dt, dt_min, dt_reset = 0.0; 4802dc7a7e3SShri Abhyankar 4812dc7a7e3SShri Abhyankar PetscFunctionBegin; 4826427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 4836427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 4846427ac75SLisandro Dalcin event = ts->event; 4852dc7a7e3SShri Abhyankar 4869566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 4879566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 4889566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 4892dc7a7e3SShri Abhyankar 4907dbe0728SLisandro Dalcin if (event->status == TSEVENT_NONE) { 4917dbe0728SLisandro Dalcin event->timestep_prev = dt; 492d294eb03SHong Zhang event->ptime_end = t; 4937dbe0728SLisandro Dalcin } 4942dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 495d294eb03SHong Zhang /* user has specified a PostEventInterval dt */ 496458122a4SShri Abhyankar dt = event->timestep_posteventinterval; 497e97c63d7SStefano Zampini if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 498e97c63d7SStefano Zampini PetscReal maxdt = ts->max_time - t; 499e97c63d7SStefano Zampini dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt); 500e97c63d7SStefano Zampini } 5019566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 5022dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 5032dc7a7e3SShri Abhyankar } 5042dc7a7e3SShri Abhyankar 5059566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(U)); 5069566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx)); 5079566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(U)); 5089e12be75SShri Abhyankar 509d294eb03SHong Zhang /* Detect the events */ 5109566063dSJacob Faibussowitsch PetscCall(TSEventDetection(ts)); 511d294eb03SHong Zhang 512d294eb03SHong Zhang /* Locate the events */ 513d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) { 514d294eb03SHong Zhang /* Approach the zero crosing by setting a new step size */ 5159566063dSJacob Faibussowitsch PetscCall(TSEventLocation(ts, &dt)); 516d294eb03SHong Zhang /* Roll back when new events are detected */ 517d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL) { 5189566063dSJacob Faibussowitsch PetscCall(TSRollBack(ts)); 5199566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); 520d294eb03SHong Zhang event->iterctr++; 521006e6a18SShri Abhyankar } 5221c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&dt, &dt_min, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 523ea3dac1cSHong Zhang if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset; 5249566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt_min)); 525d294eb03SHong Zhang /* Found the zero crossing */ 5269e12be75SShri Abhyankar if (event->status == TSEVENT_ZERO) { 5279566063dSJacob Faibussowitsch PetscCall(TSPostEvent(ts, t, U)); 5289e12be75SShri Abhyankar 529e2cdd850SShri Abhyankar dt = event->ptime_end - t; 530ad508104SStefano Zampini if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */ 531ad508104SStefano Zampini dt = event->timestep_prev; 532ad508104SStefano Zampini event->status = TSEVENT_NONE; 533ad508104SStefano Zampini } 534d294eb03SHong Zhang if (event->timestep_postevent) { /* user has specified a PostEvent dt*/ 535d294eb03SHong Zhang dt = event->timestep_postevent; 536d294eb03SHong Zhang } 537e97c63d7SStefano Zampini if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 538e97c63d7SStefano Zampini PetscReal maxdt = ts->max_time - t; 539e97c63d7SStefano Zampini dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt); 540e97c63d7SStefano Zampini } 5419566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 54238bf2713SShri Abhyankar event->iterctr = 0; 5439e12be75SShri Abhyankar } 544d294eb03SHong Zhang /* Have not found the zero crosing yet */ 545d294eb03SHong Zhang if (event->status == TSEVENT_PROCESSING) { 54648a46eb9SPierre 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)); 547d294eb03SHong Zhang event->iterctr++; 548d294eb03SHong Zhang } 549d294eb03SHong Zhang } 550d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */ 5512dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 552e2cdd850SShri Abhyankar event->ptime_right = t; 5532dc7a7e3SShri Abhyankar } else { 554d294eb03SHong Zhang for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 55538bf2713SShri Abhyankar event->ptime_prev = t; 55638bf2713SShri Abhyankar } 5572dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 5582dc7a7e3SShri Abhyankar } 5592dc7a7e3SShri Abhyankar 560d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointEventHandler(TS ts) 561d71ae5a4SJacob Faibussowitsch { 5626427ac75SLisandro Dalcin TSEvent event; 563d0578d90SShri Abhyankar PetscReal t; 564d0578d90SShri Abhyankar Vec U; 565d0578d90SShri Abhyankar PetscInt ctr; 566d0578d90SShri Abhyankar PetscBool forwardsolve = PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 567d0578d90SShri Abhyankar 568d0578d90SShri Abhyankar PetscFunctionBegin; 5696427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5706427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 5716427ac75SLisandro Dalcin event = ts->event; 572d0578d90SShri Abhyankar 5739566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 5749566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 575d0578d90SShri Abhyankar 576d0578d90SShri Abhyankar ctr = event->recorder.ctr - 1; 577bcbf8bb3SShri Abhyankar if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 578d0578d90SShri Abhyankar /* Call the user postevent function */ 579d0578d90SShri Abhyankar if (event->postevent) { 5809566063dSJacob Faibussowitsch PetscCall((*event->postevent)(ts, event->recorder.nevents[ctr], event->recorder.eventidx[ctr], t, U, forwardsolve, event->ctx)); 581d0578d90SShri Abhyankar event->recorder.ctr--; 582d0578d90SShri Abhyankar } 583d0578d90SShri Abhyankar } 584d0578d90SShri Abhyankar 585d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 586d0578d90SShri Abhyankar PetscFunctionReturn(0); 587d0578d90SShri Abhyankar } 5881ea83e56SMiguel 5891ea83e56SMiguel /*@ 5901ea83e56SMiguel TSGetNumEvents - Get the numbers of events set 5911ea83e56SMiguel 5921ea83e56SMiguel Logically Collective 5931ea83e56SMiguel 5944165533cSJose E. Roman Input Parameter: 595*bcf0153eSBarry Smith . ts - the `TS` context 5961ea83e56SMiguel 5974165533cSJose E. Roman Output Parameter: 5981ea83e56SMiguel . nevents - number of events 5991ea83e56SMiguel 6001ea83e56SMiguel Level: intermediate 6011ea83e56SMiguel 602*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSEvent`, `TSSetEventHandler()` 6031ea83e56SMiguel @*/ 604d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetNumEvents(TS ts, PetscInt *nevents) 605d71ae5a4SJacob Faibussowitsch { 6061ea83e56SMiguel PetscFunctionBegin; 6071ea83e56SMiguel *nevents = ts->event->nevents; 6081ea83e56SMiguel PetscFunctionReturn(0); 6091ea83e56SMiguel } 610