12dc7a7e3SShri Abhyankar 22dc7a7e3SShri Abhyankar #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 32dc7a7e3SShri Abhyankar 42dc7a7e3SShri Abhyankar #undef __FUNCT__ 52dc7a7e3SShri Abhyankar #define __FUNCT__ "TSSetEventMonitor" 62dc7a7e3SShri Abhyankar /*@C 72dc7a7e3SShri Abhyankar TSSetEventMonitor - Sets a monitoring function used for detecting events 82dc7a7e3SShri Abhyankar 92dc7a7e3SShri Abhyankar Logically Collective on TS 102dc7a7e3SShri Abhyankar 112dc7a7e3SShri Abhyankar Input Parameters: 122dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 132dc7a7e3SShri Abhyankar . nevents - number of local events 14d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 15d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 16d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 17d94325d3SShri Abhyankar event is detected (one for each event) 182dc7a7e3SShri Abhyankar . eventmonitor - event monitoring routine 192dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 202dc7a7e3SShri Abhyankar - mectx - [optional] user-defined context for private data for the 212dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 222dc7a7e3SShri Abhyankar context is desired) 232dc7a7e3SShri Abhyankar 242dc7a7e3SShri Abhyankar Calling sequence of eventmonitor: 25d94325d3SShri Abhyankar PetscErrorCode EventMonitor(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void* mectx) 262dc7a7e3SShri Abhyankar 272dc7a7e3SShri Abhyankar Input Parameters: 282dc7a7e3SShri Abhyankar + ts - the TS context 292dc7a7e3SShri Abhyankar . t - current time 302dc7a7e3SShri Abhyankar . U - current iterate 312dc7a7e3SShri Abhyankar - ctx - [optional] context passed with eventmonitor 322dc7a7e3SShri Abhyankar 332dc7a7e3SShri Abhyankar Output parameters: 34d94325d3SShri Abhyankar . fvalue - function value of events at time t 352dc7a7e3SShri Abhyankar 362dc7a7e3SShri Abhyankar Calling sequence of postevent: 37031fbad4SShri Abhyankar PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero, PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 382dc7a7e3SShri Abhyankar 392dc7a7e3SShri Abhyankar Input Parameters: 402dc7a7e3SShri Abhyankar + ts - the TS context 412dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 422dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 432dc7a7e3SShri Abhyankar . t - current time 442dc7a7e3SShri Abhyankar . U - current solution 45031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 462dc7a7e3SShri Abhyankar - ctx - the context passed with eventmonitor 472dc7a7e3SShri Abhyankar 482dc7a7e3SShri Abhyankar Level: intermediate 492dc7a7e3SShri Abhyankar 502dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor 512dc7a7e3SShri Abhyankar 522dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 532dc7a7e3SShri Abhyankar @*/ 54031fbad4SShri Abhyankar PetscErrorCode TSSetEventMonitor(TS ts,PetscInt nevents,PetscInt *direction,PetscBool *terminate,PetscErrorCode (*eventmonitor)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void *mectx) 552dc7a7e3SShri Abhyankar { 562dc7a7e3SShri Abhyankar PetscErrorCode ierr; 572dc7a7e3SShri Abhyankar PetscReal t; 582dc7a7e3SShri Abhyankar Vec U; 592dc7a7e3SShri Abhyankar TSEvent event; 60d94325d3SShri Abhyankar PetscInt i; 612dc7a7e3SShri Abhyankar 622dc7a7e3SShri Abhyankar PetscFunctionBegin; 6342ea6711SShri ierr = PetscNew(&event);CHKERRQ(ierr); 64854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 65854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 66854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 67854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 68d94325d3SShri Abhyankar for (i=0; i < nevents; i++) { 69d94325d3SShri Abhyankar event->direction[i] = direction[i]; 70d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 71d94325d3SShri Abhyankar } 72854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 732dc7a7e3SShri Abhyankar event->monitor = eventmonitor; 742dc7a7e3SShri Abhyankar event->postevent = postevent; 752dc7a7e3SShri Abhyankar event->monitorcontext = (void*)mectx; 762dc7a7e3SShri Abhyankar event->nevents = nevents; 772dc7a7e3SShri Abhyankar 782dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 792dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&event->initial_timestep);CHKERRQ(ierr); 802dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 812dc7a7e3SShri Abhyankar event->ptime_prev = t; 82d94325d3SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,mectx);CHKERRQ(ierr); 832dc7a7e3SShri Abhyankar ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"TS Event options","");CHKERRQ(ierr); 842dc7a7e3SShri Abhyankar { 852dc7a7e3SShri Abhyankar event->tol = 1.0e-6; 862dc7a7e3SShri Abhyankar ierr = PetscOptionsReal("-ts_event_tol","","",event->tol,&event->tol,NULL);CHKERRQ(ierr); 872dc7a7e3SShri Abhyankar } 882dc7a7e3SShri Abhyankar ierr = PetscOptionsEnd();CHKERRQ(ierr); 89d94325d3SShri Abhyankar ts->event = event; 902dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 912dc7a7e3SShri Abhyankar } 922dc7a7e3SShri Abhyankar 932dc7a7e3SShri Abhyankar #undef __FUNCT__ 942dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent" 95031fbad4SShri Abhyankar /* 96031fbad4SShri Abhyankar TSPostEvent - Does post event processing by calling the user-defined postevent function 97031fbad4SShri Abhyankar 98031fbad4SShri Abhyankar Logically Collective on TS 99031fbad4SShri Abhyankar 100031fbad4SShri Abhyankar Input Parameters: 101031fbad4SShri Abhyankar + ts - the TS context 102031fbad4SShri Abhyankar . nevents_zero - number of local events whose event function is zero 103031fbad4SShri Abhyankar . events_zero - indices of local events which have reached zero 104031fbad4SShri Abhyankar . t - current time 105031fbad4SShri Abhyankar . U - current solution 106031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 107031fbad4SShri Abhyankar - ctx - the context passed with eventmonitor 108031fbad4SShri Abhyankar 109031fbad4SShri Abhyankar Level: intermediate 110031fbad4SShri Abhyankar 111031fbad4SShri Abhyankar .keywords: TS, event, set, monitor 112031fbad4SShri Abhyankar 113031fbad4SShri Abhyankar .seealso: TSSetEventMonitor(),TSEvent 114031fbad4SShri Abhyankar */ 115031fbad4SShri Abhyankar #undef __FUNCT__ 116031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent" 117031fbad4SShri Abhyankar PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void *ctx) 1182dc7a7e3SShri Abhyankar { 1192dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1202dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1212dc7a7e3SShri Abhyankar PetscBool terminate=PETSC_FALSE; 1222dc7a7e3SShri Abhyankar PetscInt i; 1232dc7a7e3SShri Abhyankar PetscBool ts_terminate; 1242dc7a7e3SShri Abhyankar 1252dc7a7e3SShri Abhyankar PetscFunctionBegin; 1262dc7a7e3SShri Abhyankar if (event->postevent) { 127031fbad4SShri Abhyankar ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,forwardsolve,ctx);CHKERRQ(ierr); 1282dc7a7e3SShri Abhyankar } 1292dc7a7e3SShri Abhyankar for(i = 0; i < nevents_zero;i++) { 130e105d053SSatish Balay terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]); 1312dc7a7e3SShri Abhyankar } 132b4549700SJed Brown ierr = MPI_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 1332dc7a7e3SShri Abhyankar if (terminate) { 1342dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 1352dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 1362dc7a7e3SShri Abhyankar } else { 1372dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 1382dc7a7e3SShri Abhyankar } 139*73967516SShri Abhyankar /* Reset the event residual functions as states might get changed by the postevent callback */ 140*73967516SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,event->monitorcontext);CHKERRQ(ierr); 141*73967516SShri Abhyankar event->ptime_prev = t; 142*73967516SShri Abhyankar 1432dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1442dc7a7e3SShri Abhyankar } 1452dc7a7e3SShri Abhyankar 1462dc7a7e3SShri Abhyankar #undef __FUNCT__ 1472dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitorDestroy" 1482dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitorDestroy(TSEvent *event) 1492dc7a7e3SShri Abhyankar { 1502dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1512dc7a7e3SShri Abhyankar 1522dc7a7e3SShri Abhyankar PetscFunctionBegin; 1532dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 1542dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 1552dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 156d94325d3SShri Abhyankar ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 1572dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 1582dc7a7e3SShri Abhyankar ierr = PetscFree(*event);CHKERRQ(ierr); 1592dc7a7e3SShri Abhyankar *event = NULL; 1602dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1612dc7a7e3SShri Abhyankar } 1622dc7a7e3SShri Abhyankar 1632dc7a7e3SShri Abhyankar #undef __FUNCT__ 1642dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitor" 1652dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitor(TS ts) 1662dc7a7e3SShri Abhyankar { 1672dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1682dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1692dc7a7e3SShri Abhyankar PetscReal t; 1702dc7a7e3SShri Abhyankar Vec U; 1712dc7a7e3SShri Abhyankar PetscInt i; 1722dc7a7e3SShri Abhyankar PetscReal dt; 173b4549700SJed Brown TSEventStatus status = event->status; 1742dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 175031fbad4SShri Abhyankar PetscBool forwardsolve=PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 1762dc7a7e3SShri Abhyankar 1772dc7a7e3SShri Abhyankar PetscFunctionBegin; 1782dc7a7e3SShri Abhyankar 1792dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 1802dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 1812dc7a7e3SShri Abhyankar 1822dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 1832dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 1842dc7a7e3SShri Abhyankar /* Take initial time step */ 1852dc7a7e3SShri Abhyankar dt = event->initial_timestep; 1862dc7a7e3SShri Abhyankar ts->time_step = dt; 1872dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 1882dc7a7e3SShri Abhyankar } 1892dc7a7e3SShri Abhyankar 1902dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 1912dc7a7e3SShri Abhyankar event->tstepend = t; 1922dc7a7e3SShri Abhyankar } 1932dc7a7e3SShri Abhyankar 1942dc7a7e3SShri Abhyankar event->nevents_zero = 0; 1952dc7a7e3SShri Abhyankar 196d94325d3SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue,event->monitorcontext);CHKERRQ(ierr); 1972dc7a7e3SShri Abhyankar if (event->status != TSEVENT_NONE) { 1982dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 199e105d053SSatish Balay if (PetscAbsScalar(event->fvalue[i]) < event->tol) { 2002dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 2012dc7a7e3SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 2022dc7a7e3SShri Abhyankar } 2032dc7a7e3SShri Abhyankar } 2042dc7a7e3SShri Abhyankar } 2052dc7a7e3SShri Abhyankar 2062dc7a7e3SShri Abhyankar status = event->status; 207b4549700SJed Brown ierr = MPI_Allreduce((PetscEnum*)&status,(PetscEnum*)&event->status,1,MPIU_ENUM,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 2082dc7a7e3SShri Abhyankar 2092dc7a7e3SShri Abhyankar if (event->status == TSEVENT_ZERO) { 2102dc7a7e3SShri Abhyankar dt = event->tstepend-t; 2112dc7a7e3SShri Abhyankar ts->time_step = dt; 212031fbad4SShri Abhyankar ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->monitorcontext);CHKERRQ(ierr); 2132dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2142dc7a7e3SShri Abhyankar } 2152dc7a7e3SShri Abhyankar 2162dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 217d94325d3SShri Abhyankar PetscInt fvalue_sign,fvalueprev_sign; 218d94325d3SShri Abhyankar fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 219d94325d3SShri Abhyankar fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 2207a728e9fSShri Abhyankar if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) { 221d94325d3SShri Abhyankar switch (event->direction[i]) { 222d94325d3SShri Abhyankar case -1: 223d94325d3SShri Abhyankar if (fvalue_sign < 0) { 2242dc7a7e3SShri Abhyankar rollback = 1; 2252dc7a7e3SShri Abhyankar /* Compute linearly interpolated new time step */ 226e105d053SSatish Balay dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 2272dc7a7e3SShri Abhyankar } 228d94325d3SShri Abhyankar break; 229d94325d3SShri Abhyankar case 1: 230d94325d3SShri Abhyankar if (fvalue_sign > 0) { 231d94325d3SShri Abhyankar rollback = 1; 232d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 233d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 2342dc7a7e3SShri Abhyankar } 235d94325d3SShri Abhyankar break; 236d94325d3SShri Abhyankar case 0: 237d94325d3SShri Abhyankar rollback = 1; 238d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 239d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 240d94325d3SShri Abhyankar break; 241d94325d3SShri Abhyankar } 242d94325d3SShri Abhyankar } 243d94325d3SShri Abhyankar } 244d94325d3SShri Abhyankar if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 245d94325d3SShri Abhyankar 2462dc7a7e3SShri Abhyankar in[0] = event->status; 2472dc7a7e3SShri Abhyankar in[1] = rollback; 2482dc7a7e3SShri Abhyankar ierr = MPI_Allreduce(in,out,2,MPIU_INT,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 2492dc7a7e3SShri Abhyankar 2502dc7a7e3SShri Abhyankar rollback = out[1]; 2512dc7a7e3SShri Abhyankar if (rollback) { 2522dc7a7e3SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 2532dc7a7e3SShri Abhyankar } 2542dc7a7e3SShri Abhyankar 2552dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 2562dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 2572dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 2582dc7a7e3SShri Abhyankar } else { 2592dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 2602dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 2612dc7a7e3SShri Abhyankar } 2622dc7a7e3SShri Abhyankar event->ptime_prev = t; 2632dc7a7e3SShri Abhyankar if (event->status == TSEVENT_PROCESSING) { 2642dc7a7e3SShri Abhyankar dt = event->tstepend - event->ptime_prev; 2652dc7a7e3SShri Abhyankar } 2662dc7a7e3SShri Abhyankar } 267e6b5be14SSatish Balay ierr = MPI_Allreduce(&dt,&(ts->time_step),1,MPIU_REAL,MPI_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 2682dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2692dc7a7e3SShri Abhyankar } 2702dc7a7e3SShri Abhyankar 271