1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 22dc7a7e3SShri Abhyankar 36fea3669SShri Abhyankar /* 46427ac75SLisandro Dalcin TSEventInitialize - Initializes TSEvent for TSSolve 56fea3669SShri Abhyankar */ 6*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventInitialize(TSEvent event, TS ts, PetscReal t, Vec U) 7*d71ae5a4SJacob 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 19*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventDestroy(TSEvent *event) 20*d71ae5a4SJacob 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 61458122a4SShri Abhyankar Options Database Keys: 62458122a4SShri Abhyankar . -ts_event_post_eventinterval_step <dt> time-step after event interval 63458122a4SShri Abhyankar 64458122a4SShri Abhyankar Notes: 65ac6a796dSBarry Smith TSSetPostEventIntervalStep allows one to set a time-step that is used immediately following an event interval. 66458122a4SShri Abhyankar 67458122a4SShri Abhyankar This function should be called from the postevent function set with TSSetEventHandler(). 68458122a4SShri Abhyankar 69458122a4SShri Abhyankar The post event interval time-step should be selected based on the dynamics following the event. 70458122a4SShri Abhyankar If the dynamics are stiff, a conservative (small) step should be used. 71458122a4SShri Abhyankar If not, then a larger time-step can be used. 72458122a4SShri Abhyankar 73458122a4SShri Abhyankar Level: Advanced 74db781477SPatrick Sanan .seealso: `TS`, `TSEvent`, `TSSetEventHandler()` 75458122a4SShri Abhyankar @*/ 76*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt) 77*d71ae5a4SJacob Faibussowitsch { 78458122a4SShri Abhyankar PetscFunctionBegin; 79458122a4SShri Abhyankar ts->event->timestep_posteventinterval = dt; 80458122a4SShri Abhyankar PetscFunctionReturn(0); 81458122a4SShri Abhyankar } 82458122a4SShri Abhyankar 83458122a4SShri Abhyankar /*@ 84e3005195SShri Abhyankar TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 85e3005195SShri Abhyankar 86e3005195SShri Abhyankar Logically Collective 87e3005195SShri Abhyankar 884165533cSJose E. Roman Input Parameters: 89e3005195SShri Abhyankar + ts - time integration context 90e3005195SShri Abhyankar . tol - scalar tolerance, PETSC_DECIDE to leave current value 91e3005195SShri Abhyankar - vtol - array of tolerances or NULL, used in preference to tol if present 92e3005195SShri Abhyankar 932589fa24SLisandro Dalcin Options Database Keys: 94147403d9SBarry Smith . -ts_event_tol <tol> - tolerance for event zero crossing 95e3005195SShri Abhyankar 96e3005195SShri Abhyankar Notes: 97f25fe674SLisandro Dalcin Must call TSSetEventHandler() before setting the tolerances. 98e3005195SShri Abhyankar 99e3005195SShri Abhyankar The size of vtol is equal to the number of events. 100e3005195SShri Abhyankar 101e3005195SShri Abhyankar Level: beginner 102e3005195SShri Abhyankar 103db781477SPatrick Sanan .seealso: `TS`, `TSEvent`, `TSSetEventHandler()` 104e3005195SShri Abhyankar @*/ 105*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetEventTolerances(TS ts, PetscReal tol, PetscReal vtol[]) 106*d71ae5a4SJacob Faibussowitsch { 1076427ac75SLisandro Dalcin TSEvent event; 108e3005195SShri Abhyankar PetscInt i; 109e3005195SShri Abhyankar 110e3005195SShri Abhyankar PetscFunctionBegin; 1116427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1126427ac75SLisandro Dalcin if (vtol) PetscValidRealPointer(vtol, 3); 1133c633725SBarry Smith PetscCheck(ts->event, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set the events first by calling TSSetEventHandler()"); 1146427ac75SLisandro Dalcin 1156427ac75SLisandro Dalcin event = ts->event; 116e3005195SShri Abhyankar if (vtol) { 117e3005195SShri Abhyankar for (i = 0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 118e3005195SShri Abhyankar } else { 119e3005195SShri Abhyankar if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 120e3005195SShri Abhyankar for (i = 0; i < event->nevents; i++) event->vtol[i] = tol; 121e3005195SShri Abhyankar } 122e3005195SShri Abhyankar } 123e3005195SShri Abhyankar PetscFunctionReturn(0); 124e3005195SShri Abhyankar } 125e3005195SShri Abhyankar 1262dc7a7e3SShri Abhyankar /*@C 127ac6a796dSBarry Smith TSSetEventHandler - Sets a function used for detecting events 1282dc7a7e3SShri Abhyankar 1292dc7a7e3SShri Abhyankar Logically Collective on TS 1302dc7a7e3SShri Abhyankar 1312dc7a7e3SShri Abhyankar Input Parameters: 1322dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 1332dc7a7e3SShri Abhyankar . nevents - number of local events 134d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 135d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 136d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 137d94325d3SShri Abhyankar event is detected (one for each event) 1386427ac75SLisandro Dalcin . eventhandler - event monitoring routine 1392dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 1402589fa24SLisandro Dalcin - ctx - [optional] user-defined context for private data for the 1412dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 1422dc7a7e3SShri Abhyankar context is desired) 1432dc7a7e3SShri Abhyankar 1446427ac75SLisandro Dalcin Calling sequence of eventhandler: 1453a88037aSBarry Smith PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 1462dc7a7e3SShri Abhyankar 1472dc7a7e3SShri Abhyankar Input Parameters: 1482dc7a7e3SShri Abhyankar + ts - the TS context 1492dc7a7e3SShri Abhyankar . t - current time 1502dc7a7e3SShri Abhyankar . U - current iterate 1516427ac75SLisandro Dalcin - ctx - [optional] context passed with eventhandler 1522dc7a7e3SShri Abhyankar 1532dc7a7e3SShri Abhyankar Output parameters: 154d94325d3SShri Abhyankar . fvalue - function value of events at time t 1552dc7a7e3SShri Abhyankar 1562dc7a7e3SShri Abhyankar Calling sequence of postevent: 1572589fa24SLisandro Dalcin PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 1582dc7a7e3SShri Abhyankar 1592dc7a7e3SShri Abhyankar Input Parameters: 1602dc7a7e3SShri Abhyankar + ts - the TS context 1612dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 1622dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 1632dc7a7e3SShri Abhyankar . t - current time 1642dc7a7e3SShri Abhyankar . U - current solution 165031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 1666427ac75SLisandro Dalcin - ctx - the context passed with eventhandler 1672dc7a7e3SShri Abhyankar 1682dc7a7e3SShri Abhyankar Level: intermediate 1692dc7a7e3SShri Abhyankar 170db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetTimeStep()`, `TSSetConvergedReason()` 1712dc7a7e3SShri Abhyankar @*/ 172*d71ae5a4SJacob 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) 173*d71ae5a4SJacob Faibussowitsch { 174d294eb03SHong Zhang TSAdapt adapt; 175d294eb03SHong Zhang PetscReal hmin; 1762dc7a7e3SShri Abhyankar TSEvent event; 177d94325d3SShri Abhyankar PetscInt i; 178006e6a18SShri Abhyankar PetscBool flg; 179a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE 180a6c783d2SShri Abhyankar PetscReal tol = 1e-4; 181a6c783d2SShri Abhyankar #else 182d569cc17SSatish Balay PetscReal tol = 1e-6; 183a6c783d2SShri Abhyankar #endif 1842dc7a7e3SShri Abhyankar 1852dc7a7e3SShri Abhyankar PetscFunctionBegin; 1866427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1870a82b154SShri if (nevents) { 188064a246eSJacob Faibussowitsch PetscValidIntPointer(direction, 3); 189064a246eSJacob Faibussowitsch PetscValidBoolPointer(terminate, 4); 1900a82b154SShri } 1916427ac75SLisandro Dalcin 1924dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&event)); 1939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue)); 1949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue_prev)); 1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->fvalue_right)); 1969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->zerocrossing)); 1979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->side)); 1989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->direction)); 1999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->terminate)); 2009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->vtol)); 201d94325d3SShri Abhyankar for (i = 0; i < nevents; i++) { 202d94325d3SShri Abhyankar event->direction[i] = direction[i]; 203d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 204e2cdd850SShri Abhyankar event->zerocrossing[i] = PETSC_FALSE; 205e2cdd850SShri Abhyankar event->side[i] = 0; 206d94325d3SShri Abhyankar } 2079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents, &event->events_zero)); 2082589fa24SLisandro Dalcin event->nevents = nevents; 2096427ac75SLisandro Dalcin event->eventhandler = eventhandler; 2102dc7a7e3SShri Abhyankar event->postevent = postevent; 2116427ac75SLisandro Dalcin event->ctx = ctx; 212458122a4SShri Abhyankar event->timestep_posteventinterval = ts->time_step; 2139566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 2149566063dSJacob Faibussowitsch PetscCall(TSAdaptGetStepLimits(adapt, &hmin, NULL)); 215d294eb03SHong Zhang event->timestep_min = hmin; 2162dc7a7e3SShri Abhyankar 21702749585SLisandro Dalcin event->recsize = 8; /* Initial size of the recorder */ 218d0609cedSBarry Smith PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "TS Event options", "TS"); 2192dc7a7e3SShri Abhyankar { 2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_tol", "Scalar event tolerance for zero crossing check", "TSSetEventTolerances", tol, &tol, NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_event_monitor", "Print choices made by event handler", "", &flg)); 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size", "Initial size of event recorder", "", event->recsize, &event->recsize, NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_post_eventinterval_step", "Time step after event interval", "", event->timestep_posteventinterval, &event->timestep_posteventinterval, NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_post_event_step", "Time step after event", "", event->timestep_postevent, &event->timestep_postevent, NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_dt_min", "Minimum time step considered for TSEvent", "", event->timestep_min, &event->timestep_min, NULL)); 2262dc7a7e3SShri Abhyankar } 227d0609cedSBarry Smith PetscOptionsEnd(); 228a4ffd976SShri Abhyankar 2299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.time)); 2309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.stepnum)); 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.nevents)); 2329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize, &event->recorder.eventidx)); 23348a46eb9SPierre Jolivet for (i = 0; i < event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &event->recorder.eventidx[i])); 234e7069c78SShri /* Initialize the event recorder */ 235e7069c78SShri event->recorder.ctr = 0; 236a4ffd976SShri Abhyankar 237e3005195SShri Abhyankar for (i = 0; i < event->nevents; i++) event->vtol[i] = tol; 2389566063dSJacob Faibussowitsch if (flg) PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF, "stdout", &event->monitor)); 2399e12be75SShri Abhyankar 2409566063dSJacob Faibussowitsch PetscCall(TSEventDestroy(&ts->event)); 241d94325d3SShri Abhyankar ts->event = event; 242e7069c78SShri ts->event->refct = 1; 2432dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2442dc7a7e3SShri Abhyankar } 2452dc7a7e3SShri Abhyankar 246a4ffd976SShri Abhyankar /* 247a4ffd976SShri Abhyankar TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize) 248a4ffd976SShri Abhyankar is reached. 249a4ffd976SShri Abhyankar */ 250*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventRecorderResize(TSEvent event) 251*d71ae5a4SJacob Faibussowitsch { 252a4ffd976SShri Abhyankar PetscReal *time; 253a4ffd976SShri Abhyankar PetscInt *stepnum; 254a4ffd976SShri Abhyankar PetscInt *nevents; 255a4ffd976SShri Abhyankar PetscInt **eventidx; 256a4ffd976SShri Abhyankar PetscInt i, fact = 2; 257a4ffd976SShri Abhyankar 258a4ffd976SShri Abhyankar PetscFunctionBegin; 259a4ffd976SShri Abhyankar 260a4ffd976SShri Abhyankar /* Create large arrays */ 2619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &time)); 2629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &stepnum)); 2639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &nevents)); 2649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(fact * event->recsize, &eventidx)); 26548a46eb9SPierre Jolivet for (i = 0; i < fact * event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &eventidx[i])); 266a4ffd976SShri Abhyankar 267a4ffd976SShri Abhyankar /* Copy over data */ 2689566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(time, event->recorder.time, event->recsize)); 2699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(stepnum, event->recorder.stepnum, event->recsize)); 2709566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nevents, event->recorder.nevents, event->recsize)); 27148a46eb9SPierre Jolivet for (i = 0; i < event->recsize; i++) PetscCall(PetscArraycpy(eventidx[i], event->recorder.eventidx[i], event->recorder.nevents[i])); 272a4ffd976SShri Abhyankar 273a4ffd976SShri Abhyankar /* Destroy old arrays */ 27448a46eb9SPierre Jolivet for (i = 0; i < event->recsize; i++) PetscCall(PetscFree(event->recorder.eventidx[i])); 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.eventidx)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.nevents)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.stepnum)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.time)); 279a4ffd976SShri Abhyankar 280a4ffd976SShri Abhyankar /* Set pointers */ 281a4ffd976SShri Abhyankar event->recorder.time = time; 282a4ffd976SShri Abhyankar event->recorder.stepnum = stepnum; 283a4ffd976SShri Abhyankar event->recorder.nevents = nevents; 284a4ffd976SShri Abhyankar event->recorder.eventidx = eventidx; 285a4ffd976SShri Abhyankar 286a4ffd976SShri Abhyankar /* Double size */ 287a4ffd976SShri Abhyankar event->recsize *= fact; 288a4ffd976SShri Abhyankar 289a4ffd976SShri Abhyankar PetscFunctionReturn(0); 290a4ffd976SShri Abhyankar } 291a4ffd976SShri Abhyankar 292031fbad4SShri Abhyankar /* 293ac6a796dSBarry Smith Helper routine to handle user postevents and recording 294031fbad4SShri Abhyankar */ 295*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSPostEvent(TS ts, PetscReal t, Vec U) 296*d71ae5a4SJacob Faibussowitsch { 2972dc7a7e3SShri Abhyankar TSEvent event = ts->event; 2982dc7a7e3SShri Abhyankar PetscBool terminate = PETSC_FALSE; 29928d5b5d6SLisandro Dalcin PetscBool restart = PETSC_FALSE; 300d0578d90SShri Abhyankar PetscInt i, ctr, stepnum; 3017324a0ffSLisandro Dalcin PetscBool inflag[2], outflag[2]; 3024597913aSLisandro Dalcin PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 3032dc7a7e3SShri Abhyankar 3042dc7a7e3SShri Abhyankar PetscFunctionBegin; 3052dc7a7e3SShri Abhyankar if (event->postevent) { 30628d5b5d6SLisandro Dalcin PetscObjectState state_prev, state_post; 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &state_prev)); 3089566063dSJacob Faibussowitsch PetscCall((*event->postevent)(ts, event->nevents_zero, event->events_zero, t, U, forwardsolve, event->ctx)); 3099566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U, &state_post)); 31028d5b5d6SLisandro Dalcin if (state_prev != state_post) restart = PETSC_TRUE; 3112dc7a7e3SShri Abhyankar } 3124597913aSLisandro Dalcin 31328d5b5d6SLisandro Dalcin /* Handle termination events and step restart */ 3149371c9d4SSatish Balay for (i = 0; i < event->nevents_zero; i++) 3159371c9d4SSatish Balay if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE; 3169371c9d4SSatish Balay inflag[0] = restart; 3179371c9d4SSatish Balay inflag[1] = terminate; 3181c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(inflag, outflag, 2, MPIU_BOOL, MPI_LOR, ((PetscObject)ts)->comm)); 3199371c9d4SSatish Balay restart = outflag[0]; 3209371c9d4SSatish Balay terminate = outflag[1]; 3219566063dSJacob Faibussowitsch if (restart) PetscCall(TSRestartStep(ts)); 3229566063dSJacob Faibussowitsch if (terminate) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_EVENT)); 3237324a0ffSLisandro Dalcin event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP; 324f7aea88cSShri Abhyankar 3254597913aSLisandro Dalcin /* Reset event residual functions as states might get changed by the postevent callback */ 326f443add6SLisandro Dalcin if (event->postevent) { 3279566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(U)); 3289566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx)); 3299566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(U)); 330f443add6SLisandro Dalcin } 331f443add6SLisandro Dalcin 332f443add6SLisandro Dalcin /* Cache current time and event residual functions */ 333f443add6SLisandro Dalcin event->ptime_prev = t; 3349371c9d4SSatish Balay for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 3354597913aSLisandro Dalcin 336d0578d90SShri Abhyankar /* Record the event in the event recorder */ 3379566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &stepnum)); 338f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 33948a46eb9SPierre Jolivet if (ctr == event->recsize) PetscCall(TSEventRecorderResize(event)); 340f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 341d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 3424597913aSLisandro Dalcin event->recorder.nevents[ctr] = event->nevents_zero; 3434597913aSLisandro Dalcin for (i = 0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 344f7aea88cSShri Abhyankar event->recorder.ctr++; 3452dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3462dc7a7e3SShri Abhyankar } 3472dc7a7e3SShri Abhyankar 34802749585SLisandro Dalcin /* Uses Anderson-Bjorck variant of regula falsi method */ 349*d71ae5a4SJacob Faibussowitsch static inline PetscReal TSEventComputeStepSize(PetscReal tleft, PetscReal t, PetscReal tright, PetscScalar fleft, PetscScalar f, PetscScalar fright, PetscInt side, PetscReal dt) 350*d71ae5a4SJacob Faibussowitsch { 35102749585SLisandro Dalcin PetscReal new_dt, scal = 1.0; 352e2cdd850SShri Abhyankar if (PetscRealPart(fleft) * PetscRealPart(f) < 0) { 353e2cdd850SShri Abhyankar if (side == 1) { 354a6c783d2SShri Abhyankar scal = (PetscRealPart(fright) - PetscRealPart(f)) / PetscRealPart(fright); 355a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 356e2cdd850SShri Abhyankar } 35702749585SLisandro Dalcin new_dt = (scal * PetscRealPart(fleft) * t - PetscRealPart(f) * tleft) / (scal * PetscRealPart(fleft) - PetscRealPart(f)) - tleft; 358e2cdd850SShri Abhyankar } else { 359e2cdd850SShri Abhyankar if (side == -1) { 360a6c783d2SShri Abhyankar scal = (PetscRealPart(fleft) - PetscRealPart(f)) / PetscRealPart(fleft); 361a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 362e2cdd850SShri Abhyankar } 36302749585SLisandro Dalcin new_dt = (PetscRealPart(f) * tright - scal * PetscRealPart(fright) * t) / (PetscRealPart(f) - scal * PetscRealPart(fright)) - t; 364e2cdd850SShri Abhyankar } 36502749585SLisandro Dalcin return PetscMin(dt, new_dt); 36638bf2713SShri Abhyankar } 367e2cdd850SShri Abhyankar 368*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventDetection(TS ts) 369*d71ae5a4SJacob Faibussowitsch { 370d294eb03SHong Zhang TSEvent event = ts->event; 371d294eb03SHong Zhang PetscReal t; 372d294eb03SHong Zhang PetscInt i; 373d294eb03SHong Zhang PetscInt fvalue_sign, fvalueprev_sign; 374ea3dac1cSHong Zhang PetscInt in, out; 375d294eb03SHong Zhang 376d294eb03SHong Zhang PetscFunctionBegin; 3779566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 378d294eb03SHong Zhang for (i = 0; i < event->nevents; i++) { 379d294eb03SHong Zhang if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 380d294eb03SHong Zhang if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 381d294eb03SHong Zhang event->status = TSEVENT_LOCATED_INTERVAL; 382d294eb03SHong Zhang if (event->monitor) { 38363a3b9bcSJacob 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)); 384d294eb03SHong Zhang } 385d294eb03SHong Zhang continue; 386d294eb03SHong Zhang } 387ea3dac1cSHong Zhang if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */ 388d294eb03SHong Zhang fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 389d294eb03SHong Zhang fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 390d294eb03SHong Zhang if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) { 391d294eb03SHong Zhang if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 392d294eb03SHong Zhang event->status = TSEVENT_LOCATED_INTERVAL; 39348a46eb9SPierre 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)); 394d294eb03SHong Zhang } 395d294eb03SHong Zhang } 396d5f990dbSBarry Smith in = (PetscInt)event->status; 3971c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&in, &out, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts))); 398ea3dac1cSHong Zhang event->status = (TSEventStatus)out; 399d294eb03SHong Zhang PetscFunctionReturn(0); 400d294eb03SHong Zhang } 401d294eb03SHong Zhang 402*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEventLocation(TS ts, PetscReal *dt) 403*d71ae5a4SJacob Faibussowitsch { 404d294eb03SHong Zhang TSEvent event = ts->event; 405d294eb03SHong Zhang PetscInt i; 406d294eb03SHong Zhang PetscReal t; 407ea3dac1cSHong Zhang PetscInt fvalue_sign, fvalueprev_sign; 408ea3dac1cSHong Zhang PetscInt rollback = 0, in[2], out[2]; 409d294eb03SHong Zhang 410d294eb03SHong Zhang PetscFunctionBegin; 4119566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 412d294eb03SHong Zhang event->nevents_zero = 0; 413d294eb03SHong Zhang for (i = 0; i < event->nevents; i++) { 414d294eb03SHong Zhang if (event->zerocrossing[i]) { 415d294eb03SHong 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 */ 416d294eb03SHong Zhang event->status = TSEVENT_ZERO; 417d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 418d294eb03SHong Zhang continue; 419d294eb03SHong Zhang } 420d294eb03SHong Zhang /* Compute new time step */ 421d294eb03SHong 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); 422d294eb03SHong Zhang fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 423ea3dac1cSHong Zhang fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 424d294eb03SHong Zhang switch (event->direction[i]) { 425d294eb03SHong Zhang case -1: 426d294eb03SHong Zhang if (fvalue_sign < 0) { 427ea3dac1cSHong Zhang rollback = 1; 428d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 429d294eb03SHong Zhang event->side[i] = 1; 430d294eb03SHong Zhang } 431d294eb03SHong Zhang break; 432d294eb03SHong Zhang case 1: 433d294eb03SHong Zhang if (fvalue_sign > 0) { 434ea3dac1cSHong Zhang rollback = 1; 435d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 436d294eb03SHong Zhang event->side[i] = 1; 437d294eb03SHong Zhang } 438d294eb03SHong Zhang break; 439d294eb03SHong Zhang case 0: 440ea3dac1cSHong Zhang if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */ 441ea3dac1cSHong Zhang rollback = 1; 442d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 443d294eb03SHong Zhang event->side[i] = 1; 444ea3dac1cSHong Zhang } 445d294eb03SHong Zhang break; 446d294eb03SHong Zhang } 447ea3dac1cSHong Zhang if (event->status == TSEVENT_PROCESSING) event->side[i] = -1; 448d294eb03SHong Zhang } 449d294eb03SHong Zhang } 4509371c9d4SSatish Balay in[0] = (PetscInt)event->status; 4519371c9d4SSatish Balay in[1] = rollback; 4521c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts))); 4539371c9d4SSatish Balay event->status = (TSEventStatus)out[0]; 4549371c9d4SSatish Balay rollback = out[1]; 455ea3dac1cSHong 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 */ 456ea3dac1cSHong Zhang if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 457d294eb03SHong Zhang if (event->status == TSEVENT_ZERO) { 458ea3dac1cSHong Zhang for (i = 0; i < event->nevents; i++) { 459ea3dac1cSHong Zhang if (event->zerocrossing[i]) { 460ea3dac1cSHong 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 */ 461ea3dac1cSHong Zhang event->events_zero[event->nevents_zero++] = i; 46248a46eb9SPierre 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)); 463ea3dac1cSHong Zhang event->zerocrossing[i] = PETSC_FALSE; 464ea3dac1cSHong Zhang } 465ea3dac1cSHong Zhang } 466ea3dac1cSHong Zhang event->side[i] = 0; 467ea3dac1cSHong Zhang } 468d294eb03SHong Zhang } 469d294eb03SHong Zhang PetscFunctionReturn(0); 470d294eb03SHong Zhang } 47138bf2713SShri Abhyankar 472*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSEventHandler(TS ts) 473*d71ae5a4SJacob Faibussowitsch { 4746427ac75SLisandro Dalcin TSEvent event; 4752dc7a7e3SShri Abhyankar PetscReal t; 4762dc7a7e3SShri Abhyankar Vec U; 4772dc7a7e3SShri Abhyankar PetscInt i; 478ea3dac1cSHong Zhang PetscReal dt, dt_min, dt_reset = 0.0; 4792dc7a7e3SShri Abhyankar 4802dc7a7e3SShri Abhyankar PetscFunctionBegin; 4816427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 4826427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 4836427ac75SLisandro Dalcin event = ts->event; 4842dc7a7e3SShri Abhyankar 4859566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 4869566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 4879566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 4882dc7a7e3SShri Abhyankar 4897dbe0728SLisandro Dalcin if (event->status == TSEVENT_NONE) { 4907dbe0728SLisandro Dalcin event->timestep_prev = dt; 491d294eb03SHong Zhang event->ptime_end = t; 4927dbe0728SLisandro Dalcin } 4932dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 494d294eb03SHong Zhang /* user has specified a PostEventInterval dt */ 495458122a4SShri Abhyankar dt = event->timestep_posteventinterval; 496e97c63d7SStefano Zampini if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 497e97c63d7SStefano Zampini PetscReal maxdt = ts->max_time - t; 498e97c63d7SStefano Zampini dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt); 499e97c63d7SStefano Zampini } 5009566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 5012dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 5022dc7a7e3SShri Abhyankar } 5032dc7a7e3SShri Abhyankar 5049566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(U)); 5059566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts, t, U, event->fvalue, event->ctx)); 5069566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(U)); 5079e12be75SShri Abhyankar 508d294eb03SHong Zhang /* Detect the events */ 5099566063dSJacob Faibussowitsch PetscCall(TSEventDetection(ts)); 510d294eb03SHong Zhang 511d294eb03SHong Zhang /* Locate the events */ 512d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) { 513d294eb03SHong Zhang /* Approach the zero crosing by setting a new step size */ 5149566063dSJacob Faibussowitsch PetscCall(TSEventLocation(ts, &dt)); 515d294eb03SHong Zhang /* Roll back when new events are detected */ 516d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL) { 5179566063dSJacob Faibussowitsch PetscCall(TSRollBack(ts)); 5189566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); 519d294eb03SHong Zhang event->iterctr++; 520006e6a18SShri Abhyankar } 5211c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&dt, &dt_min, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 522ea3dac1cSHong Zhang if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset; 5239566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt_min)); 524d294eb03SHong Zhang /* Found the zero crossing */ 5259e12be75SShri Abhyankar if (event->status == TSEVENT_ZERO) { 5269566063dSJacob Faibussowitsch PetscCall(TSPostEvent(ts, t, U)); 5279e12be75SShri Abhyankar 528e2cdd850SShri Abhyankar dt = event->ptime_end - t; 529ad508104SStefano Zampini if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */ 530ad508104SStefano Zampini dt = event->timestep_prev; 531ad508104SStefano Zampini event->status = TSEVENT_NONE; 532ad508104SStefano Zampini } 533d294eb03SHong Zhang if (event->timestep_postevent) { /* user has specified a PostEvent dt*/ 534d294eb03SHong Zhang dt = event->timestep_postevent; 535d294eb03SHong Zhang } 536e97c63d7SStefano Zampini if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 537e97c63d7SStefano Zampini PetscReal maxdt = ts->max_time - t; 538e97c63d7SStefano Zampini dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt); 539e97c63d7SStefano Zampini } 5409566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 54138bf2713SShri Abhyankar event->iterctr = 0; 5429e12be75SShri Abhyankar } 543d294eb03SHong Zhang /* Have not found the zero crosing yet */ 544d294eb03SHong Zhang if (event->status == TSEVENT_PROCESSING) { 54548a46eb9SPierre 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)); 546d294eb03SHong Zhang event->iterctr++; 547d294eb03SHong Zhang } 548d294eb03SHong Zhang } 549d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */ 5502dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 551e2cdd850SShri Abhyankar event->ptime_right = t; 5522dc7a7e3SShri Abhyankar } else { 553d294eb03SHong Zhang for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 55438bf2713SShri Abhyankar event->ptime_prev = t; 55538bf2713SShri Abhyankar } 5562dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 5572dc7a7e3SShri Abhyankar } 5582dc7a7e3SShri Abhyankar 559*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointEventHandler(TS ts) 560*d71ae5a4SJacob Faibussowitsch { 5616427ac75SLisandro Dalcin TSEvent event; 562d0578d90SShri Abhyankar PetscReal t; 563d0578d90SShri Abhyankar Vec U; 564d0578d90SShri Abhyankar PetscInt ctr; 565d0578d90SShri Abhyankar PetscBool forwardsolve = PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 566d0578d90SShri Abhyankar 567d0578d90SShri Abhyankar PetscFunctionBegin; 5686427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5696427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 5706427ac75SLisandro Dalcin event = ts->event; 571d0578d90SShri Abhyankar 5729566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 5739566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 574d0578d90SShri Abhyankar 575d0578d90SShri Abhyankar ctr = event->recorder.ctr - 1; 576bcbf8bb3SShri Abhyankar if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 577d0578d90SShri Abhyankar /* Call the user postevent function */ 578d0578d90SShri Abhyankar if (event->postevent) { 5799566063dSJacob Faibussowitsch PetscCall((*event->postevent)(ts, event->recorder.nevents[ctr], event->recorder.eventidx[ctr], t, U, forwardsolve, event->ctx)); 580d0578d90SShri Abhyankar event->recorder.ctr--; 581d0578d90SShri Abhyankar } 582d0578d90SShri Abhyankar } 583d0578d90SShri Abhyankar 584d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 585d0578d90SShri Abhyankar PetscFunctionReturn(0); 586d0578d90SShri Abhyankar } 5871ea83e56SMiguel 5881ea83e56SMiguel /*@ 5891ea83e56SMiguel TSGetNumEvents - Get the numbers of events set 5901ea83e56SMiguel 5911ea83e56SMiguel Logically Collective 5921ea83e56SMiguel 5934165533cSJose E. Roman Input Parameter: 59499c90e12SSatish Balay . ts - the TS context 5951ea83e56SMiguel 5964165533cSJose E. Roman Output Parameter: 5971ea83e56SMiguel . nevents - number of events 5981ea83e56SMiguel 5991ea83e56SMiguel Level: intermediate 6001ea83e56SMiguel 601db781477SPatrick Sanan .seealso: `TSSetEventHandler()` 6021ea83e56SMiguel 6031ea83e56SMiguel @*/ 604*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetNumEvents(TS ts, PetscInt *nevents) 605*d71ae5a4SJacob Faibussowitsch { 6061ea83e56SMiguel PetscFunctionBegin; 6071ea83e56SMiguel *nevents = ts->event->nevents; 6081ea83e56SMiguel PetscFunctionReturn(0); 6091ea83e56SMiguel } 610