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 78*f7aea88cSShri Abhyankar /* Initialize the event recorder */ 79*f7aea88cSShri Abhyankar event->recorder.ctr = 0; 80*f7aea88cSShri Abhyankar for(i=0; i < MAXEVENTRECORDERS; i++) { 81*f7aea88cSShri Abhyankar ierr = PetscMalloc1(nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 82*f7aea88cSShri Abhyankar } 83*f7aea88cSShri Abhyankar 842dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 852dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&event->initial_timestep);CHKERRQ(ierr); 862dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 872dc7a7e3SShri Abhyankar event->ptime_prev = t; 88d94325d3SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,mectx);CHKERRQ(ierr); 892dc7a7e3SShri Abhyankar ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"TS Event options","");CHKERRQ(ierr); 902dc7a7e3SShri Abhyankar { 912dc7a7e3SShri Abhyankar event->tol = 1.0e-6; 922dc7a7e3SShri Abhyankar ierr = PetscOptionsReal("-ts_event_tol","","",event->tol,&event->tol,NULL);CHKERRQ(ierr); 932dc7a7e3SShri Abhyankar } 942dc7a7e3SShri Abhyankar ierr = PetscOptionsEnd();CHKERRQ(ierr); 95d94325d3SShri Abhyankar ts->event = event; 962dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 972dc7a7e3SShri Abhyankar } 982dc7a7e3SShri Abhyankar 992dc7a7e3SShri Abhyankar #undef __FUNCT__ 1002dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent" 101031fbad4SShri Abhyankar /* 102031fbad4SShri Abhyankar TSPostEvent - Does post event processing by calling the user-defined postevent function 103031fbad4SShri Abhyankar 104031fbad4SShri Abhyankar Logically Collective on TS 105031fbad4SShri Abhyankar 106031fbad4SShri Abhyankar Input Parameters: 107031fbad4SShri Abhyankar + ts - the TS context 108031fbad4SShri Abhyankar . nevents_zero - number of local events whose event function is zero 109031fbad4SShri Abhyankar . events_zero - indices of local events which have reached zero 110031fbad4SShri Abhyankar . t - current time 111031fbad4SShri Abhyankar . U - current solution 112031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 113031fbad4SShri Abhyankar - ctx - the context passed with eventmonitor 114031fbad4SShri Abhyankar 115031fbad4SShri Abhyankar Level: intermediate 116031fbad4SShri Abhyankar 117031fbad4SShri Abhyankar .keywords: TS, event, set, monitor 118031fbad4SShri Abhyankar 119031fbad4SShri Abhyankar .seealso: TSSetEventMonitor(),TSEvent 120031fbad4SShri Abhyankar */ 121031fbad4SShri Abhyankar #undef __FUNCT__ 122031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent" 123031fbad4SShri Abhyankar PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void *ctx) 1242dc7a7e3SShri Abhyankar { 1252dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1262dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1272dc7a7e3SShri Abhyankar PetscBool terminate=PETSC_FALSE; 128*f7aea88cSShri Abhyankar PetscInt i,ctr; 1292dc7a7e3SShri Abhyankar PetscBool ts_terminate; 1302dc7a7e3SShri Abhyankar 1312dc7a7e3SShri Abhyankar PetscFunctionBegin; 1322dc7a7e3SShri Abhyankar if (event->postevent) { 133031fbad4SShri Abhyankar ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,forwardsolve,ctx);CHKERRQ(ierr); 1342dc7a7e3SShri Abhyankar } 1352dc7a7e3SShri Abhyankar for(i = 0; i < nevents_zero;i++) { 136e105d053SSatish Balay terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]); 1372dc7a7e3SShri Abhyankar } 138b4549700SJed Brown ierr = MPI_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 1392dc7a7e3SShri Abhyankar if (terminate) { 1402dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 1412dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 1422dc7a7e3SShri Abhyankar } else { 1432dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 1442dc7a7e3SShri Abhyankar } 145*f7aea88cSShri Abhyankar 146*f7aea88cSShri Abhyankar /* Record the events in the event recorder */ 147*f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 148*f7aea88cSShri Abhyankar if (ctr == MAXEVENTRECORDERS) { 149*f7aea88cSShri Abhyankar SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS); 150*f7aea88cSShri Abhyankar } 151*f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 152*f7aea88cSShri Abhyankar event->recorder.nevents[ctr] = nevents_zero; 153*f7aea88cSShri Abhyankar for(i=0; i < nevents_zero; i++) event->recorder.eventidx[ctr][i] = events_zero[i]; 154*f7aea88cSShri Abhyankar event->recorder.ctr++; 155*f7aea88cSShri Abhyankar 15673967516SShri Abhyankar /* Reset the event residual functions as states might get changed by the postevent callback */ 15773967516SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,event->monitorcontext);CHKERRQ(ierr); 15873967516SShri Abhyankar event->ptime_prev = t; 15973967516SShri Abhyankar 1602dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1612dc7a7e3SShri Abhyankar } 1622dc7a7e3SShri Abhyankar 1632dc7a7e3SShri Abhyankar #undef __FUNCT__ 1642dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitorDestroy" 1652dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitorDestroy(TSEvent *event) 1662dc7a7e3SShri Abhyankar { 1672dc7a7e3SShri Abhyankar PetscErrorCode ierr; 168*f7aea88cSShri Abhyankar PetscInt i; 1692dc7a7e3SShri Abhyankar 1702dc7a7e3SShri Abhyankar PetscFunctionBegin; 1712dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 1722dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 1732dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 174d94325d3SShri Abhyankar ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 1752dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 176*f7aea88cSShri Abhyankar for(i=0; i < MAXEVENTRECORDERS; i++) { 177*f7aea88cSShri Abhyankar ierr = PetscFree((*event)->recorder.eventidx[i]); 178*f7aea88cSShri Abhyankar } 1792dc7a7e3SShri Abhyankar ierr = PetscFree(*event);CHKERRQ(ierr); 1802dc7a7e3SShri Abhyankar *event = NULL; 1812dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1822dc7a7e3SShri Abhyankar } 1832dc7a7e3SShri Abhyankar 1842dc7a7e3SShri Abhyankar #undef __FUNCT__ 1852dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitor" 1862dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitor(TS ts) 1872dc7a7e3SShri Abhyankar { 1882dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1892dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1902dc7a7e3SShri Abhyankar PetscReal t; 1912dc7a7e3SShri Abhyankar Vec U; 1922dc7a7e3SShri Abhyankar PetscInt i; 1932dc7a7e3SShri Abhyankar PetscReal dt; 194b4549700SJed Brown TSEventStatus status = event->status; 1952dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 196031fbad4SShri Abhyankar PetscBool forwardsolve=PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 1972dc7a7e3SShri Abhyankar 1982dc7a7e3SShri Abhyankar PetscFunctionBegin; 1992dc7a7e3SShri Abhyankar 2002dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 2012dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 2022dc7a7e3SShri Abhyankar 2032dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 2042dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 2052dc7a7e3SShri Abhyankar /* Take initial time step */ 2062dc7a7e3SShri Abhyankar dt = event->initial_timestep; 2072dc7a7e3SShri Abhyankar ts->time_step = dt; 2082dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2092dc7a7e3SShri Abhyankar } 2102dc7a7e3SShri Abhyankar 2112dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 2122dc7a7e3SShri Abhyankar event->tstepend = t; 2132dc7a7e3SShri Abhyankar } 2142dc7a7e3SShri Abhyankar 2152dc7a7e3SShri Abhyankar event->nevents_zero = 0; 2162dc7a7e3SShri Abhyankar 217d94325d3SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue,event->monitorcontext);CHKERRQ(ierr); 2182dc7a7e3SShri Abhyankar if (event->status != TSEVENT_NONE) { 2192dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 220e105d053SSatish Balay if (PetscAbsScalar(event->fvalue[i]) < event->tol) { 2212dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 2222dc7a7e3SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 2232dc7a7e3SShri Abhyankar } 2242dc7a7e3SShri Abhyankar } 2252dc7a7e3SShri Abhyankar } 2262dc7a7e3SShri Abhyankar 2272dc7a7e3SShri Abhyankar status = event->status; 228b4549700SJed Brown ierr = MPI_Allreduce((PetscEnum*)&status,(PetscEnum*)&event->status,1,MPIU_ENUM,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 2292dc7a7e3SShri Abhyankar 2302dc7a7e3SShri Abhyankar if (event->status == TSEVENT_ZERO) { 231031fbad4SShri Abhyankar ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->monitorcontext);CHKERRQ(ierr); 232c540466cSShri Abhyankar dt = event->tstepend-t; 233c540466cSShri Abhyankar if(PetscAbsScalar(dt) < PETSC_SMALL) dt += event->initial_timestep; 234c540466cSShri Abhyankar ts->time_step = dt; 2352dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2362dc7a7e3SShri Abhyankar } 2372dc7a7e3SShri Abhyankar 2382dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 239d94325d3SShri Abhyankar PetscInt fvalue_sign,fvalueprev_sign; 240d94325d3SShri Abhyankar fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 241d94325d3SShri Abhyankar fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 2427a728e9fSShri Abhyankar if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) { 243d94325d3SShri Abhyankar switch (event->direction[i]) { 244d94325d3SShri Abhyankar case -1: 245d94325d3SShri Abhyankar if (fvalue_sign < 0) { 2462dc7a7e3SShri Abhyankar rollback = 1; 2472dc7a7e3SShri Abhyankar /* Compute linearly interpolated new time step */ 248e105d053SSatish Balay dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 2492dc7a7e3SShri Abhyankar } 250d94325d3SShri Abhyankar break; 251d94325d3SShri Abhyankar case 1: 252d94325d3SShri Abhyankar if (fvalue_sign > 0) { 253d94325d3SShri Abhyankar rollback = 1; 254d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 255d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 2562dc7a7e3SShri Abhyankar } 257d94325d3SShri Abhyankar break; 258d94325d3SShri Abhyankar case 0: 259d94325d3SShri Abhyankar rollback = 1; 260d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 261d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 262d94325d3SShri Abhyankar break; 263d94325d3SShri Abhyankar } 264d94325d3SShri Abhyankar } 265d94325d3SShri Abhyankar } 266d94325d3SShri Abhyankar if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 267d94325d3SShri Abhyankar 2682dc7a7e3SShri Abhyankar in[0] = event->status; 2692dc7a7e3SShri Abhyankar in[1] = rollback; 2702dc7a7e3SShri Abhyankar ierr = MPI_Allreduce(in,out,2,MPIU_INT,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 2712dc7a7e3SShri Abhyankar 2722dc7a7e3SShri Abhyankar rollback = out[1]; 2732dc7a7e3SShri Abhyankar if (rollback) { 2742dc7a7e3SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 2752dc7a7e3SShri Abhyankar } 2762dc7a7e3SShri Abhyankar 2772dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 2782dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 279c540466cSShri Abhyankar ts->steps--; 280c540466cSShri Abhyankar ts->total_steps--; 2812dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 2822dc7a7e3SShri Abhyankar } else { 2832dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 2842dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 2852dc7a7e3SShri Abhyankar } 2862dc7a7e3SShri Abhyankar event->ptime_prev = t; 2872dc7a7e3SShri Abhyankar if (event->status == TSEVENT_PROCESSING) { 2882dc7a7e3SShri Abhyankar dt = event->tstepend - event->ptime_prev; 2892dc7a7e3SShri Abhyankar } 2902dc7a7e3SShri Abhyankar } 291e6b5be14SSatish Balay ierr = MPI_Allreduce(&dt,&(ts->time_step),1,MPIU_REAL,MPI_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 2922dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2932dc7a7e3SShri Abhyankar } 2942dc7a7e3SShri Abhyankar 295