1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 22dc7a7e3SShri Abhyankar 32dc7a7e3SShri Abhyankar #undef __FUNCT__ 46427ac75SLisandro Dalcin #define __FUNCT__ "TSEventInitialize" 56fea3669SShri Abhyankar /* 66427ac75SLisandro Dalcin TSEventInitialize - Initializes TSEvent for TSSolve 76fea3669SShri Abhyankar */ 86427ac75SLisandro Dalcin PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U) 96fea3669SShri Abhyankar { 106fea3669SShri Abhyankar PetscErrorCode ierr; 116fea3669SShri Abhyankar 126fea3669SShri Abhyankar PetscFunctionBegin; 136427ac75SLisandro Dalcin if (!event) PetscFunctionReturn(0); 146427ac75SLisandro Dalcin PetscValidPointer(event,1); 156427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 166427ac75SLisandro Dalcin PetscValidHeaderSpecific(U,VEC_CLASSID,4); 176fea3669SShri Abhyankar event->ptime_prev = t; 1838bf2713SShri Abhyankar event->iterctr = 0; 196427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr); 206fea3669SShri Abhyankar /* Initialize the event recorder */ 216fea3669SShri Abhyankar event->recorder.ctr = 0; 226fea3669SShri Abhyankar PetscFunctionReturn(0); 236fea3669SShri Abhyankar } 246fea3669SShri Abhyankar 256fea3669SShri Abhyankar #undef __FUNCT__ 267dbe0728SLisandro Dalcin #define __FUNCT__ "TSEventDestroy" 277dbe0728SLisandro Dalcin PetscErrorCode TSEventDestroy(TSEvent *event) 287dbe0728SLisandro Dalcin { 297dbe0728SLisandro Dalcin PetscErrorCode ierr; 307dbe0728SLisandro Dalcin PetscInt i; 317dbe0728SLisandro Dalcin 327dbe0728SLisandro Dalcin PetscFunctionBegin; 337dbe0728SLisandro Dalcin PetscValidPointer(event,1); 347dbe0728SLisandro Dalcin if (!*event) PetscFunctionReturn(0); 357dbe0728SLisandro Dalcin ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 367dbe0728SLisandro Dalcin ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 37e2cdd850SShri Abhyankar ierr = PetscFree((*event)->fvalue_right);CHKERRQ(ierr); 38e2cdd850SShri Abhyankar ierr = PetscFree((*event)->zerocrossing);CHKERRQ(ierr); 39e2cdd850SShri Abhyankar ierr = PetscFree((*event)->side);CHKERRQ(ierr); 407dbe0728SLisandro Dalcin ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 417dbe0728SLisandro Dalcin ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 427dbe0728SLisandro Dalcin ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 437dbe0728SLisandro Dalcin ierr = PetscFree((*event)->vtol);CHKERRQ(ierr); 44a4ffd976SShri Abhyankar 45a4ffd976SShri Abhyankar for (i=0; i < (*event)->recsize; i++) { 467dbe0728SLisandro Dalcin ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr); 477dbe0728SLisandro Dalcin } 48a4ffd976SShri Abhyankar ierr = PetscFree((*event)->recorder.eventidx);CHKERRQ(ierr); 49a4ffd976SShri Abhyankar ierr = PetscFree((*event)->recorder.nevents);CHKERRQ(ierr); 50a4ffd976SShri Abhyankar ierr = PetscFree((*event)->recorder.stepnum);CHKERRQ(ierr); 51a4ffd976SShri Abhyankar ierr = PetscFree((*event)->recorder.time);CHKERRQ(ierr); 52a4ffd976SShri Abhyankar 537dbe0728SLisandro Dalcin ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr); 547dbe0728SLisandro Dalcin ierr = PetscFree(*event);CHKERRQ(ierr); 557dbe0728SLisandro Dalcin PetscFunctionReturn(0); 567dbe0728SLisandro Dalcin } 577dbe0728SLisandro Dalcin 587dbe0728SLisandro Dalcin #undef __FUNCT__ 59e3005195SShri Abhyankar #define __FUNCT__ "TSSetEventTolerances" 60e3005195SShri Abhyankar /*@ 61e3005195SShri Abhyankar TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 62e3005195SShri Abhyankar 63e3005195SShri Abhyankar Logically Collective 64e3005195SShri Abhyankar 65e3005195SShri Abhyankar Input Arguments: 66e3005195SShri Abhyankar + ts - time integration context 67e3005195SShri Abhyankar . tol - scalar tolerance, PETSC_DECIDE to leave current value 68e3005195SShri Abhyankar - vtol - array of tolerances or NULL, used in preference to tol if present 69e3005195SShri Abhyankar 702589fa24SLisandro Dalcin Options Database Keys: 712589fa24SLisandro Dalcin . -ts_event_tol <tol> tolerance for event zero crossing 72e3005195SShri Abhyankar 73e3005195SShri Abhyankar Notes: 74f25fe674SLisandro Dalcin Must call TSSetEventHandler() before setting the tolerances. 75e3005195SShri Abhyankar 76e3005195SShri Abhyankar The size of vtol is equal to the number of events. 77e3005195SShri Abhyankar 78e3005195SShri Abhyankar Level: beginner 79e3005195SShri Abhyankar 80f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler() 81e3005195SShri Abhyankar @*/ 826427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[]) 83e3005195SShri Abhyankar { 846427ac75SLisandro Dalcin TSEvent event; 85e3005195SShri Abhyankar PetscInt i; 86e3005195SShri Abhyankar 87e3005195SShri Abhyankar PetscFunctionBegin; 886427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 896427ac75SLisandro Dalcin if (vtol) PetscValidRealPointer(vtol,3); 90f25fe674SLisandro Dalcin if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()"); 916427ac75SLisandro Dalcin 926427ac75SLisandro Dalcin event = ts->event; 93e3005195SShri Abhyankar if (vtol) { 94e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 95e3005195SShri Abhyankar } else { 96e3005195SShri Abhyankar if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 97e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 98e3005195SShri Abhyankar } 99e3005195SShri Abhyankar } 100e3005195SShri Abhyankar PetscFunctionReturn(0); 101e3005195SShri Abhyankar } 102e3005195SShri Abhyankar 103e3005195SShri Abhyankar #undef __FUNCT__ 104f25fe674SLisandro Dalcin #define __FUNCT__ "TSSetEventHandler" 1052dc7a7e3SShri Abhyankar /*@C 106f25fe674SLisandro Dalcin TSSetEventHandler - Sets a monitoring function used for detecting events 1072dc7a7e3SShri Abhyankar 1082dc7a7e3SShri Abhyankar Logically Collective on TS 1092dc7a7e3SShri Abhyankar 1102dc7a7e3SShri Abhyankar Input Parameters: 1112dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 1122dc7a7e3SShri Abhyankar . nevents - number of local events 113d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 114d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 115d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 116d94325d3SShri Abhyankar event is detected (one for each event) 1176427ac75SLisandro Dalcin . eventhandler - event monitoring routine 1182dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 1192589fa24SLisandro Dalcin - ctx - [optional] user-defined context for private data for the 1202dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 1212dc7a7e3SShri Abhyankar context is desired) 1222dc7a7e3SShri Abhyankar 1236427ac75SLisandro Dalcin Calling sequence of eventhandler: 1242589fa24SLisandro Dalcin PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 1252dc7a7e3SShri Abhyankar 1262dc7a7e3SShri Abhyankar Input Parameters: 1272dc7a7e3SShri Abhyankar + ts - the TS context 1282dc7a7e3SShri Abhyankar . t - current time 1292dc7a7e3SShri Abhyankar . U - current iterate 1306427ac75SLisandro Dalcin - ctx - [optional] context passed with eventhandler 1312dc7a7e3SShri Abhyankar 1322dc7a7e3SShri Abhyankar Output parameters: 133d94325d3SShri Abhyankar . fvalue - function value of events at time t 1342dc7a7e3SShri Abhyankar 1352dc7a7e3SShri Abhyankar Calling sequence of postevent: 1362589fa24SLisandro Dalcin PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 1372dc7a7e3SShri Abhyankar 1382dc7a7e3SShri Abhyankar Input Parameters: 1392dc7a7e3SShri Abhyankar + ts - the TS context 1402dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 1412dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 1422dc7a7e3SShri Abhyankar . t - current time 1432dc7a7e3SShri Abhyankar . U - current solution 144031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 1456427ac75SLisandro Dalcin - ctx - the context passed with eventhandler 1462dc7a7e3SShri Abhyankar 1472dc7a7e3SShri Abhyankar Level: intermediate 1482dc7a7e3SShri Abhyankar 14902749585SLisandro Dalcin .keywords: TS, event, set 1502dc7a7e3SShri Abhyankar 1512dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 1522dc7a7e3SShri Abhyankar @*/ 1532589fa24SLisandro 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) 1542dc7a7e3SShri Abhyankar { 1552dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1562dc7a7e3SShri Abhyankar TSEvent event; 157d94325d3SShri Abhyankar PetscInt i; 158006e6a18SShri Abhyankar PetscBool flg; 159a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE 160a6c783d2SShri Abhyankar PetscReal tol=1e-4; 161a6c783d2SShri Abhyankar #else 162d569cc17SSatish Balay PetscReal tol=1e-6; 163a6c783d2SShri Abhyankar #endif 1642dc7a7e3SShri Abhyankar 1652dc7a7e3SShri Abhyankar PetscFunctionBegin; 1666427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1676427ac75SLisandro Dalcin PetscValidIntPointer(direction,2); 1686427ac75SLisandro Dalcin PetscValidIntPointer(terminate,3); 1696427ac75SLisandro Dalcin 1706427ac75SLisandro Dalcin ierr = PetscNewLog(ts,&event);CHKERRQ(ierr); 171854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 172854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 173e2cdd850SShri Abhyankar ierr = PetscMalloc1(nevents,&event->fvalue_right);CHKERRQ(ierr); 174e2cdd850SShri Abhyankar ierr = PetscMalloc1(nevents,&event->zerocrossing);CHKERRQ(ierr); 175e2cdd850SShri Abhyankar ierr = PetscMalloc1(nevents,&event->side);CHKERRQ(ierr); 176854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 177854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 178e3005195SShri Abhyankar ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr); 179d94325d3SShri Abhyankar for (i=0; i < nevents; i++) { 180d94325d3SShri Abhyankar event->direction[i] = direction[i]; 181d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 182e2cdd850SShri Abhyankar event->zerocrossing[i] = PETSC_FALSE; 183e2cdd850SShri Abhyankar event->side[i] = 0; 184d94325d3SShri Abhyankar } 185854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 1862589fa24SLisandro Dalcin event->nevents = nevents; 1876427ac75SLisandro Dalcin event->eventhandler = eventhandler; 1882dc7a7e3SShri Abhyankar event->postevent = postevent; 1896427ac75SLisandro Dalcin event->ctx = ctx; 1902dc7a7e3SShri Abhyankar 19102749585SLisandro Dalcin event->recsize = 8; /* Initial size of the recorder */ 192a9514d71SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr); 1932dc7a7e3SShri Abhyankar { 19402749585SLisandro Dalcin ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL);CHKERRQ(ierr); 195006e6a18SShri Abhyankar ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr); 196a4ffd976SShri Abhyankar ierr = PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);CHKERRQ(ierr); 1972dc7a7e3SShri Abhyankar } 1989e12be75SShri Abhyankar PetscOptionsEnd(); 199a4ffd976SShri Abhyankar 200a4ffd976SShri Abhyankar ierr = PetscMalloc1(event->recsize,&event->recorder.time);CHKERRQ(ierr); 201a4ffd976SShri Abhyankar ierr = PetscMalloc1(event->recsize,&event->recorder.stepnum);CHKERRQ(ierr); 202a4ffd976SShri Abhyankar ierr = PetscMalloc1(event->recsize,&event->recorder.nevents);CHKERRQ(ierr); 203a4ffd976SShri Abhyankar ierr = PetscMalloc1(event->recsize,&event->recorder.eventidx);CHKERRQ(ierr); 204a4ffd976SShri Abhyankar for (i=0; i < event->recsize; i++) { 205a4ffd976SShri Abhyankar ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 206a4ffd976SShri Abhyankar } 207a4ffd976SShri Abhyankar 208e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 2096427ac75SLisandro Dalcin if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);} 2109e12be75SShri Abhyankar 2116427ac75SLisandro Dalcin ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr); 212d94325d3SShri Abhyankar ts->event = event; 2132dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2142dc7a7e3SShri Abhyankar } 2152dc7a7e3SShri Abhyankar 216a4ffd976SShri Abhyankar #undef __FUNCT__ 217a4ffd976SShri Abhyankar #define __FUNCT__ "TSEventRecorderResize" 218a4ffd976SShri Abhyankar /* 219a4ffd976SShri Abhyankar TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize) 220a4ffd976SShri Abhyankar is reached. 221a4ffd976SShri Abhyankar */ 22202749585SLisandro Dalcin static PetscErrorCode TSEventRecorderResize(TSEvent event) 223a4ffd976SShri Abhyankar { 224a4ffd976SShri Abhyankar PetscErrorCode ierr; 225a4ffd976SShri Abhyankar PetscReal *time; 226a4ffd976SShri Abhyankar PetscInt *stepnum; 227a4ffd976SShri Abhyankar PetscInt *nevents; 228a4ffd976SShri Abhyankar PetscInt **eventidx; 229a4ffd976SShri Abhyankar PetscInt i,fact=2; 230a4ffd976SShri Abhyankar 231a4ffd976SShri Abhyankar PetscFunctionBegin; 232a4ffd976SShri Abhyankar 233a4ffd976SShri Abhyankar /* Create large arrays */ 234a4ffd976SShri Abhyankar ierr = PetscMalloc1(fact*event->recsize,&time);CHKERRQ(ierr); 235a4ffd976SShri Abhyankar ierr = PetscMalloc1(fact*event->recsize,&stepnum);CHKERRQ(ierr); 236a4ffd976SShri Abhyankar ierr = PetscMalloc1(fact*event->recsize,&nevents);CHKERRQ(ierr); 237a4ffd976SShri Abhyankar ierr = PetscMalloc1(fact*event->recsize,&eventidx);CHKERRQ(ierr); 238a4ffd976SShri Abhyankar for (i=0; i < fact*event->recsize; i++) { 239a4ffd976SShri Abhyankar ierr = PetscMalloc1(event->nevents,&eventidx[i]);CHKERRQ(ierr); 240a4ffd976SShri Abhyankar } 241a4ffd976SShri Abhyankar 242a4ffd976SShri Abhyankar /* Copy over data */ 243a4ffd976SShri Abhyankar ierr = PetscMemcpy(time,event->recorder.time,event->recsize*sizeof(PetscReal));CHKERRQ(ierr); 244a4ffd976SShri Abhyankar ierr = PetscMemcpy(stepnum,event->recorder.stepnum,event->recsize*sizeof(PetscInt));CHKERRQ(ierr); 245a4ffd976SShri Abhyankar ierr = PetscMemcpy(nevents,event->recorder.nevents,event->recsize*sizeof(PetscInt));CHKERRQ(ierr); 246a4ffd976SShri Abhyankar for (i=0; i < event->recsize; i++) { 247a4ffd976SShri Abhyankar ierr = PetscMemcpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]*sizeof(PetscInt));CHKERRQ(ierr); 248a4ffd976SShri Abhyankar } 249a4ffd976SShri Abhyankar 250a4ffd976SShri Abhyankar /* Destroy old arrays */ 251a4ffd976SShri Abhyankar for (i=0; i < event->recsize; i++) { 252a4ffd976SShri Abhyankar ierr = PetscFree(event->recorder.eventidx[i]);CHKERRQ(ierr); 253a4ffd976SShri Abhyankar } 254a4ffd976SShri Abhyankar ierr = PetscFree(event->recorder.eventidx);CHKERRQ(ierr); 255a4ffd976SShri Abhyankar ierr = PetscFree(event->recorder.nevents);CHKERRQ(ierr); 256a4ffd976SShri Abhyankar ierr = PetscFree(event->recorder.stepnum);CHKERRQ(ierr); 257a4ffd976SShri Abhyankar ierr = PetscFree(event->recorder.time);CHKERRQ(ierr); 258a4ffd976SShri Abhyankar 259a4ffd976SShri Abhyankar /* Set pointers */ 260a4ffd976SShri Abhyankar event->recorder.time = time; 261a4ffd976SShri Abhyankar event->recorder.stepnum = stepnum; 262a4ffd976SShri Abhyankar event->recorder.nevents = nevents; 263a4ffd976SShri Abhyankar event->recorder.eventidx = eventidx; 264a4ffd976SShri Abhyankar 265a4ffd976SShri Abhyankar /* Double size */ 266a4ffd976SShri Abhyankar event->recsize *= fact; 267a4ffd976SShri Abhyankar 268a4ffd976SShri Abhyankar PetscFunctionReturn(0); 269a4ffd976SShri Abhyankar } 270a4ffd976SShri Abhyankar 271031fbad4SShri Abhyankar /* 2724597913aSLisandro Dalcin Helper rutine to handle user postenvents and recording 273031fbad4SShri Abhyankar */ 274031fbad4SShri Abhyankar #undef __FUNCT__ 275031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent" 2764597913aSLisandro Dalcin static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U) 2772dc7a7e3SShri Abhyankar { 2782dc7a7e3SShri Abhyankar PetscErrorCode ierr; 2792dc7a7e3SShri Abhyankar TSEvent event = ts->event; 2802dc7a7e3SShri Abhyankar PetscBool terminate = PETSC_FALSE; 281d0578d90SShri Abhyankar PetscInt i,ctr,stepnum; 2822dc7a7e3SShri Abhyankar PetscBool ts_terminate; 2834597913aSLisandro Dalcin PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 2842dc7a7e3SShri Abhyankar 2852dc7a7e3SShri Abhyankar PetscFunctionBegin; 2862dc7a7e3SShri Abhyankar if (event->postevent) { 2874597913aSLisandro Dalcin ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 2882dc7a7e3SShri Abhyankar } 2894597913aSLisandro Dalcin 2904597913aSLisandro Dalcin /* Handle termination events */ 2914597913aSLisandro Dalcin for (i=0; i < event->nevents_zero; i++) { 2924597913aSLisandro Dalcin terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]); 2932dc7a7e3SShri Abhyankar } 294b2566f29SBarry Smith ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 2959e12be75SShri Abhyankar if (ts_terminate) { 2962dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 2972dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2982dc7a7e3SShri Abhyankar } else { 2992dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 3002dc7a7e3SShri Abhyankar } 301f7aea88cSShri Abhyankar 3024597913aSLisandro Dalcin event->ptime_prev = t; 3034597913aSLisandro Dalcin /* Reset event residual functions as states might get changed by the postevent callback */ 3044597913aSLisandro Dalcin if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);} 3054597913aSLisandro Dalcin /* Cache current event residual functions */ 3064597913aSLisandro Dalcin for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 3074597913aSLisandro Dalcin 308d0578d90SShri Abhyankar /* Record the event in the event recorder */ 309d0578d90SShri Abhyankar ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr); 310f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 311a4ffd976SShri Abhyankar if (ctr == event->recsize) { 312a4ffd976SShri Abhyankar ierr = TSEventRecorderResize(event);CHKERRQ(ierr); 313a4ffd976SShri Abhyankar } 314f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 315d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 3164597913aSLisandro Dalcin event->recorder.nevents[ctr] = event->nevents_zero; 3174597913aSLisandro Dalcin for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 318f7aea88cSShri Abhyankar event->recorder.ctr++; 3192dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3202dc7a7e3SShri Abhyankar } 3212dc7a7e3SShri Abhyankar 32202749585SLisandro Dalcin /* Uses Anderson-Bjorck variant of regula falsi method */ 323a3a8645aSShri Abhyankar PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt) 32438bf2713SShri Abhyankar { 32502749585SLisandro Dalcin PetscReal new_dt, scal = 1.0; 326e2cdd850SShri Abhyankar if (PetscRealPart(fleft)*PetscRealPart(f) < 0) { 327e2cdd850SShri Abhyankar if (side == 1) { 328a6c783d2SShri Abhyankar scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright); 329a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 330e2cdd850SShri Abhyankar } 33102749585SLisandro Dalcin new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft; 332e2cdd850SShri Abhyankar } else { 333e2cdd850SShri Abhyankar if (side == -1) { 334a6c783d2SShri Abhyankar scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft); 335a6c783d2SShri Abhyankar if (scal < PETSC_SMALL) scal = 0.5; 336e2cdd850SShri Abhyankar } 33702749585SLisandro Dalcin new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t; 338e2cdd850SShri Abhyankar } 33902749585SLisandro Dalcin return PetscMin(dt,new_dt); 34038bf2713SShri Abhyankar } 341e2cdd850SShri Abhyankar 34238bf2713SShri Abhyankar 3432dc7a7e3SShri Abhyankar #undef __FUNCT__ 3446427ac75SLisandro Dalcin #define __FUNCT__ "TSEventHandler" 3456427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts) 3462dc7a7e3SShri Abhyankar { 3472dc7a7e3SShri Abhyankar PetscErrorCode ierr; 3486427ac75SLisandro Dalcin TSEvent event; 3492dc7a7e3SShri Abhyankar PetscReal t; 3502dc7a7e3SShri Abhyankar Vec U; 3512dc7a7e3SShri Abhyankar PetscInt i; 3527dbe0728SLisandro Dalcin PetscReal dt,dt_min; 3532dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 3549e12be75SShri Abhyankar PetscInt fvalue_sign,fvalueprev_sign; 3552dc7a7e3SShri Abhyankar 3562dc7a7e3SShri Abhyankar PetscFunctionBegin; 3576427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3586427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 3596427ac75SLisandro Dalcin event = ts->event; 3602dc7a7e3SShri Abhyankar 3612dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 3627dbe0728SLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 3632dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 3642dc7a7e3SShri Abhyankar 3657dbe0728SLisandro Dalcin if (event->status == TSEVENT_NONE) { 366*2d5bd56dSLisandro Dalcin if (ts->steps == 1) /* After first successful step */ 367*2d5bd56dSLisandro Dalcin event->timestep_orig = ts->ptime - ts->ptime_prev; 3687dbe0728SLisandro Dalcin event->timestep_prev = dt; 3697dbe0728SLisandro Dalcin } 3707dbe0728SLisandro Dalcin 3712dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 372*2d5bd56dSLisandro Dalcin /* Restore time step */ 373*2d5bd56dSLisandro Dalcin dt = PetscMin(event->timestep_orig,event->timestep_prev); 3747dbe0728SLisandro Dalcin ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 3752dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 3762dc7a7e3SShri Abhyankar } 3772dc7a7e3SShri Abhyankar 3782dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 3797dbe0728SLisandro Dalcin event->ptime_end = t; 3802dc7a7e3SShri Abhyankar } 3812dc7a7e3SShri Abhyankar 3826427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); 3839e12be75SShri Abhyankar 3842dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 385e3005195SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 3862dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 387e2cdd850SShri Abhyankar event->fvalue_right[i] = event->fvalue[i]; 3889e12be75SShri Abhyankar continue; 389006e6a18SShri Abhyankar } 390d94325d3SShri Abhyankar fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 391d94325d3SShri Abhyankar fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 392e3005195SShri Abhyankar if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) { 393d94325d3SShri Abhyankar switch (event->direction[i]) { 394d94325d3SShri Abhyankar case -1: 395d94325d3SShri Abhyankar if (fvalue_sign < 0) { 3962dc7a7e3SShri Abhyankar rollback = 1; 397e2cdd850SShri Abhyankar 398e2cdd850SShri Abhyankar /* Compute new time step */ 399e2cdd850SShri Abhyankar dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); 400e2cdd850SShri Abhyankar 4016427ac75SLisandro Dalcin if (event->monitor) { 40238bf2713SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 403006e6a18SShri Abhyankar } 404e2cdd850SShri Abhyankar event->fvalue_right[i] = event->fvalue[i]; 405e2cdd850SShri Abhyankar event->side[i] = 1; 406e2cdd850SShri Abhyankar 407e2cdd850SShri Abhyankar if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 408e2cdd850SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 4092dc7a7e3SShri Abhyankar } 410d94325d3SShri Abhyankar break; 411d94325d3SShri Abhyankar case 1: 412d94325d3SShri Abhyankar if (fvalue_sign > 0) { 413d94325d3SShri Abhyankar rollback = 1; 414e2cdd850SShri Abhyankar 415e2cdd850SShri Abhyankar /* Compute new time step */ 416e2cdd850SShri Abhyankar dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); 417e2cdd850SShri Abhyankar 4186427ac75SLisandro Dalcin if (event->monitor) { 41938bf2713SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 420006e6a18SShri Abhyankar } 421e2cdd850SShri Abhyankar event->fvalue_right[i] = event->fvalue[i]; 422e2cdd850SShri Abhyankar event->side[i] = 1; 423e2cdd850SShri Abhyankar 424e2cdd850SShri Abhyankar if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 425e2cdd850SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 4262dc7a7e3SShri Abhyankar } 427d94325d3SShri Abhyankar break; 428d94325d3SShri Abhyankar case 0: 429d94325d3SShri Abhyankar rollback = 1; 430e2cdd850SShri Abhyankar 431e2cdd850SShri Abhyankar /* Compute new time step */ 432e2cdd850SShri Abhyankar dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); 433e2cdd850SShri Abhyankar 4346427ac75SLisandro Dalcin if (event->monitor) { 43538bf2713SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 436006e6a18SShri Abhyankar } 437e2cdd850SShri Abhyankar event->fvalue_right[i] = event->fvalue[i]; 438e2cdd850SShri Abhyankar event->side[i] = 1; 439e2cdd850SShri Abhyankar 440e2cdd850SShri Abhyankar if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 441e2cdd850SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 442e2cdd850SShri Abhyankar 443d94325d3SShri Abhyankar break; 444d94325d3SShri Abhyankar } 445d94325d3SShri Abhyankar } 446d94325d3SShri Abhyankar } 4477dbe0728SLisandro Dalcin 4482589fa24SLisandro Dalcin in[0] = event->status; in[1] = rollback; 4497dbe0728SLisandro Dalcin ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 4502589fa24SLisandro Dalcin event->status = (TSEventStatus)out[0]; rollback = out[1]; 4512589fa24SLisandro Dalcin if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 4522dc7a7e3SShri Abhyankar 4537dbe0728SLisandro Dalcin event->nevents_zero = 0; 4549e12be75SShri Abhyankar if (event->status == TSEVENT_ZERO) { 4559e12be75SShri Abhyankar for (i=0; i < event->nevents; i++) { 4569e12be75SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 4579e12be75SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 4586427ac75SLisandro Dalcin if (event->monitor) { 45938bf2713SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);CHKERRQ(ierr); 4609e12be75SShri Abhyankar } 461e2cdd850SShri Abhyankar event->zerocrossing[i] = PETSC_FALSE; 4629e12be75SShri Abhyankar } 463e2cdd850SShri Abhyankar event->side[i] = 0; 4649e12be75SShri Abhyankar } 4654597913aSLisandro Dalcin ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr); 4669e12be75SShri Abhyankar 467e2cdd850SShri Abhyankar dt = event->ptime_end - t; 468*2d5bd56dSLisandro Dalcin if (PetscAbsReal(dt) < PETSC_SMALL) dt += PetscMin(event->timestep_orig,event->timestep_prev); /* XXX Should be done better */ 4699e12be75SShri Abhyankar ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 47038bf2713SShri Abhyankar event->iterctr = 0; 4719e12be75SShri Abhyankar PetscFunctionReturn(0); 4729e12be75SShri Abhyankar } 4739e12be75SShri Abhyankar 4742dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 4752dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 4762589fa24SLisandro Dalcin ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr); 4772dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 478e2cdd850SShri Abhyankar event->ptime_right = t; 4792dc7a7e3SShri Abhyankar } else { 4802dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 481e2cdd850SShri Abhyankar if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) { 482e2cdd850SShri Abhyankar /* Compute new time step */ 483e2cdd850SShri Abhyankar dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); 484e2cdd850SShri Abhyankar event->side[i] = -1; 485e2cdd850SShri Abhyankar } 4862dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 4872dc7a7e3SShri Abhyankar } 488e2cdd850SShri Abhyankar if (event->monitor && event->status == TSEVENT_PROCESSING) { 48938bf2713SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Stepping forward as no event detected in interval [%g - %g]\n",event->iterctr,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 4902dc7a7e3SShri Abhyankar } 49138bf2713SShri Abhyankar event->ptime_prev = t; 49238bf2713SShri Abhyankar } 49338bf2713SShri Abhyankar 49438bf2713SShri Abhyankar if (event->status == TSEVENT_PROCESSING) event->iterctr++; 4959e12be75SShri Abhyankar 4967dbe0728SLisandro Dalcin ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 4977dbe0728SLisandro Dalcin ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr); 4982dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 4992dc7a7e3SShri Abhyankar } 5002dc7a7e3SShri Abhyankar 501d0578d90SShri Abhyankar #undef __FUNCT__ 5026427ac75SLisandro Dalcin #define __FUNCT__ "TSAdjointEventHandler" 5036427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts) 504d0578d90SShri Abhyankar { 505d0578d90SShri Abhyankar PetscErrorCode ierr; 5066427ac75SLisandro Dalcin TSEvent event; 507d0578d90SShri Abhyankar PetscReal t; 508d0578d90SShri Abhyankar Vec U; 509d0578d90SShri Abhyankar PetscInt ctr; 510d0578d90SShri Abhyankar PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 511d0578d90SShri Abhyankar 512d0578d90SShri Abhyankar PetscFunctionBegin; 5136427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 5146427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 5156427ac75SLisandro Dalcin event = ts->event; 516d0578d90SShri Abhyankar 517d0578d90SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 518d0578d90SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 519d0578d90SShri Abhyankar 520d0578d90SShri Abhyankar ctr = event->recorder.ctr-1; 521bcbf8bb3SShri Abhyankar if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 522d0578d90SShri Abhyankar /* Call the user postevent function */ 523d0578d90SShri Abhyankar if (event->postevent) { 5246427ac75SLisandro Dalcin ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 525d0578d90SShri Abhyankar event->recorder.ctr--; 526d0578d90SShri Abhyankar } 527d0578d90SShri Abhyankar } 528d0578d90SShri Abhyankar 529d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 530d0578d90SShri Abhyankar PetscFunctionReturn(0); 531d0578d90SShri Abhyankar } 532