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 142dc7a7e3SShri Abhyankar . eventmonitor - event monitoring routine 152dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 162dc7a7e3SShri Abhyankar - mectx - [optional] user-defined context for private data for the 172dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 182dc7a7e3SShri Abhyankar context is desired) 192dc7a7e3SShri Abhyankar 202dc7a7e3SShri Abhyankar Calling sequence of eventmonitor: 212dc7a7e3SShri Abhyankar PetscErrorCode EventMonitor(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,PetscInt *direction,PetscBool *terminate,void* mectx) 222dc7a7e3SShri Abhyankar 232dc7a7e3SShri Abhyankar Input Parameters: 242dc7a7e3SShri Abhyankar + ts - the TS context 252dc7a7e3SShri Abhyankar . t - current time 262dc7a7e3SShri Abhyankar . U - current iterate 272dc7a7e3SShri Abhyankar - ctx - [optional] context passed with eventmonitor 282dc7a7e3SShri Abhyankar 292dc7a7e3SShri Abhyankar Output parameters: 302dc7a7e3SShri Abhyankar + fvalue - function value of events at time t 312dc7a7e3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 322dc7a7e3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways 332dc7a7e3SShri Abhyankar - terminate - terminate time stepping after event is detected. 342dc7a7e3SShri Abhyankar 352dc7a7e3SShri Abhyankar Calling sequence of postevent: 362dc7a7e3SShri Abhyankar PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero, PetscReal t,Vec U,void* ctx) 372dc7a7e3SShri Abhyankar 382dc7a7e3SShri Abhyankar Input Parameters: 392dc7a7e3SShri Abhyankar + ts - the TS context 402dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 412dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 422dc7a7e3SShri Abhyankar . t - current time 432dc7a7e3SShri Abhyankar . U - current solution 442dc7a7e3SShri Abhyankar - ctx - the context passed with eventmonitor 452dc7a7e3SShri Abhyankar 462dc7a7e3SShri Abhyankar Level: intermediate 472dc7a7e3SShri Abhyankar 482dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor 492dc7a7e3SShri Abhyankar 502dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 512dc7a7e3SShri Abhyankar @*/ 522dc7a7e3SShri Abhyankar PetscErrorCode TSSetEventMonitor(TS ts,PetscInt nevents,PetscErrorCode (*eventmonitor)(TS,PetscReal,Vec,PetscScalar*,PetscInt*,PetscBool*,void*),PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,void*),void *mectx) 532dc7a7e3SShri Abhyankar { 542dc7a7e3SShri Abhyankar PetscErrorCode ierr; 552dc7a7e3SShri Abhyankar PetscReal t; 562dc7a7e3SShri Abhyankar Vec U; 572dc7a7e3SShri Abhyankar TSEvent event; 582dc7a7e3SShri Abhyankar 592dc7a7e3SShri Abhyankar PetscFunctionBegin; 602dc7a7e3SShri Abhyankar ierr = PetscNew(struct _p_TSEvent,&ts->event);CHKERRQ(ierr); 612dc7a7e3SShri Abhyankar event = ts->event; 622dc7a7e3SShri Abhyankar ierr = PetscMalloc(nevents*sizeof(PetscScalar),&event->fvalue);CHKERRQ(ierr); 632dc7a7e3SShri Abhyankar ierr = PetscMalloc(nevents*sizeof(PetscScalar),&event->fvalue_prev);CHKERRQ(ierr); 642dc7a7e3SShri Abhyankar ierr = PetscMalloc(nevents*sizeof(PetscBool),&event->terminate);CHKERRQ(ierr); 652dc7a7e3SShri Abhyankar ierr = PetscMalloc(nevents*sizeof(PetscInt),&event->direction);CHKERRQ(ierr); 662dc7a7e3SShri Abhyankar ierr = PetscMalloc(nevents*sizeof(PetscInt),&event->events_zero);CHKERRQ(ierr); 672dc7a7e3SShri Abhyankar event->monitor = eventmonitor; 682dc7a7e3SShri Abhyankar event->postevent = postevent; 692dc7a7e3SShri Abhyankar event->monitorcontext = (void*)mectx; 702dc7a7e3SShri Abhyankar event->nevents = nevents; 712dc7a7e3SShri Abhyankar 722dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 732dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&event->initial_timestep);CHKERRQ(ierr); 742dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 752dc7a7e3SShri Abhyankar event->ptime_prev = t; 762dc7a7e3SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,event->direction,event->terminate,mectx);CHKERRQ(ierr); 772dc7a7e3SShri Abhyankar ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"TS Event options","");CHKERRQ(ierr); 782dc7a7e3SShri Abhyankar { 792dc7a7e3SShri Abhyankar event->tol = 1.0e-6; 802dc7a7e3SShri Abhyankar ierr = PetscOptionsReal("-ts_event_tol","","",event->tol,&event->tol,NULL);CHKERRQ(ierr); 812dc7a7e3SShri Abhyankar } 822dc7a7e3SShri Abhyankar ierr = PetscOptionsEnd();CHKERRQ(ierr); 832dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 842dc7a7e3SShri Abhyankar } 852dc7a7e3SShri Abhyankar 862dc7a7e3SShri Abhyankar #undef __FUNCT__ 872dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent" 882dc7a7e3SShri Abhyankar PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,void *ctx) 892dc7a7e3SShri Abhyankar { 902dc7a7e3SShri Abhyankar PetscErrorCode ierr; 912dc7a7e3SShri Abhyankar TSEvent event=ts->event; 922dc7a7e3SShri Abhyankar PetscBool terminate=PETSC_FALSE; 932dc7a7e3SShri Abhyankar PetscInt i; 942dc7a7e3SShri Abhyankar PetscBool ts_terminate; 952dc7a7e3SShri Abhyankar 962dc7a7e3SShri Abhyankar PetscFunctionBegin; 972dc7a7e3SShri Abhyankar if (event->postevent) { 982dc7a7e3SShri Abhyankar ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,ctx);CHKERRQ(ierr); 992dc7a7e3SShri Abhyankar } 1002dc7a7e3SShri Abhyankar for(i = 0; i < nevents_zero;i++) { 101e105d053SSatish Balay terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]); 1022dc7a7e3SShri Abhyankar } 103*b4549700SJed Brown ierr = MPI_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 1042dc7a7e3SShri Abhyankar if (terminate) { 1052dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 1062dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 1072dc7a7e3SShri Abhyankar } else { 1082dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 1092dc7a7e3SShri Abhyankar } 1102dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1112dc7a7e3SShri Abhyankar } 1122dc7a7e3SShri Abhyankar 1132dc7a7e3SShri Abhyankar #undef __FUNCT__ 1142dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitorDestroy" 1152dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitorDestroy(TSEvent *event) 1162dc7a7e3SShri Abhyankar { 1172dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1182dc7a7e3SShri Abhyankar 1192dc7a7e3SShri Abhyankar PetscFunctionBegin; 1202dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 1212dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 1222dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 1232dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 1242dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 1252dc7a7e3SShri Abhyankar ierr = PetscFree(*event);CHKERRQ(ierr); 1262dc7a7e3SShri Abhyankar *event = NULL; 1272dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1282dc7a7e3SShri Abhyankar } 1292dc7a7e3SShri Abhyankar 1302dc7a7e3SShri Abhyankar #undef __FUNCT__ 1312dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitor" 1322dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitor(TS ts) 1332dc7a7e3SShri Abhyankar { 1342dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1352dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1362dc7a7e3SShri Abhyankar PetscReal t; 1372dc7a7e3SShri Abhyankar Vec U; 1382dc7a7e3SShri Abhyankar PetscInt i; 1392dc7a7e3SShri Abhyankar PetscReal dt; 140*b4549700SJed Brown TSEventStatus status = event->status; 1412dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 1422dc7a7e3SShri Abhyankar 1432dc7a7e3SShri Abhyankar PetscFunctionBegin; 1442dc7a7e3SShri Abhyankar 1452dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 1462dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 1472dc7a7e3SShri Abhyankar 1482dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 1492dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 1502dc7a7e3SShri Abhyankar /* Take initial time step */ 1512dc7a7e3SShri Abhyankar dt = event->initial_timestep; 1522dc7a7e3SShri Abhyankar ts->time_step = dt; 1532dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 1542dc7a7e3SShri Abhyankar } 1552dc7a7e3SShri Abhyankar 1562dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 1572dc7a7e3SShri Abhyankar event->tstepend = t; 1582dc7a7e3SShri Abhyankar } 1592dc7a7e3SShri Abhyankar 1602dc7a7e3SShri Abhyankar event->nevents_zero = 0; 1612dc7a7e3SShri Abhyankar 1622dc7a7e3SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue,event->direction,event->terminate,event->monitorcontext);CHKERRQ(ierr); 1632dc7a7e3SShri Abhyankar if (event->status != TSEVENT_NONE) { 1642dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 165e105d053SSatish Balay if (PetscAbsScalar(event->fvalue[i]) < event->tol) { 1662dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 1672dc7a7e3SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 1682dc7a7e3SShri Abhyankar } 1692dc7a7e3SShri Abhyankar } 1702dc7a7e3SShri Abhyankar } 1712dc7a7e3SShri Abhyankar 1722dc7a7e3SShri Abhyankar status = event->status; 173*b4549700SJed Brown ierr = MPI_Allreduce((PetscEnum*)&status,(PetscEnum*)&event->status,1,MPIU_ENUM,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 1742dc7a7e3SShri Abhyankar 1752dc7a7e3SShri Abhyankar if (event->status == TSEVENT_ZERO) { 1762dc7a7e3SShri Abhyankar dt = event->tstepend-t; 1772dc7a7e3SShri Abhyankar ts->time_step = dt; 1782dc7a7e3SShri Abhyankar ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,event->monitorcontext);CHKERRQ(ierr); 1792dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 1802dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 1812dc7a7e3SShri Abhyankar } 1822dc7a7e3SShri Abhyankar event->ptime_prev = t; 1832dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1842dc7a7e3SShri Abhyankar } 1852dc7a7e3SShri Abhyankar 1862dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 187e105d053SSatish Balay if ((event->direction[i] < 0 && PetscSign(PetscRealPart(event->fvalue[i])) <= 0 && PetscSign(PetscRealPart(event->fvalue_prev[i])) >= 0) || \ 188e105d053SSatish Balay (event->direction[i] > 0 && PetscSign(PetscRealPart(event->fvalue[i])) >= 0 && PetscSign(PetscRealPart(event->fvalue_prev[i])) <= 0) || \ 189e105d053SSatish Balay (event->direction[i] == 0 && PetscSign(PetscRealPart(event->fvalue[i]))*PetscSign(PetscRealPart(event->fvalue_prev[i])) <= 0)) { 1902dc7a7e3SShri Abhyankar 1912dc7a7e3SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 1922dc7a7e3SShri Abhyankar rollback = 1; 1932dc7a7e3SShri Abhyankar /* Compute linearly interpolated new time step */ 194e105d053SSatish Balay dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 1952dc7a7e3SShri Abhyankar } 1962dc7a7e3SShri Abhyankar } 1972dc7a7e3SShri Abhyankar in[0] = event->status; 1982dc7a7e3SShri Abhyankar in[1] = rollback; 1992dc7a7e3SShri Abhyankar ierr = MPI_Allreduce(in,out,2,MPIU_INT,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 2002dc7a7e3SShri Abhyankar 2012dc7a7e3SShri Abhyankar rollback = out[1]; 2022dc7a7e3SShri Abhyankar if (rollback) { 2032dc7a7e3SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 2042dc7a7e3SShri Abhyankar } 2052dc7a7e3SShri Abhyankar 2062dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 2072dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 2082dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 2092dc7a7e3SShri Abhyankar } else { 2102dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 2112dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 2122dc7a7e3SShri Abhyankar } 2132dc7a7e3SShri Abhyankar event->ptime_prev = t; 2142dc7a7e3SShri Abhyankar if (event->status == TSEVENT_PROCESSING) { 2152dc7a7e3SShri Abhyankar dt = event->tstepend - event->ptime_prev; 2162dc7a7e3SShri Abhyankar } 2172dc7a7e3SShri Abhyankar } 218e6b5be14SSatish Balay ierr = MPI_Allreduce(&dt,&(ts->time_step),1,MPIU_REAL,MPI_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 2192dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2202dc7a7e3SShri Abhyankar } 2212dc7a7e3SShri Abhyankar 222