1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 22dc7a7e3SShri Abhyankar 36fea3669SShri Abhyankar /* 46427ac75SLisandro Dalcin TSEventInitialize - Initializes TSEvent for TSSolve 56fea3669SShri Abhyankar */ 66427ac75SLisandro Dalcin PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U) 76fea3669SShri Abhyankar { 86fea3669SShri Abhyankar PetscFunctionBegin; 96427ac75SLisandro Dalcin if (!event) PetscFunctionReturn(0); 106427ac75SLisandro Dalcin PetscValidPointer(event,1); 116427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 126427ac75SLisandro Dalcin PetscValidHeaderSpecific(U,VEC_CLASSID,4); 136fea3669SShri Abhyankar event->ptime_prev = t; 1438bf2713SShri Abhyankar event->iterctr = 0; 159566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx)); 166fea3669SShri Abhyankar PetscFunctionReturn(0); 176fea3669SShri Abhyankar } 186fea3669SShri Abhyankar 197dbe0728SLisandro Dalcin PetscErrorCode TSEventDestroy(TSEvent *event) 207dbe0728SLisandro Dalcin { 217dbe0728SLisandro Dalcin PetscInt i; 227dbe0728SLisandro Dalcin 237dbe0728SLisandro Dalcin PetscFunctionBegin; 247dbe0728SLisandro Dalcin PetscValidPointer(event,1); 257dbe0728SLisandro Dalcin if (!*event) PetscFunctionReturn(0); 26c793f718SLisandro Dalcin if (--(*event)->refct > 0) {*event = NULL; PetscFunctionReturn(0);} 27e7069c78SShri 289566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->fvalue)); 299566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->fvalue_prev)); 309566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->fvalue_right)); 319566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->zerocrossing)); 329566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->side)); 339566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->direction)); 349566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->terminate)); 359566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->events_zero)); 369566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->vtol)); 37a4ffd976SShri Abhyankar 38a4ffd976SShri Abhyankar for (i=0; i < (*event)->recsize; i++) { 399566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.eventidx[i])); 407dbe0728SLisandro Dalcin } 419566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.eventidx)); 429566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.nevents)); 439566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.stepnum)); 449566063dSJacob Faibussowitsch PetscCall(PetscFree((*event)->recorder.time)); 45a4ffd976SShri Abhyankar 469566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*event)->monitor)); 479566063dSJacob Faibussowitsch PetscCall(PetscFree(*event)); 487dbe0728SLisandro Dalcin PetscFunctionReturn(0); 497dbe0728SLisandro Dalcin } 507dbe0728SLisandro Dalcin 51e3005195SShri Abhyankar /*@ 52ac6a796dSBarry Smith TSSetPostEventIntervalStep - Set the time-step used immediately following the event interval 53458122a4SShri Abhyankar 54458122a4SShri Abhyankar Logically Collective 55458122a4SShri Abhyankar 564165533cSJose E. Roman Input Parameters: 57458122a4SShri Abhyankar + ts - time integration context 58458122a4SShri Abhyankar - dt - post event interval step 59458122a4SShri Abhyankar 60458122a4SShri Abhyankar Options Database Keys: 61458122a4SShri Abhyankar . -ts_event_post_eventinterval_step <dt> time-step after event interval 62458122a4SShri Abhyankar 63458122a4SShri Abhyankar Notes: 64ac6a796dSBarry Smith TSSetPostEventIntervalStep allows one to set a time-step that is used immediately following an event interval. 65458122a4SShri Abhyankar 66458122a4SShri Abhyankar This function should be called from the postevent function set with TSSetEventHandler(). 67458122a4SShri Abhyankar 68458122a4SShri Abhyankar The post event interval time-step should be selected based on the dynamics following the event. 69458122a4SShri Abhyankar If the dynamics are stiff, a conservative (small) step should be used. 70458122a4SShri Abhyankar If not, then a larger time-step can be used. 71458122a4SShri Abhyankar 72458122a4SShri Abhyankar Level: Advanced 73458122a4SShri Abhyankar .seealso: TS, TSEvent, TSSetEventHandler() 74458122a4SShri Abhyankar @*/ 75458122a4SShri Abhyankar PetscErrorCode TSSetPostEventIntervalStep(TS ts,PetscReal dt) 76458122a4SShri Abhyankar { 77458122a4SShri Abhyankar PetscFunctionBegin; 78458122a4SShri Abhyankar ts->event->timestep_posteventinterval = dt; 79458122a4SShri Abhyankar PetscFunctionReturn(0); 80458122a4SShri Abhyankar } 81458122a4SShri Abhyankar 82458122a4SShri Abhyankar /*@ 83e3005195SShri Abhyankar TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 84e3005195SShri Abhyankar 85e3005195SShri Abhyankar Logically Collective 86e3005195SShri Abhyankar 874165533cSJose E. Roman Input Parameters: 88e3005195SShri Abhyankar + ts - time integration context 89e3005195SShri Abhyankar . tol - scalar tolerance, PETSC_DECIDE to leave current value 90e3005195SShri Abhyankar - vtol - array of tolerances or NULL, used in preference to tol if present 91e3005195SShri Abhyankar 922589fa24SLisandro Dalcin Options Database Keys: 93147403d9SBarry Smith . -ts_event_tol <tol> - tolerance for event zero crossing 94e3005195SShri Abhyankar 95e3005195SShri Abhyankar Notes: 96f25fe674SLisandro Dalcin Must call TSSetEventHandler() before setting the tolerances. 97e3005195SShri Abhyankar 98e3005195SShri Abhyankar The size of vtol is equal to the number of events. 99e3005195SShri Abhyankar 100e3005195SShri Abhyankar Level: beginner 101e3005195SShri Abhyankar 102f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler() 103e3005195SShri Abhyankar @*/ 1046427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[]) 105e3005195SShri Abhyankar { 1066427ac75SLisandro Dalcin TSEvent event; 107e3005195SShri Abhyankar PetscInt i; 108e3005195SShri Abhyankar 109e3005195SShri Abhyankar PetscFunctionBegin; 1106427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1116427ac75SLisandro Dalcin if (vtol) PetscValidRealPointer(vtol,3); 1123c633725SBarry Smith PetscCheck(ts->event,PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()"); 1136427ac75SLisandro Dalcin 1146427ac75SLisandro Dalcin event = ts->event; 115e3005195SShri Abhyankar if (vtol) { 116e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 117e3005195SShri Abhyankar } else { 118e3005195SShri Abhyankar if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 119e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 120e3005195SShri Abhyankar } 121e3005195SShri Abhyankar } 122e3005195SShri Abhyankar PetscFunctionReturn(0); 123e3005195SShri Abhyankar } 124e3005195SShri Abhyankar 1252dc7a7e3SShri Abhyankar /*@C 126ac6a796dSBarry Smith TSSetEventHandler - Sets a function used for detecting events 1272dc7a7e3SShri Abhyankar 1282dc7a7e3SShri Abhyankar Logically Collective on TS 1292dc7a7e3SShri Abhyankar 1302dc7a7e3SShri Abhyankar Input Parameters: 1312dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 1322dc7a7e3SShri Abhyankar . nevents - number of local events 133d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 134d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 135d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 136d94325d3SShri Abhyankar event is detected (one for each event) 1376427ac75SLisandro Dalcin . eventhandler - event monitoring routine 1382dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 1392589fa24SLisandro Dalcin - ctx - [optional] user-defined context for private data for the 1402dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 1412dc7a7e3SShri Abhyankar context is desired) 1422dc7a7e3SShri Abhyankar 1436427ac75SLisandro Dalcin Calling sequence of eventhandler: 1443a88037aSBarry Smith PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 1452dc7a7e3SShri Abhyankar 1462dc7a7e3SShri Abhyankar Input Parameters: 1472dc7a7e3SShri Abhyankar + ts - the TS context 1482dc7a7e3SShri Abhyankar . t - current time 1492dc7a7e3SShri Abhyankar . U - current iterate 1506427ac75SLisandro Dalcin - ctx - [optional] context passed with eventhandler 1512dc7a7e3SShri Abhyankar 1522dc7a7e3SShri Abhyankar Output parameters: 153d94325d3SShri Abhyankar . fvalue - function value of events at time t 1542dc7a7e3SShri Abhyankar 1552dc7a7e3SShri Abhyankar Calling sequence of postevent: 1562589fa24SLisandro Dalcin PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 1572dc7a7e3SShri Abhyankar 1582dc7a7e3SShri Abhyankar Input Parameters: 1592dc7a7e3SShri Abhyankar + ts - the TS context 1602dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 1612dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 1622dc7a7e3SShri Abhyankar . t - current time 1632dc7a7e3SShri Abhyankar . U - current solution 164031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 1656427ac75SLisandro Dalcin - ctx - the context passed with eventhandler 1662dc7a7e3SShri Abhyankar 1672dc7a7e3SShri Abhyankar Level: intermediate 1682dc7a7e3SShri Abhyankar 1692dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 1702dc7a7e3SShri Abhyankar @*/ 1712589fa24SLisandro Dalcin PetscErrorCode TSSetEventHandler(TS ts,PetscInt nevents,PetscInt direction[],PetscBool terminate[],PetscErrorCode (*eventhandler)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void *ctx) 1722dc7a7e3SShri Abhyankar { 173d294eb03SHong Zhang TSAdapt adapt; 174d294eb03SHong Zhang PetscReal hmin; 1752dc7a7e3SShri Abhyankar TSEvent event; 176d94325d3SShri Abhyankar PetscInt i; 177006e6a18SShri Abhyankar PetscBool flg; 178a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE 179a6c783d2SShri Abhyankar PetscReal tol=1e-4; 180a6c783d2SShri Abhyankar #else 181d569cc17SSatish Balay PetscReal tol=1e-6; 182a6c783d2SShri Abhyankar #endif 1832dc7a7e3SShri Abhyankar 1842dc7a7e3SShri Abhyankar PetscFunctionBegin; 1856427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1860a82b154SShri if (nevents) { 187064a246eSJacob Faibussowitsch PetscValidIntPointer(direction,3); 188064a246eSJacob Faibussowitsch PetscValidBoolPointer(terminate,4); 1890a82b154SShri } 1906427ac75SLisandro Dalcin 1919566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&event)); 1929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->fvalue)); 1939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->fvalue_prev)); 1949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->fvalue_right)); 1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->zerocrossing)); 1969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->side)); 1979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->direction)); 1989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->terminate)); 1999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->vtol)); 200d94325d3SShri Abhyankar for (i=0; i < nevents; i++) { 201d94325d3SShri Abhyankar event->direction[i] = direction[i]; 202d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 203e2cdd850SShri Abhyankar event->zerocrossing[i] = PETSC_FALSE; 204e2cdd850SShri Abhyankar event->side[i] = 0; 205d94325d3SShri Abhyankar } 2069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nevents,&event->events_zero)); 2072589fa24SLisandro Dalcin event->nevents = nevents; 2086427ac75SLisandro Dalcin event->eventhandler = eventhandler; 2092dc7a7e3SShri Abhyankar event->postevent = postevent; 2106427ac75SLisandro Dalcin event->ctx = ctx; 211458122a4SShri Abhyankar event->timestep_posteventinterval = ts->time_step; 2129566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 2139566063dSJacob Faibussowitsch PetscCall(TSAdaptGetStepLimits(adapt,&hmin,NULL)); 214d294eb03SHong Zhang event->timestep_min = hmin; 2152dc7a7e3SShri Abhyankar 21602749585SLisandro Dalcin event->recsize = 8; /* Initial size of the recorder */ 217d0609cedSBarry Smith PetscOptionsBegin(((PetscObject)ts)->comm,((PetscObject)ts)->prefix,"TS Event options","TS"); 2182dc7a7e3SShri Abhyankar { 2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg)); 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_post_eventinterval_step","Time step after event interval","",event->timestep_posteventinterval,&event->timestep_posteventinterval,NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_post_event_step","Time step after event","",event->timestep_postevent,&event->timestep_postevent,NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_event_dt_min","Minimum time step considered for TSEvent","",event->timestep_min,&event->timestep_min,NULL)); 2252dc7a7e3SShri Abhyankar } 226d0609cedSBarry Smith PetscOptionsEnd(); 227a4ffd976SShri Abhyankar 2289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize,&event->recorder.time)); 2299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize,&event->recorder.stepnum)); 2309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize,&event->recorder.nevents)); 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->recsize,&event->recorder.eventidx)); 232a4ffd976SShri Abhyankar for (i=0; i < event->recsize; i++) { 2339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->nevents,&event->recorder.eventidx[i])); 234a4ffd976SShri Abhyankar } 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 */ 25102749585SLisandro Dalcin static PetscErrorCode TSEventRecorderResize(TSEvent event) 252a4ffd976SShri Abhyankar { 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)); 266a4ffd976SShri Abhyankar for (i=0; i < fact*event->recsize; i++) { 2679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(event->nevents,&eventidx[i])); 268a4ffd976SShri Abhyankar } 269a4ffd976SShri Abhyankar 270a4ffd976SShri Abhyankar /* Copy over data */ 2719566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(time,event->recorder.time,event->recsize)); 2729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(stepnum,event->recorder.stepnum,event->recsize)); 2739566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(nevents,event->recorder.nevents,event->recsize)); 274a4ffd976SShri Abhyankar for (i=0; i < event->recsize; i++) { 2759566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i])); 276a4ffd976SShri Abhyankar } 277a4ffd976SShri Abhyankar 278a4ffd976SShri Abhyankar /* Destroy old arrays */ 279a4ffd976SShri Abhyankar for (i=0; i < event->recsize; i++) { 2809566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.eventidx[i])); 281a4ffd976SShri Abhyankar } 2829566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.eventidx)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.nevents)); 2849566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.stepnum)); 2859566063dSJacob Faibussowitsch PetscCall(PetscFree(event->recorder.time)); 286a4ffd976SShri Abhyankar 287a4ffd976SShri Abhyankar /* Set pointers */ 288a4ffd976SShri Abhyankar event->recorder.time = time; 289a4ffd976SShri Abhyankar event->recorder.stepnum = stepnum; 290a4ffd976SShri Abhyankar event->recorder.nevents = nevents; 291a4ffd976SShri Abhyankar event->recorder.eventidx = eventidx; 292a4ffd976SShri Abhyankar 293a4ffd976SShri Abhyankar /* Double size */ 294a4ffd976SShri Abhyankar event->recsize *= fact; 295a4ffd976SShri Abhyankar 296a4ffd976SShri Abhyankar PetscFunctionReturn(0); 297a4ffd976SShri Abhyankar } 298a4ffd976SShri Abhyankar 299031fbad4SShri Abhyankar /* 300ac6a796dSBarry Smith Helper routine to handle user postevents and recording 301031fbad4SShri Abhyankar */ 3024597913aSLisandro Dalcin static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U) 3032dc7a7e3SShri Abhyankar { 3042dc7a7e3SShri Abhyankar TSEvent event = ts->event; 3052dc7a7e3SShri Abhyankar PetscBool terminate = PETSC_FALSE; 30628d5b5d6SLisandro Dalcin PetscBool restart = PETSC_FALSE; 307d0578d90SShri Abhyankar PetscInt i,ctr,stepnum; 3087324a0ffSLisandro Dalcin PetscBool inflag[2],outflag[2]; 3094597913aSLisandro Dalcin PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 3102dc7a7e3SShri Abhyankar 3112dc7a7e3SShri Abhyankar PetscFunctionBegin; 3122dc7a7e3SShri Abhyankar if (event->postevent) { 31328d5b5d6SLisandro Dalcin PetscObjectState state_prev,state_post; 3149566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U,&state_prev)); 3159566063dSJacob Faibussowitsch PetscCall((*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx)); 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)U,&state_post)); 31728d5b5d6SLisandro Dalcin if (state_prev != state_post) restart = PETSC_TRUE; 3182dc7a7e3SShri Abhyankar } 3194597913aSLisandro Dalcin 32028d5b5d6SLisandro Dalcin /* Handle termination events and step restart */ 32128d5b5d6SLisandro Dalcin for (i=0; i<event->nevents_zero; i++) if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE; 3227324a0ffSLisandro Dalcin inflag[0] = restart; inflag[1] = terminate; 3231c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(inflag,outflag,2,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm)); 3247324a0ffSLisandro Dalcin restart = outflag[0]; terminate = outflag[1]; 3259566063dSJacob Faibussowitsch if (restart) PetscCall(TSRestartStep(ts)); 3269566063dSJacob Faibussowitsch if (terminate) PetscCall(TSSetConvergedReason(ts,TS_CONVERGED_EVENT)); 3277324a0ffSLisandro Dalcin event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP; 328f7aea88cSShri Abhyankar 3294597913aSLisandro Dalcin /* Reset event residual functions as states might get changed by the postevent callback */ 330f443add6SLisandro Dalcin if (event->postevent) { 3319566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(U)); 3329566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts,t,U,event->fvalue,event->ctx)); 3339566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(U)); 334f443add6SLisandro Dalcin } 335f443add6SLisandro Dalcin 336f443add6SLisandro Dalcin /* Cache current time and event residual functions */ 337f443add6SLisandro Dalcin event->ptime_prev = t; 338f443add6SLisandro Dalcin for (i=0; i<event->nevents; i++) 339f443add6SLisandro Dalcin event->fvalue_prev[i] = event->fvalue[i]; 3404597913aSLisandro Dalcin 341d0578d90SShri Abhyankar /* Record the event in the event recorder */ 3429566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&stepnum)); 343f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 344a4ffd976SShri Abhyankar if (ctr == event->recsize) { 3459566063dSJacob Faibussowitsch PetscCall(TSEventRecorderResize(event)); 346a4ffd976SShri Abhyankar } 347f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 348d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 3494597913aSLisandro Dalcin event->recorder.nevents[ctr] = event->nevents_zero; 3504597913aSLisandro Dalcin for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 351f7aea88cSShri Abhyankar event->recorder.ctr++; 3522dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3532dc7a7e3SShri Abhyankar } 3542dc7a7e3SShri Abhyankar 35502749585SLisandro Dalcin /* Uses Anderson-Bjorck variant of regula falsi method */ 3569fbee547SJacob Faibussowitsch static inline PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt) 35738bf2713SShri Abhyankar { 35802749585SLisandro Dalcin PetscReal new_dt, scal = 1.0; 359e2cdd850SShri Abhyankar if (PetscRealPart(fleft)*PetscRealPart(f) < 0) { 360e2cdd850SShri Abhyankar if (side == 1) { 361a6c783d2SShri Abhyankar scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright); 362a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 363e2cdd850SShri Abhyankar } 36402749585SLisandro Dalcin new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft; 365e2cdd850SShri Abhyankar } else { 366e2cdd850SShri Abhyankar if (side == -1) { 367a6c783d2SShri Abhyankar scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft); 368a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 369e2cdd850SShri Abhyankar } 37002749585SLisandro Dalcin new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t; 371e2cdd850SShri Abhyankar } 37202749585SLisandro Dalcin return PetscMin(dt,new_dt); 37338bf2713SShri Abhyankar } 374e2cdd850SShri Abhyankar 375d294eb03SHong Zhang static PetscErrorCode TSEventDetection(TS ts) 376d294eb03SHong Zhang { 377d294eb03SHong Zhang TSEvent event = ts->event; 378d294eb03SHong Zhang PetscReal t; 379d294eb03SHong Zhang PetscInt i; 380d294eb03SHong Zhang PetscInt fvalue_sign,fvalueprev_sign; 381ea3dac1cSHong Zhang PetscInt in,out; 382d294eb03SHong Zhang 383d294eb03SHong Zhang PetscFunctionBegin; 3849566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&t)); 385d294eb03SHong Zhang for (i=0; i < event->nevents; i++) { 386d294eb03SHong Zhang if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 387d294eb03SHong Zhang if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 388d294eb03SHong Zhang event->status = TSEVENT_LOCATED_INTERVAL; 389d294eb03SHong Zhang if (event->monitor) { 390*63a3b9bcSJacob 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)); 391d294eb03SHong Zhang } 392d294eb03SHong Zhang continue; 393d294eb03SHong Zhang } 394ea3dac1cSHong Zhang if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */ 395d294eb03SHong Zhang fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 396d294eb03SHong Zhang fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 397d294eb03SHong Zhang if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) { 398d294eb03SHong Zhang if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 399d294eb03SHong Zhang event->status = TSEVENT_LOCATED_INTERVAL; 400d294eb03SHong Zhang if (event->monitor) { 401*63a3b9bcSJacob Faibussowitsch 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)); 402d294eb03SHong Zhang } 403d294eb03SHong Zhang } 404d294eb03SHong Zhang } 405d5f990dbSBarry Smith in = (PetscInt)event->status; 4061c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&in,&out,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts))); 407ea3dac1cSHong Zhang event->status = (TSEventStatus)out; 408d294eb03SHong Zhang PetscFunctionReturn(0); 409d294eb03SHong Zhang } 410d294eb03SHong Zhang 411d294eb03SHong Zhang static PetscErrorCode TSEventLocation(TS ts,PetscReal *dt) 412d294eb03SHong Zhang { 413d294eb03SHong Zhang TSEvent event = ts->event; 414d294eb03SHong Zhang PetscInt i; 415d294eb03SHong Zhang PetscReal t; 416ea3dac1cSHong Zhang PetscInt fvalue_sign,fvalueprev_sign; 417ea3dac1cSHong Zhang PetscInt rollback=0,in[2],out[2]; 418d294eb03SHong Zhang 419d294eb03SHong Zhang PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&t)); 421d294eb03SHong Zhang event->nevents_zero = 0; 422d294eb03SHong Zhang for (i=0; i < event->nevents; i++) { 423d294eb03SHong Zhang if (event->zerocrossing[i]) { 424d294eb03SHong 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 */ 425d294eb03SHong Zhang event->status = TSEVENT_ZERO; 426d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 427d294eb03SHong Zhang continue; 428d294eb03SHong Zhang } 429d294eb03SHong Zhang /* Compute new time step */ 430d294eb03SHong 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); 431d294eb03SHong Zhang fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 432ea3dac1cSHong Zhang fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 433d294eb03SHong Zhang switch (event->direction[i]) { 434d294eb03SHong Zhang case -1: 435d294eb03SHong Zhang if (fvalue_sign < 0) { 436ea3dac1cSHong Zhang rollback = 1; 437d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 438d294eb03SHong Zhang event->side[i] = 1; 439d294eb03SHong Zhang } 440d294eb03SHong Zhang break; 441d294eb03SHong Zhang case 1: 442d294eb03SHong Zhang if (fvalue_sign > 0) { 443ea3dac1cSHong Zhang rollback = 1; 444d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 445d294eb03SHong Zhang event->side[i] = 1; 446d294eb03SHong Zhang } 447d294eb03SHong Zhang break; 448d294eb03SHong Zhang case 0: 449ea3dac1cSHong Zhang if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */ 450ea3dac1cSHong Zhang rollback = 1; 451d294eb03SHong Zhang event->fvalue_right[i] = event->fvalue[i]; 452d294eb03SHong Zhang event->side[i] = 1; 453ea3dac1cSHong Zhang } 454d294eb03SHong Zhang break; 455d294eb03SHong Zhang } 456ea3dac1cSHong Zhang if (event->status == TSEVENT_PROCESSING) event->side[i] = -1; 457d294eb03SHong Zhang } 458d294eb03SHong Zhang } 459d5f990dbSBarry Smith in[0] = (PetscInt)event->status; in[1] = rollback; 4601c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts))); 461ea3dac1cSHong Zhang event->status = (TSEventStatus)out[0]; rollback = out[1]; 462ea3dac1cSHong 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 */ 463ea3dac1cSHong Zhang if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 464d294eb03SHong Zhang if (event->status == TSEVENT_ZERO) { 465ea3dac1cSHong Zhang for (i=0; i < event->nevents; i++) { 466ea3dac1cSHong Zhang if (event->zerocrossing[i]) { 467ea3dac1cSHong 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 */ 468ea3dac1cSHong Zhang event->events_zero[event->nevents_zero++] = i; 469ea3dac1cSHong Zhang if (event->monitor) { 470*63a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " zero crossing located at time %g\n",event->iterctr,i,(double)t)); 471ea3dac1cSHong Zhang } 472ea3dac1cSHong Zhang event->zerocrossing[i] = PETSC_FALSE; 473ea3dac1cSHong Zhang } 474ea3dac1cSHong Zhang } 475ea3dac1cSHong Zhang event->side[i] = 0; 476ea3dac1cSHong Zhang } 477d294eb03SHong Zhang } 478d294eb03SHong Zhang PetscFunctionReturn(0); 479d294eb03SHong Zhang } 48038bf2713SShri Abhyankar 4816427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts) 4822dc7a7e3SShri Abhyankar { 4836427ac75SLisandro Dalcin TSEvent event; 4842dc7a7e3SShri Abhyankar PetscReal t; 4852dc7a7e3SShri Abhyankar Vec U; 4862dc7a7e3SShri Abhyankar PetscInt i; 487ea3dac1cSHong Zhang PetscReal dt,dt_min,dt_reset = 0.0; 4882dc7a7e3SShri Abhyankar 4892dc7a7e3SShri Abhyankar PetscFunctionBegin; 4906427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4916427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 4926427ac75SLisandro Dalcin event = ts->event; 4932dc7a7e3SShri Abhyankar 4949566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&t)); 4959566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 4969566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&U)); 4972dc7a7e3SShri Abhyankar 4987dbe0728SLisandro Dalcin if (event->status == TSEVENT_NONE) { 4997dbe0728SLisandro Dalcin event->timestep_prev = dt; 500d294eb03SHong Zhang event->ptime_end = t; 5017dbe0728SLisandro Dalcin } 5022dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 503d294eb03SHong Zhang /* user has specified a PostEventInterval dt */ 504458122a4SShri Abhyankar dt = event->timestep_posteventinterval; 505e97c63d7SStefano Zampini if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 506e97c63d7SStefano Zampini PetscReal maxdt = ts->max_time-t; 507e97c63d7SStefano Zampini dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt); 508e97c63d7SStefano Zampini } 5099566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 5102dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 5112dc7a7e3SShri Abhyankar } 5122dc7a7e3SShri Abhyankar 5139566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(U)); 5149566063dSJacob Faibussowitsch PetscCall((*event->eventhandler)(ts,t,U,event->fvalue,event->ctx)); 5159566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(U)); 5169e12be75SShri Abhyankar 517d294eb03SHong Zhang /* Detect the events */ 5189566063dSJacob Faibussowitsch PetscCall(TSEventDetection(ts)); 519d294eb03SHong Zhang 520d294eb03SHong Zhang /* Locate the events */ 521d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) { 522d294eb03SHong Zhang /* Approach the zero crosing by setting a new step size */ 5239566063dSJacob Faibussowitsch PetscCall(TSEventLocation(ts,&dt)); 524d294eb03SHong Zhang /* Roll back when new events are detected */ 525d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL) { 5269566063dSJacob Faibussowitsch PetscCall(TSRollBack(ts)); 5279566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts,TS_CONVERGED_ITERATING)); 528d294eb03SHong Zhang event->iterctr++; 529006e6a18SShri Abhyankar } 5301c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts))); 531ea3dac1cSHong Zhang if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset; 5329566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt_min)); 533d294eb03SHong Zhang /* Found the zero crossing */ 5349e12be75SShri Abhyankar if (event->status == TSEVENT_ZERO) { 5359566063dSJacob Faibussowitsch PetscCall(TSPostEvent(ts,t,U)); 5369e12be75SShri Abhyankar 537e2cdd850SShri Abhyankar dt = event->ptime_end - t; 538ad508104SStefano Zampini if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */ 539ad508104SStefano Zampini dt = event->timestep_prev; 540ad508104SStefano Zampini event->status = TSEVENT_NONE; 541ad508104SStefano Zampini } 542d294eb03SHong Zhang if (event->timestep_postevent) { /* user has specified a PostEvent dt*/ 543d294eb03SHong Zhang dt = event->timestep_postevent; 544d294eb03SHong Zhang } 545e97c63d7SStefano Zampini if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 546e97c63d7SStefano Zampini PetscReal maxdt = ts->max_time-t; 547e97c63d7SStefano Zampini dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt); 548e97c63d7SStefano Zampini } 5499566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 55038bf2713SShri Abhyankar event->iterctr = 0; 5519e12be75SShri Abhyankar } 552d294eb03SHong Zhang /* Have not found the zero crosing yet */ 553d294eb03SHong Zhang if (event->status == TSEVENT_PROCESSING) { 554d294eb03SHong Zhang if (event->monitor) { 555*63a3b9bcSJacob Faibussowitsch 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)); 556d294eb03SHong Zhang } 557d294eb03SHong Zhang event->iterctr++; 558d294eb03SHong Zhang } 559d294eb03SHong Zhang } 560d294eb03SHong Zhang if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */ 5612dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 562e2cdd850SShri Abhyankar event->ptime_right = t; 5632dc7a7e3SShri Abhyankar } else { 564d294eb03SHong Zhang for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 56538bf2713SShri Abhyankar event->ptime_prev = t; 56638bf2713SShri Abhyankar } 5672dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 5682dc7a7e3SShri Abhyankar } 5692dc7a7e3SShri Abhyankar 5706427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts) 571d0578d90SShri Abhyankar { 5726427ac75SLisandro Dalcin TSEvent event; 573d0578d90SShri Abhyankar PetscReal t; 574d0578d90SShri Abhyankar Vec U; 575d0578d90SShri Abhyankar PetscInt ctr; 576d0578d90SShri Abhyankar PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 577d0578d90SShri Abhyankar 578d0578d90SShri Abhyankar PetscFunctionBegin; 5796427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5806427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 5816427ac75SLisandro Dalcin event = ts->event; 582d0578d90SShri Abhyankar 5839566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&t)); 5849566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&U)); 585d0578d90SShri Abhyankar 586d0578d90SShri Abhyankar ctr = event->recorder.ctr-1; 587bcbf8bb3SShri Abhyankar if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 588d0578d90SShri Abhyankar /* Call the user postevent function */ 589d0578d90SShri Abhyankar if (event->postevent) { 5909566063dSJacob Faibussowitsch PetscCall((*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx)); 591d0578d90SShri Abhyankar event->recorder.ctr--; 592d0578d90SShri Abhyankar } 593d0578d90SShri Abhyankar } 594d0578d90SShri Abhyankar 595d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 596d0578d90SShri Abhyankar PetscFunctionReturn(0); 597d0578d90SShri Abhyankar } 5981ea83e56SMiguel 5991ea83e56SMiguel /*@ 6001ea83e56SMiguel TSGetNumEvents - Get the numbers of events set 6011ea83e56SMiguel 6021ea83e56SMiguel Logically Collective 6031ea83e56SMiguel 6044165533cSJose E. Roman Input Parameter: 60599c90e12SSatish Balay . ts - the TS context 6061ea83e56SMiguel 6074165533cSJose E. Roman Output Parameter: 6081ea83e56SMiguel . nevents - number of events 6091ea83e56SMiguel 6101ea83e56SMiguel Level: intermediate 6111ea83e56SMiguel 6121ea83e56SMiguel .seealso: TSSetEventHandler() 6131ea83e56SMiguel 6141ea83e56SMiguel @*/ 6151ea83e56SMiguel PetscErrorCode TSGetNumEvents(TS ts,PetscInt * nevents) 6161ea83e56SMiguel { 6171ea83e56SMiguel PetscFunctionBegin; 6181ea83e56SMiguel *nevents = ts->event->nevents; 6191ea83e56SMiguel PetscFunctionReturn(0); 6201ea83e56SMiguel } 621