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 186fea3669SShri Abhyankar ierr = TSGetTimeStep(ts,&event->initial_timestep);CHKERRQ(ierr); 196fea3669SShri Abhyankar event->ptime_prev = t; 206427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr); 216fea3669SShri Abhyankar /* Initialize the event recorder */ 226fea3669SShri Abhyankar event->recorder.ctr = 0; 236fea3669SShri Abhyankar PetscFunctionReturn(0); 246fea3669SShri Abhyankar } 256fea3669SShri Abhyankar 266fea3669SShri Abhyankar #undef __FUNCT__ 27e3005195SShri Abhyankar #define __FUNCT__ "TSSetEventTolerances" 28e3005195SShri Abhyankar /*@ 29e3005195SShri Abhyankar TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 30e3005195SShri Abhyankar 31e3005195SShri Abhyankar Logically Collective 32e3005195SShri Abhyankar 33e3005195SShri Abhyankar Input Arguments: 34e3005195SShri Abhyankar + ts - time integration context 35e3005195SShri Abhyankar . tol - scalar tolerance, PETSC_DECIDE to leave current value 36e3005195SShri Abhyankar - vtol - array of tolerances or NULL, used in preference to tol if present 37e3005195SShri Abhyankar 38*2589fa24SLisandro Dalcin Options Database Keys: 39*2589fa24SLisandro Dalcin . -ts_event_tol <tol> tolerance for event zero crossing 40e3005195SShri Abhyankar 41e3005195SShri Abhyankar Notes: 42f25fe674SLisandro Dalcin Must call TSSetEventHandler() before setting the tolerances. 43e3005195SShri Abhyankar 44e3005195SShri Abhyankar The size of vtol is equal to the number of events. 45e3005195SShri Abhyankar 46e3005195SShri Abhyankar Level: beginner 47e3005195SShri Abhyankar 48f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler() 49e3005195SShri Abhyankar @*/ 506427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[]) 51e3005195SShri Abhyankar { 526427ac75SLisandro Dalcin TSEvent event; 53e3005195SShri Abhyankar PetscInt i; 54e3005195SShri Abhyankar 55e3005195SShri Abhyankar PetscFunctionBegin; 566427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 576427ac75SLisandro Dalcin if (vtol) PetscValidRealPointer(vtol,3); 58f25fe674SLisandro Dalcin if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()"); 596427ac75SLisandro Dalcin 606427ac75SLisandro Dalcin event = ts->event; 61e3005195SShri Abhyankar if (vtol) { 62e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 63e3005195SShri Abhyankar } else { 64e3005195SShri Abhyankar if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 65e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 66e3005195SShri Abhyankar } 67e3005195SShri Abhyankar } 68e3005195SShri Abhyankar PetscFunctionReturn(0); 69e3005195SShri Abhyankar } 70e3005195SShri Abhyankar 71e3005195SShri Abhyankar #undef __FUNCT__ 72f25fe674SLisandro Dalcin #define __FUNCT__ "TSSetEventHandler" 732dc7a7e3SShri Abhyankar /*@C 74f25fe674SLisandro Dalcin TSSetEventHandler - Sets a monitoring function used for detecting events 752dc7a7e3SShri Abhyankar 762dc7a7e3SShri Abhyankar Logically Collective on TS 772dc7a7e3SShri Abhyankar 782dc7a7e3SShri Abhyankar Input Parameters: 792dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 802dc7a7e3SShri Abhyankar . nevents - number of local events 81d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 82d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 83d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 84d94325d3SShri Abhyankar event is detected (one for each event) 856427ac75SLisandro Dalcin . eventhandler - event monitoring routine 862dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 87*2589fa24SLisandro Dalcin - ctx - [optional] user-defined context for private data for the 882dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 892dc7a7e3SShri Abhyankar context is desired) 902dc7a7e3SShri Abhyankar 916427ac75SLisandro Dalcin Calling sequence of eventhandler: 92*2589fa24SLisandro Dalcin PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 932dc7a7e3SShri Abhyankar 942dc7a7e3SShri Abhyankar Input Parameters: 952dc7a7e3SShri Abhyankar + ts - the TS context 962dc7a7e3SShri Abhyankar . t - current time 972dc7a7e3SShri Abhyankar . U - current iterate 986427ac75SLisandro Dalcin - ctx - [optional] context passed with eventhandler 992dc7a7e3SShri Abhyankar 1002dc7a7e3SShri Abhyankar Output parameters: 101d94325d3SShri Abhyankar . fvalue - function value of events at time t 1022dc7a7e3SShri Abhyankar 1032dc7a7e3SShri Abhyankar Calling sequence of postevent: 104*2589fa24SLisandro Dalcin PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero[], PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 1052dc7a7e3SShri Abhyankar 1062dc7a7e3SShri Abhyankar Input Parameters: 1072dc7a7e3SShri Abhyankar + ts - the TS context 1082dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 1092dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 1102dc7a7e3SShri Abhyankar . t - current time 1112dc7a7e3SShri Abhyankar . U - current solution 112031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 1136427ac75SLisandro Dalcin - ctx - the context passed with eventhandler 1142dc7a7e3SShri Abhyankar 1152dc7a7e3SShri Abhyankar Level: intermediate 1162dc7a7e3SShri Abhyankar 1172dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor 1182dc7a7e3SShri Abhyankar 1192dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 1202dc7a7e3SShri Abhyankar @*/ 121*2589fa24SLisandro 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) 1222dc7a7e3SShri Abhyankar { 1232dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1242dc7a7e3SShri Abhyankar TSEvent event; 125d94325d3SShri Abhyankar PetscInt i; 126006e6a18SShri Abhyankar PetscBool flg; 127e3005195SShri Abhyankar PetscReal tol=1e-6; 1282dc7a7e3SShri Abhyankar 1292dc7a7e3SShri Abhyankar PetscFunctionBegin; 1306427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1316427ac75SLisandro Dalcin PetscValidIntPointer(direction,2); 1326427ac75SLisandro Dalcin PetscValidIntPointer(terminate,3); 1336427ac75SLisandro Dalcin 1346427ac75SLisandro Dalcin ierr = PetscNewLog(ts,&event);CHKERRQ(ierr); 135854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 136854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 137854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 138854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 139e3005195SShri Abhyankar ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr); 140d94325d3SShri Abhyankar for (i=0; i < nevents; i++) { 141d94325d3SShri Abhyankar event->direction[i] = direction[i]; 142d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 143d94325d3SShri Abhyankar } 144854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 145*2589fa24SLisandro Dalcin event->nevents = nevents; 1466427ac75SLisandro Dalcin event->eventhandler = eventhandler; 1472dc7a7e3SShri Abhyankar event->postevent = postevent; 1486427ac75SLisandro Dalcin event->ctx = ctx; 1492dc7a7e3SShri Abhyankar 150f7aea88cSShri Abhyankar for (i=0; i < MAXEVENTRECORDERS; i++) { 151*2589fa24SLisandro Dalcin ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 152f7aea88cSShri Abhyankar } 153f7aea88cSShri Abhyankar 154a9514d71SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr); 1552dc7a7e3SShri Abhyankar { 156e3005195SShri Abhyankar ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr); 157006e6a18SShri Abhyankar ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr); 1582dc7a7e3SShri Abhyankar } 1599e12be75SShri Abhyankar PetscOptionsEnd(); 160e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 1616427ac75SLisandro Dalcin if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);} 1629e12be75SShri Abhyankar 1636427ac75SLisandro Dalcin ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr); 164d94325d3SShri Abhyankar ts->event = event; 1652dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1662dc7a7e3SShri Abhyankar } 1672dc7a7e3SShri Abhyankar 1682dc7a7e3SShri Abhyankar #undef __FUNCT__ 1692dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent" 170031fbad4SShri Abhyankar /* 171031fbad4SShri Abhyankar TSPostEvent - Does post event processing by calling the user-defined postevent function 172031fbad4SShri Abhyankar 173031fbad4SShri Abhyankar Logically Collective on TS 174031fbad4SShri Abhyankar 175031fbad4SShri Abhyankar Input Parameters: 176031fbad4SShri Abhyankar + ts - the TS context 177031fbad4SShri Abhyankar . nevents_zero - number of local events whose event function is zero 178031fbad4SShri Abhyankar . events_zero - indices of local events which have reached zero 179031fbad4SShri Abhyankar . t - current time 180031fbad4SShri Abhyankar . U - current solution 1816427ac75SLisandro Dalcin - forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 182031fbad4SShri Abhyankar 183031fbad4SShri Abhyankar Level: intermediate 184031fbad4SShri Abhyankar 185031fbad4SShri Abhyankar .keywords: TS, event, set, monitor 186031fbad4SShri Abhyankar 187f25fe674SLisandro Dalcin .seealso: TSSetEventHandler(), TSEvent 188031fbad4SShri Abhyankar */ 189031fbad4SShri Abhyankar #undef __FUNCT__ 190031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent" 1916427ac75SLisandro Dalcin PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve) 1922dc7a7e3SShri Abhyankar { 1932dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1942dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1952dc7a7e3SShri Abhyankar PetscBool terminate=PETSC_FALSE; 196d0578d90SShri Abhyankar PetscInt i,ctr,stepnum; 1972dc7a7e3SShri Abhyankar PetscBool ts_terminate; 1982dc7a7e3SShri Abhyankar 1992dc7a7e3SShri Abhyankar PetscFunctionBegin; 2002dc7a7e3SShri Abhyankar if (event->postevent) { 2016427ac75SLisandro Dalcin ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 2022dc7a7e3SShri Abhyankar } 2032dc7a7e3SShri Abhyankar for (i=0; i < nevents_zero;i++) { 204e105d053SSatish Balay terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]); 2052dc7a7e3SShri Abhyankar } 206b2566f29SBarry Smith ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 2079e12be75SShri Abhyankar if (ts_terminate) { 2082dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 2092dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2102dc7a7e3SShri Abhyankar } else { 2112dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 2122dc7a7e3SShri Abhyankar } 213f7aea88cSShri Abhyankar 214d0578d90SShri Abhyankar /* Record the event in the event recorder */ 215d0578d90SShri Abhyankar ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr); 216f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 2176c4ed002SBarry Smith if (ctr == MAXEVENTRECORDERS) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS); 218f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 219d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 220f7aea88cSShri Abhyankar event->recorder.nevents[ctr] = nevents_zero; 221f7aea88cSShri Abhyankar for (i=0; i < nevents_zero; i++) event->recorder.eventidx[ctr][i] = events_zero[i]; 222f7aea88cSShri Abhyankar event->recorder.ctr++; 223f7aea88cSShri Abhyankar 22473967516SShri Abhyankar /* Reset the event residual functions as states might get changed by the postevent callback */ 2256427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr); 22673967516SShri Abhyankar event->ptime_prev = t; 2272dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2282dc7a7e3SShri Abhyankar } 2292dc7a7e3SShri Abhyankar 2302dc7a7e3SShri Abhyankar #undef __FUNCT__ 2316427ac75SLisandro Dalcin #define __FUNCT__ "TSEventDestroy" 2326427ac75SLisandro Dalcin PetscErrorCode TSEventDestroy(TSEvent *event) 2332dc7a7e3SShri Abhyankar { 2342dc7a7e3SShri Abhyankar PetscErrorCode ierr; 235f7aea88cSShri Abhyankar PetscInt i; 2362dc7a7e3SShri Abhyankar 2372dc7a7e3SShri Abhyankar PetscFunctionBegin; 2386427ac75SLisandro Dalcin PetscValidPointer(event,1); 2396427ac75SLisandro Dalcin if (!*event) PetscFunctionReturn(0); 2402dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 2412dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 2422dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 243d94325d3SShri Abhyankar ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 2442dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 245e3005195SShri Abhyankar ierr = PetscFree((*event)->vtol);CHKERRQ(ierr); 246f7aea88cSShri Abhyankar for (i=0; i < MAXEVENTRECORDERS; i++) { 247302440fdSBarry Smith ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr); 248f7aea88cSShri Abhyankar } 2496427ac75SLisandro Dalcin ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr); 2502dc7a7e3SShri Abhyankar ierr = PetscFree(*event);CHKERRQ(ierr); 2512dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2522dc7a7e3SShri Abhyankar } 2532dc7a7e3SShri Abhyankar 2542dc7a7e3SShri Abhyankar #undef __FUNCT__ 2556427ac75SLisandro Dalcin #define __FUNCT__ "TSEventHandler" 2566427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts) 2572dc7a7e3SShri Abhyankar { 2582dc7a7e3SShri Abhyankar PetscErrorCode ierr; 2596427ac75SLisandro Dalcin TSEvent event; 2602dc7a7e3SShri Abhyankar PetscReal t; 2612dc7a7e3SShri Abhyankar Vec U; 2622dc7a7e3SShri Abhyankar PetscInt i; 2632dc7a7e3SShri Abhyankar PetscReal dt; 2642dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 265031fbad4SShri Abhyankar PetscBool forwardsolve=PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 2669e12be75SShri Abhyankar PetscInt fvalue_sign,fvalueprev_sign; 2672dc7a7e3SShri Abhyankar 2682dc7a7e3SShri Abhyankar PetscFunctionBegin; 2696427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2706427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 2716427ac75SLisandro Dalcin event = ts->event; 2722dc7a7e3SShri Abhyankar 2732dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 2742dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 2752dc7a7e3SShri Abhyankar 2762dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 2772dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 2782dc7a7e3SShri Abhyankar /* Take initial time step */ 2792dc7a7e3SShri Abhyankar dt = event->initial_timestep; 2802dc7a7e3SShri Abhyankar ts->time_step = dt; 2812dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2822dc7a7e3SShri Abhyankar } 2832dc7a7e3SShri Abhyankar 2842dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 2852dc7a7e3SShri Abhyankar event->tstepend = t; 2862dc7a7e3SShri Abhyankar } 2872dc7a7e3SShri Abhyankar 2882dc7a7e3SShri Abhyankar event->nevents_zero = 0; 2892dc7a7e3SShri Abhyankar 2906427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); 2919e12be75SShri Abhyankar 2922dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 293e3005195SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 2942dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 2959e12be75SShri Abhyankar continue; 296006e6a18SShri Abhyankar } 297d94325d3SShri Abhyankar fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 298d94325d3SShri Abhyankar fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 299e3005195SShri Abhyankar if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) { 300d94325d3SShri Abhyankar switch (event->direction[i]) { 301d94325d3SShri Abhyankar case -1: 302d94325d3SShri Abhyankar if (fvalue_sign < 0) { 3032dc7a7e3SShri Abhyankar rollback = 1; 3042dc7a7e3SShri Abhyankar /* Compute linearly interpolated new time step */ 305e105d053SSatish Balay dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3066427ac75SLisandro Dalcin if (event->monitor) { 3076427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 308006e6a18SShri Abhyankar } 3092dc7a7e3SShri Abhyankar } 310d94325d3SShri Abhyankar break; 311d94325d3SShri Abhyankar case 1: 312d94325d3SShri Abhyankar if (fvalue_sign > 0) { 313d94325d3SShri Abhyankar rollback = 1; 314d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 315d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3166427ac75SLisandro Dalcin if (event->monitor) { 3176427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 318006e6a18SShri Abhyankar } 3192dc7a7e3SShri Abhyankar } 320d94325d3SShri Abhyankar break; 321d94325d3SShri Abhyankar case 0: 322d94325d3SShri Abhyankar rollback = 1; 323d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 324d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3256427ac75SLisandro Dalcin if (event->monitor) { 3266427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 327006e6a18SShri Abhyankar } 328d94325d3SShri Abhyankar break; 329d94325d3SShri Abhyankar } 330d94325d3SShri Abhyankar } 331d94325d3SShri Abhyankar } 332d94325d3SShri Abhyankar if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 333d94325d3SShri Abhyankar 334*2589fa24SLisandro Dalcin in[0] = event->status; in[1] = rollback; 335b2566f29SBarry Smith ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 3362dc7a7e3SShri Abhyankar 337*2589fa24SLisandro Dalcin event->status = (TSEventStatus)out[0]; rollback = out[1]; 338*2589fa24SLisandro Dalcin if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 3392dc7a7e3SShri Abhyankar 3409e12be75SShri Abhyankar if (event->status == TSEVENT_ZERO) { 3419e12be75SShri Abhyankar for (i=0; i < event->nevents; i++) { 3429e12be75SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 3439e12be75SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 3446427ac75SLisandro Dalcin if (event->monitor) { 3456427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g\n",i,(double)t);CHKERRQ(ierr); 3469e12be75SShri Abhyankar } 3479e12be75SShri Abhyankar } 3489e12be75SShri Abhyankar } 3499e12be75SShri Abhyankar 3506427ac75SLisandro Dalcin ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve);CHKERRQ(ierr); 3519e12be75SShri Abhyankar 3529e12be75SShri Abhyankar dt = event->tstepend - t; 3539e12be75SShri Abhyankar if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->initial_timestep; 3549e12be75SShri Abhyankar ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 3559e12be75SShri Abhyankar PetscFunctionReturn(0); 3569e12be75SShri Abhyankar } 3579e12be75SShri Abhyankar 3582dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 3592dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 360c540466cSShri Abhyankar ts->steps--; 361c540466cSShri Abhyankar ts->total_steps--; 362*2589fa24SLisandro Dalcin ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr); 3632dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 3642dc7a7e3SShri Abhyankar } else { 3652dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 3662dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 3672dc7a7e3SShri Abhyankar } 3682dc7a7e3SShri Abhyankar event->ptime_prev = t; 3692dc7a7e3SShri Abhyankar if (event->status == TSEVENT_PROCESSING) { 3702dc7a7e3SShri Abhyankar dt = event->tstepend - event->ptime_prev; 3712dc7a7e3SShri Abhyankar } 3722dc7a7e3SShri Abhyankar } 3739e12be75SShri Abhyankar 374b2566f29SBarry Smith ierr = MPIU_Allreduce(&dt,&(ts->time_step),1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 3752dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3762dc7a7e3SShri Abhyankar } 3772dc7a7e3SShri Abhyankar 378d0578d90SShri Abhyankar #undef __FUNCT__ 3796427ac75SLisandro Dalcin #define __FUNCT__ "TSAdjointEventHandler" 3806427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts) 381d0578d90SShri Abhyankar { 382d0578d90SShri Abhyankar PetscErrorCode ierr; 3836427ac75SLisandro Dalcin TSEvent event; 384d0578d90SShri Abhyankar PetscReal t; 385d0578d90SShri Abhyankar Vec U; 386d0578d90SShri Abhyankar PetscInt ctr; 387d0578d90SShri Abhyankar PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 388d0578d90SShri Abhyankar 389d0578d90SShri Abhyankar PetscFunctionBegin; 3906427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3916427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 3926427ac75SLisandro Dalcin event = ts->event; 393d0578d90SShri Abhyankar 394d0578d90SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 395d0578d90SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 396d0578d90SShri Abhyankar 397d0578d90SShri Abhyankar ctr = event->recorder.ctr-1; 398bcbf8bb3SShri Abhyankar if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 399d0578d90SShri Abhyankar /* Call the user postevent function */ 400d0578d90SShri Abhyankar if (event->postevent) { 4016427ac75SLisandro Dalcin ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 402d0578d90SShri Abhyankar event->recorder.ctr--; 403d0578d90SShri Abhyankar } 404d0578d90SShri Abhyankar } 405d0578d90SShri Abhyankar 406d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 407d0578d90SShri Abhyankar PetscFunctionReturn(0); 408d0578d90SShri Abhyankar } 409