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; 18*38bf2713SShri 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); 377dbe0728SLisandro Dalcin ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 387dbe0728SLisandro Dalcin ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 397dbe0728SLisandro Dalcin ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 407dbe0728SLisandro Dalcin ierr = PetscFree((*event)->vtol);CHKERRQ(ierr); 417dbe0728SLisandro Dalcin for (i=0; i < MAXEVENTRECORDERS; i++) { 427dbe0728SLisandro Dalcin ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr); 437dbe0728SLisandro Dalcin } 447dbe0728SLisandro Dalcin ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr); 457dbe0728SLisandro Dalcin ierr = PetscFree(*event);CHKERRQ(ierr); 467dbe0728SLisandro Dalcin PetscFunctionReturn(0); 477dbe0728SLisandro Dalcin } 487dbe0728SLisandro Dalcin 497dbe0728SLisandro Dalcin #undef __FUNCT__ 50e3005195SShri Abhyankar #define __FUNCT__ "TSSetEventTolerances" 51e3005195SShri Abhyankar /*@ 52e3005195SShri Abhyankar TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 53e3005195SShri Abhyankar 54e3005195SShri Abhyankar Logically Collective 55e3005195SShri Abhyankar 56e3005195SShri Abhyankar Input Arguments: 57e3005195SShri Abhyankar + ts - time integration context 58e3005195SShri Abhyankar . tol - scalar tolerance, PETSC_DECIDE to leave current value 59e3005195SShri Abhyankar - vtol - array of tolerances or NULL, used in preference to tol if present 60e3005195SShri Abhyankar 612589fa24SLisandro Dalcin Options Database Keys: 622589fa24SLisandro Dalcin . -ts_event_tol <tol> tolerance for event zero crossing 63e3005195SShri Abhyankar 64e3005195SShri Abhyankar Notes: 65f25fe674SLisandro Dalcin Must call TSSetEventHandler() before setting the tolerances. 66e3005195SShri Abhyankar 67e3005195SShri Abhyankar The size of vtol is equal to the number of events. 68e3005195SShri Abhyankar 69e3005195SShri Abhyankar Level: beginner 70e3005195SShri Abhyankar 71f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler() 72e3005195SShri Abhyankar @*/ 736427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[]) 74e3005195SShri Abhyankar { 756427ac75SLisandro Dalcin TSEvent event; 76e3005195SShri Abhyankar PetscInt i; 77e3005195SShri Abhyankar 78e3005195SShri Abhyankar PetscFunctionBegin; 796427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 806427ac75SLisandro Dalcin if (vtol) PetscValidRealPointer(vtol,3); 81f25fe674SLisandro Dalcin if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()"); 826427ac75SLisandro Dalcin 836427ac75SLisandro Dalcin event = ts->event; 84e3005195SShri Abhyankar if (vtol) { 85e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 86e3005195SShri Abhyankar } else { 87e3005195SShri Abhyankar if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 88e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 89e3005195SShri Abhyankar } 90e3005195SShri Abhyankar } 91e3005195SShri Abhyankar PetscFunctionReturn(0); 92e3005195SShri Abhyankar } 93e3005195SShri Abhyankar 94e3005195SShri Abhyankar #undef __FUNCT__ 95f25fe674SLisandro Dalcin #define __FUNCT__ "TSSetEventHandler" 962dc7a7e3SShri Abhyankar /*@C 97f25fe674SLisandro Dalcin TSSetEventHandler - Sets a monitoring function used for detecting events 982dc7a7e3SShri Abhyankar 992dc7a7e3SShri Abhyankar Logically Collective on TS 1002dc7a7e3SShri Abhyankar 1012dc7a7e3SShri Abhyankar Input Parameters: 1022dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 1032dc7a7e3SShri Abhyankar . nevents - number of local events 104d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 105d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 106d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 107d94325d3SShri Abhyankar event is detected (one for each event) 1086427ac75SLisandro Dalcin . eventhandler - event monitoring routine 1092dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 1102589fa24SLisandro Dalcin - ctx - [optional] user-defined context for private data for the 1112dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 1122dc7a7e3SShri Abhyankar context is desired) 1132dc7a7e3SShri Abhyankar 1146427ac75SLisandro Dalcin Calling sequence of eventhandler: 1152589fa24SLisandro Dalcin PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 1162dc7a7e3SShri Abhyankar 1172dc7a7e3SShri Abhyankar Input Parameters: 1182dc7a7e3SShri Abhyankar + ts - the TS context 1192dc7a7e3SShri Abhyankar . t - current time 1202dc7a7e3SShri Abhyankar . U - current iterate 1216427ac75SLisandro Dalcin - ctx - [optional] context passed with eventhandler 1222dc7a7e3SShri Abhyankar 1232dc7a7e3SShri Abhyankar Output parameters: 124d94325d3SShri Abhyankar . fvalue - function value of events at time t 1252dc7a7e3SShri Abhyankar 1262dc7a7e3SShri Abhyankar Calling sequence of postevent: 1272589fa24SLisandro Dalcin PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero[], PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 1282dc7a7e3SShri Abhyankar 1292dc7a7e3SShri Abhyankar Input Parameters: 1302dc7a7e3SShri Abhyankar + ts - the TS context 1312dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 1322dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 1332dc7a7e3SShri Abhyankar . t - current time 1342dc7a7e3SShri Abhyankar . U - current solution 135031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 1366427ac75SLisandro Dalcin - ctx - the context passed with eventhandler 1372dc7a7e3SShri Abhyankar 1382dc7a7e3SShri Abhyankar Level: intermediate 1392dc7a7e3SShri Abhyankar 1402dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor 1412dc7a7e3SShri Abhyankar 1422dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 1432dc7a7e3SShri Abhyankar @*/ 1442589fa24SLisandro 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) 1452dc7a7e3SShri Abhyankar { 1462dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1472dc7a7e3SShri Abhyankar TSEvent event; 148d94325d3SShri Abhyankar PetscInt i; 149006e6a18SShri Abhyankar PetscBool flg; 150e3005195SShri Abhyankar PetscReal tol=1e-6; 1512dc7a7e3SShri Abhyankar 1522dc7a7e3SShri Abhyankar PetscFunctionBegin; 1536427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1546427ac75SLisandro Dalcin PetscValidIntPointer(direction,2); 1556427ac75SLisandro Dalcin PetscValidIntPointer(terminate,3); 1566427ac75SLisandro Dalcin 1576427ac75SLisandro Dalcin ierr = PetscNewLog(ts,&event);CHKERRQ(ierr); 158854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 159854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 160854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 161854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 162e3005195SShri Abhyankar ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr); 163d94325d3SShri Abhyankar for (i=0; i < nevents; i++) { 164d94325d3SShri Abhyankar event->direction[i] = direction[i]; 165d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 166d94325d3SShri Abhyankar } 167854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 1682589fa24SLisandro Dalcin event->nevents = nevents; 1696427ac75SLisandro Dalcin event->eventhandler = eventhandler; 1702dc7a7e3SShri Abhyankar event->postevent = postevent; 1716427ac75SLisandro Dalcin event->ctx = ctx; 1722dc7a7e3SShri Abhyankar 173f7aea88cSShri Abhyankar for (i=0; i < MAXEVENTRECORDERS; i++) { 1742589fa24SLisandro Dalcin ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 175f7aea88cSShri Abhyankar } 176f7aea88cSShri Abhyankar 177a9514d71SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr); 1782dc7a7e3SShri Abhyankar { 179e3005195SShri Abhyankar ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr); 180006e6a18SShri Abhyankar ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr); 1812dc7a7e3SShri Abhyankar } 1829e12be75SShri Abhyankar PetscOptionsEnd(); 183e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 1846427ac75SLisandro Dalcin if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);} 1859e12be75SShri Abhyankar 1866427ac75SLisandro Dalcin ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr); 187d94325d3SShri Abhyankar ts->event = event; 1882dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1892dc7a7e3SShri Abhyankar } 1902dc7a7e3SShri Abhyankar 191031fbad4SShri Abhyankar /* 1924597913aSLisandro Dalcin Helper rutine to handle user postenvents and recording 193031fbad4SShri Abhyankar */ 194031fbad4SShri Abhyankar #undef __FUNCT__ 195031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent" 1964597913aSLisandro Dalcin static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U) 1972dc7a7e3SShri Abhyankar { 1982dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1992dc7a7e3SShri Abhyankar TSEvent event = ts->event; 2002dc7a7e3SShri Abhyankar PetscBool terminate = PETSC_FALSE; 201d0578d90SShri Abhyankar PetscInt i,ctr,stepnum; 2022dc7a7e3SShri Abhyankar PetscBool ts_terminate; 2034597913aSLisandro Dalcin PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 2042dc7a7e3SShri Abhyankar 2052dc7a7e3SShri Abhyankar PetscFunctionBegin; 2062dc7a7e3SShri Abhyankar if (event->postevent) { 2074597913aSLisandro Dalcin ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 2082dc7a7e3SShri Abhyankar } 2094597913aSLisandro Dalcin 2104597913aSLisandro Dalcin /* Handle termination events */ 2114597913aSLisandro Dalcin for (i=0; i < event->nevents_zero; i++) { 2124597913aSLisandro Dalcin terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]); 2132dc7a7e3SShri Abhyankar } 214b2566f29SBarry Smith ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 2159e12be75SShri Abhyankar if (ts_terminate) { 2162dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 2172dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2182dc7a7e3SShri Abhyankar } else { 2192dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 2202dc7a7e3SShri Abhyankar } 221f7aea88cSShri Abhyankar 2224597913aSLisandro Dalcin event->ptime_prev = t; 2234597913aSLisandro Dalcin /* Reset event residual functions as states might get changed by the postevent callback */ 2244597913aSLisandro Dalcin if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);} 2254597913aSLisandro Dalcin /* Cache current event residual functions */ 2264597913aSLisandro Dalcin for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 2274597913aSLisandro Dalcin 228d0578d90SShri Abhyankar /* Record the event in the event recorder */ 229d0578d90SShri Abhyankar ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr); 230f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 2316c4ed002SBarry Smith if (ctr == MAXEVENTRECORDERS) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS); 232f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 233d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 2344597913aSLisandro Dalcin event->recorder.nevents[ctr] = event->nevents_zero; 2354597913aSLisandro Dalcin for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 236f7aea88cSShri Abhyankar event->recorder.ctr++; 2372dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2382dc7a7e3SShri Abhyankar } 2392dc7a7e3SShri Abhyankar 240*38bf2713SShri Abhyankar /*PETSC_STATIC_INLINE PetscReal ComputeInterpolatedStep(TSEvent event,PetscInt i, PetscReal dt) 241*38bf2713SShri Abhyankar { 242*38bf2713SShri Abhyankar PetscReal newt; 243*38bf2713SShri Abhyankar 244*38bf2713SShri Abhyankar newt = 245*38bf2713SShri Abhyankar 246*38bf2713SShri Abhyankar return dt; 247*38bf2713SShri Abhyankar } 248*38bf2713SShri Abhyankar */ 249*38bf2713SShri Abhyankar 2502dc7a7e3SShri Abhyankar #undef __FUNCT__ 2516427ac75SLisandro Dalcin #define __FUNCT__ "TSEventHandler" 2526427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts) 2532dc7a7e3SShri Abhyankar { 2542dc7a7e3SShri Abhyankar PetscErrorCode ierr; 2556427ac75SLisandro Dalcin TSEvent event; 2562dc7a7e3SShri Abhyankar PetscReal t; 2572dc7a7e3SShri Abhyankar Vec U; 2582dc7a7e3SShri Abhyankar PetscInt i; 2597dbe0728SLisandro Dalcin PetscReal dt,dt_min; 2602dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 2619e12be75SShri Abhyankar PetscInt fvalue_sign,fvalueprev_sign; 2622dc7a7e3SShri Abhyankar 2632dc7a7e3SShri Abhyankar PetscFunctionBegin; 2646427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2656427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 2666427ac75SLisandro Dalcin event = ts->event; 2672dc7a7e3SShri Abhyankar 2682dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 2697dbe0728SLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 2702dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 2712dc7a7e3SShri Abhyankar 2727dbe0728SLisandro Dalcin if (event->status == TSEVENT_NONE) { 2737dbe0728SLisandro Dalcin event->timestep_prev = dt; 2747dbe0728SLisandro Dalcin } 2757dbe0728SLisandro Dalcin 2762dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 2777dbe0728SLisandro Dalcin /* Restore previous time step */ 2787dbe0728SLisandro Dalcin dt = event->timestep_prev; 2797dbe0728SLisandro Dalcin ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 2802dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2812dc7a7e3SShri Abhyankar } 2822dc7a7e3SShri Abhyankar 2832dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 2847dbe0728SLisandro Dalcin event->ptime_end = t; 2852dc7a7e3SShri Abhyankar } 2862dc7a7e3SShri Abhyankar 2876427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); 2889e12be75SShri Abhyankar 2892dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 290e3005195SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 2912dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 2929e12be75SShri Abhyankar continue; 293006e6a18SShri Abhyankar } 294d94325d3SShri Abhyankar fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 295d94325d3SShri Abhyankar fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 296e3005195SShri Abhyankar if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) { 297d94325d3SShri Abhyankar switch (event->direction[i]) { 298d94325d3SShri Abhyankar case -1: 299d94325d3SShri Abhyankar if (fvalue_sign < 0) { 3002dc7a7e3SShri Abhyankar rollback = 1; 3012dc7a7e3SShri Abhyankar /* Compute linearly interpolated new time step */ 302e105d053SSatish Balay dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3036427ac75SLisandro Dalcin if (event->monitor) { 304*38bf2713SShri 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); 305006e6a18SShri Abhyankar } 3062dc7a7e3SShri Abhyankar } 307d94325d3SShri Abhyankar break; 308d94325d3SShri Abhyankar case 1: 309d94325d3SShri Abhyankar if (fvalue_sign > 0) { 310d94325d3SShri Abhyankar rollback = 1; 311d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 312d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3136427ac75SLisandro Dalcin if (event->monitor) { 314*38bf2713SShri 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); 315006e6a18SShri Abhyankar } 3162dc7a7e3SShri Abhyankar } 317d94325d3SShri Abhyankar break; 318d94325d3SShri Abhyankar case 0: 319d94325d3SShri Abhyankar rollback = 1; 320d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 321d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3226427ac75SLisandro Dalcin if (event->monitor) { 323*38bf2713SShri 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); 324006e6a18SShri Abhyankar } 325d94325d3SShri Abhyankar break; 326d94325d3SShri Abhyankar } 327d94325d3SShri Abhyankar } 328d94325d3SShri Abhyankar } 3297dbe0728SLisandro Dalcin 330d94325d3SShri Abhyankar if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 3312589fa24SLisandro Dalcin in[0] = event->status; in[1] = rollback; 3327dbe0728SLisandro Dalcin ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 3332589fa24SLisandro Dalcin event->status = (TSEventStatus)out[0]; rollback = out[1]; 3342589fa24SLisandro Dalcin if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 3352dc7a7e3SShri Abhyankar 3367dbe0728SLisandro Dalcin event->nevents_zero = 0; 3379e12be75SShri Abhyankar if (event->status == TSEVENT_ZERO) { 3389e12be75SShri Abhyankar for (i=0; i < event->nevents; i++) { 3399e12be75SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 3409e12be75SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 3416427ac75SLisandro Dalcin if (event->monitor) { 342*38bf2713SShri 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); 3439e12be75SShri Abhyankar } 3449e12be75SShri Abhyankar } 3459e12be75SShri Abhyankar } 3464597913aSLisandro Dalcin ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr); 3479e12be75SShri Abhyankar 3487dbe0728SLisandro Dalcin dt = event->ptime_end - event->ptime_prev; 3497dbe0728SLisandro Dalcin if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */ 3509e12be75SShri Abhyankar ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 351*38bf2713SShri Abhyankar event->iterctr = 0; 3529e12be75SShri Abhyankar PetscFunctionReturn(0); 3539e12be75SShri Abhyankar } 3549e12be75SShri Abhyankar 3552dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 3562dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 357c540466cSShri Abhyankar ts->steps--; 358c540466cSShri Abhyankar ts->total_steps--; 3592589fa24SLisandro Dalcin ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr); 3602dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 3612dc7a7e3SShri Abhyankar } else { 3622dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 3632dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 3642dc7a7e3SShri Abhyankar } 365*38bf2713SShri Abhyankar 3662dc7a7e3SShri Abhyankar if (event->status == TSEVENT_PROCESSING) { 367*38bf2713SShri Abhyankar dt = event->ptime_end - t; 368*38bf2713SShri Abhyankar if (event->monitor) { 369*38bf2713SShri 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); 3702dc7a7e3SShri Abhyankar } 3712dc7a7e3SShri Abhyankar } 372*38bf2713SShri Abhyankar event->ptime_prev = t; 373*38bf2713SShri Abhyankar } 374*38bf2713SShri Abhyankar 375*38bf2713SShri Abhyankar if(event->status == TSEVENT_PROCESSING) event->iterctr++; 3769e12be75SShri Abhyankar 3777dbe0728SLisandro Dalcin ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 3787dbe0728SLisandro Dalcin ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr); 3792dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3802dc7a7e3SShri Abhyankar } 3812dc7a7e3SShri Abhyankar 382d0578d90SShri Abhyankar #undef __FUNCT__ 3836427ac75SLisandro Dalcin #define __FUNCT__ "TSAdjointEventHandler" 3846427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts) 385d0578d90SShri Abhyankar { 386d0578d90SShri Abhyankar PetscErrorCode ierr; 3876427ac75SLisandro Dalcin TSEvent event; 388d0578d90SShri Abhyankar PetscReal t; 389d0578d90SShri Abhyankar Vec U; 390d0578d90SShri Abhyankar PetscInt ctr; 391d0578d90SShri Abhyankar PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 392d0578d90SShri Abhyankar 393d0578d90SShri Abhyankar PetscFunctionBegin; 3946427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3956427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 3966427ac75SLisandro Dalcin event = ts->event; 397d0578d90SShri Abhyankar 398d0578d90SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 399d0578d90SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 400d0578d90SShri Abhyankar 401d0578d90SShri Abhyankar ctr = event->recorder.ctr-1; 402bcbf8bb3SShri Abhyankar if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 403d0578d90SShri Abhyankar /* Call the user postevent function */ 404d0578d90SShri Abhyankar if (event->postevent) { 4056427ac75SLisandro Dalcin ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 406d0578d90SShri Abhyankar event->recorder.ctr--; 407d0578d90SShri Abhyankar } 408d0578d90SShri Abhyankar } 409d0578d90SShri Abhyankar 410d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 411d0578d90SShri Abhyankar PetscFunctionReturn(0); 412d0578d90SShri Abhyankar } 413