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 38e3005195SShri Abhyankar - -ts_event_tol <tol> tolerance for event zero crossing 39e3005195SShri Abhyankar 40e3005195SShri Abhyankar Notes: 41*f25fe674SLisandro Dalcin Must call TSSetEventHandler() before setting the tolerances. 42e3005195SShri Abhyankar 43e3005195SShri Abhyankar The size of vtol is equal to the number of events. 44e3005195SShri Abhyankar 45e3005195SShri Abhyankar Level: beginner 46e3005195SShri Abhyankar 47*f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler() 48e3005195SShri Abhyankar @*/ 496427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[]) 50e3005195SShri Abhyankar { 516427ac75SLisandro Dalcin TSEvent event; 52e3005195SShri Abhyankar PetscInt i; 53e3005195SShri Abhyankar 54e3005195SShri Abhyankar PetscFunctionBegin; 556427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 566427ac75SLisandro Dalcin if (vtol) PetscValidRealPointer(vtol,3); 57*f25fe674SLisandro Dalcin if(!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()"); 586427ac75SLisandro Dalcin 596427ac75SLisandro Dalcin event = ts->event; 60e3005195SShri Abhyankar if (vtol) { 61e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 62e3005195SShri Abhyankar } else { 63e3005195SShri Abhyankar if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 64e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 65e3005195SShri Abhyankar } 66e3005195SShri Abhyankar } 67e3005195SShri Abhyankar PetscFunctionReturn(0); 68e3005195SShri Abhyankar } 69e3005195SShri Abhyankar 70e3005195SShri Abhyankar #undef __FUNCT__ 71*f25fe674SLisandro Dalcin #define __FUNCT__ "TSSetEventHandler" 722dc7a7e3SShri Abhyankar /*@C 73*f25fe674SLisandro Dalcin TSSetEventHandler - Sets a monitoring function used for detecting events 742dc7a7e3SShri Abhyankar 752dc7a7e3SShri Abhyankar Logically Collective on TS 762dc7a7e3SShri Abhyankar 772dc7a7e3SShri Abhyankar Input Parameters: 782dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 792dc7a7e3SShri Abhyankar . nevents - number of local events 80d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 81d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 82d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 83d94325d3SShri Abhyankar event is detected (one for each event) 846427ac75SLisandro Dalcin . eventhandler - event monitoring routine 852dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 862dc7a7e3SShri Abhyankar - mectx - [optional] user-defined context for private data for the 872dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 882dc7a7e3SShri Abhyankar context is desired) 892dc7a7e3SShri Abhyankar 906427ac75SLisandro Dalcin Calling sequence of eventhandler: 916427ac75SLisandro Dalcin PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void* mectx) 922dc7a7e3SShri Abhyankar 932dc7a7e3SShri Abhyankar Input Parameters: 942dc7a7e3SShri Abhyankar + ts - the TS context 952dc7a7e3SShri Abhyankar . t - current time 962dc7a7e3SShri Abhyankar . U - current iterate 976427ac75SLisandro Dalcin - ctx - [optional] context passed with eventhandler 982dc7a7e3SShri Abhyankar 992dc7a7e3SShri Abhyankar Output parameters: 100d94325d3SShri Abhyankar . fvalue - function value of events at time t 1012dc7a7e3SShri Abhyankar 1022dc7a7e3SShri Abhyankar Calling sequence of postevent: 103031fbad4SShri Abhyankar PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero, PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 1042dc7a7e3SShri Abhyankar 1052dc7a7e3SShri Abhyankar Input Parameters: 1062dc7a7e3SShri Abhyankar + ts - the TS context 1072dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 1082dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 1092dc7a7e3SShri Abhyankar . t - current time 1102dc7a7e3SShri Abhyankar . U - current solution 111031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 1126427ac75SLisandro Dalcin - ctx - the context passed with eventhandler 1132dc7a7e3SShri Abhyankar 1142dc7a7e3SShri Abhyankar Level: intermediate 1152dc7a7e3SShri Abhyankar 1162dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor 1172dc7a7e3SShri Abhyankar 1182dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 1192dc7a7e3SShri Abhyankar @*/ 120*f25fe674SLisandro 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) 1212dc7a7e3SShri Abhyankar { 1222dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1232dc7a7e3SShri Abhyankar TSEvent event; 124d94325d3SShri Abhyankar PetscInt i; 125006e6a18SShri Abhyankar PetscBool flg; 126e3005195SShri Abhyankar PetscReal tol=1e-6; 1272dc7a7e3SShri Abhyankar 1282dc7a7e3SShri Abhyankar PetscFunctionBegin; 1296427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1306427ac75SLisandro Dalcin PetscValidIntPointer(direction,2); 1316427ac75SLisandro Dalcin PetscValidIntPointer(terminate,3); 1326427ac75SLisandro Dalcin 1336427ac75SLisandro Dalcin ierr = PetscNewLog(ts,&event);CHKERRQ(ierr); 134854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 135854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 136854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 137854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 138e3005195SShri Abhyankar ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr); 139d94325d3SShri Abhyankar for (i=0; i < nevents; i++) { 140d94325d3SShri Abhyankar event->direction[i] = direction[i]; 141d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 142d94325d3SShri Abhyankar } 143854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 1446427ac75SLisandro Dalcin event->eventhandler = eventhandler; 1452dc7a7e3SShri Abhyankar event->postevent = postevent; 1466427ac75SLisandro Dalcin event->ctx = ctx; 1472dc7a7e3SShri Abhyankar event->nevents = nevents; 1482dc7a7e3SShri Abhyankar 149f7aea88cSShri Abhyankar for (i=0; i < MAXEVENTRECORDERS; i++) { 150f7aea88cSShri Abhyankar ierr = PetscMalloc1(nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 151f7aea88cSShri Abhyankar } 152f7aea88cSShri Abhyankar 153a9514d71SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr); 1542dc7a7e3SShri Abhyankar { 155e3005195SShri Abhyankar ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr); 156006e6a18SShri Abhyankar ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr); 1572dc7a7e3SShri Abhyankar } 1589e12be75SShri Abhyankar PetscOptionsEnd(); 159e3005195SShri Abhyankar for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 1606427ac75SLisandro Dalcin if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);} 1619e12be75SShri Abhyankar 1626427ac75SLisandro Dalcin ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr); 163d94325d3SShri Abhyankar ts->event = event; 1642dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1652dc7a7e3SShri Abhyankar } 1662dc7a7e3SShri Abhyankar 1672dc7a7e3SShri Abhyankar #undef __FUNCT__ 1682dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent" 169031fbad4SShri Abhyankar /* 170031fbad4SShri Abhyankar TSPostEvent - Does post event processing by calling the user-defined postevent function 171031fbad4SShri Abhyankar 172031fbad4SShri Abhyankar Logically Collective on TS 173031fbad4SShri Abhyankar 174031fbad4SShri Abhyankar Input Parameters: 175031fbad4SShri Abhyankar + ts - the TS context 176031fbad4SShri Abhyankar . nevents_zero - number of local events whose event function is zero 177031fbad4SShri Abhyankar . events_zero - indices of local events which have reached zero 178031fbad4SShri Abhyankar . t - current time 179031fbad4SShri Abhyankar . U - current solution 1806427ac75SLisandro Dalcin - forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 181031fbad4SShri Abhyankar 182031fbad4SShri Abhyankar Level: intermediate 183031fbad4SShri Abhyankar 184031fbad4SShri Abhyankar .keywords: TS, event, set, monitor 185031fbad4SShri Abhyankar 186*f25fe674SLisandro Dalcin .seealso: TSSetEventHandler(), TSEvent 187031fbad4SShri Abhyankar */ 188031fbad4SShri Abhyankar #undef __FUNCT__ 189031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent" 1906427ac75SLisandro Dalcin PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve) 1912dc7a7e3SShri Abhyankar { 1922dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1932dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1942dc7a7e3SShri Abhyankar PetscBool terminate=PETSC_FALSE; 195d0578d90SShri Abhyankar PetscInt i,ctr,stepnum; 1962dc7a7e3SShri Abhyankar PetscBool ts_terminate; 1972dc7a7e3SShri Abhyankar 1982dc7a7e3SShri Abhyankar PetscFunctionBegin; 1992dc7a7e3SShri Abhyankar if (event->postevent) { 2006427ac75SLisandro Dalcin ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 2012dc7a7e3SShri Abhyankar } 2022dc7a7e3SShri Abhyankar for(i = 0; i < nevents_zero;i++) { 203e105d053SSatish Balay terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]); 2042dc7a7e3SShri Abhyankar } 205b2566f29SBarry Smith ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 2069e12be75SShri Abhyankar if (ts_terminate) { 2072dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 2082dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2092dc7a7e3SShri Abhyankar } else { 2102dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 2112dc7a7e3SShri Abhyankar } 212f7aea88cSShri Abhyankar 213d0578d90SShri Abhyankar /* Record the event in the event recorder */ 214d0578d90SShri Abhyankar ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr); 215f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 2166c4ed002SBarry Smith if (ctr == MAXEVENTRECORDERS) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS); 217f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 218d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 219f7aea88cSShri Abhyankar event->recorder.nevents[ctr] = nevents_zero; 220f7aea88cSShri Abhyankar for(i=0; i < nevents_zero; i++) event->recorder.eventidx[ctr][i] = events_zero[i]; 221f7aea88cSShri Abhyankar event->recorder.ctr++; 222f7aea88cSShri Abhyankar 22373967516SShri Abhyankar /* Reset the event residual functions as states might get changed by the postevent callback */ 2246427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr); 22573967516SShri Abhyankar event->ptime_prev = t; 2262dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2272dc7a7e3SShri Abhyankar } 2282dc7a7e3SShri Abhyankar 2292dc7a7e3SShri Abhyankar #undef __FUNCT__ 2306427ac75SLisandro Dalcin #define __FUNCT__ "TSEventDestroy" 2316427ac75SLisandro Dalcin PetscErrorCode TSEventDestroy(TSEvent *event) 2322dc7a7e3SShri Abhyankar { 2332dc7a7e3SShri Abhyankar PetscErrorCode ierr; 234f7aea88cSShri Abhyankar PetscInt i; 2352dc7a7e3SShri Abhyankar 2362dc7a7e3SShri Abhyankar PetscFunctionBegin; 2376427ac75SLisandro Dalcin PetscValidPointer(event,1); 2386427ac75SLisandro Dalcin if (!*event) PetscFunctionReturn(0); 2392dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 2402dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 2412dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 242d94325d3SShri Abhyankar ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 2432dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 244e3005195SShri Abhyankar ierr = PetscFree((*event)->vtol);CHKERRQ(ierr); 245f7aea88cSShri Abhyankar for(i=0; i < MAXEVENTRECORDERS; i++) { 246302440fdSBarry Smith ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr); 247f7aea88cSShri Abhyankar } 2486427ac75SLisandro Dalcin ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr); 2492dc7a7e3SShri Abhyankar ierr = PetscFree(*event);CHKERRQ(ierr); 2502dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2512dc7a7e3SShri Abhyankar } 2522dc7a7e3SShri Abhyankar 2532dc7a7e3SShri Abhyankar #undef __FUNCT__ 2546427ac75SLisandro Dalcin #define __FUNCT__ "TSEventHandler" 2556427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts) 2562dc7a7e3SShri Abhyankar { 2572dc7a7e3SShri Abhyankar PetscErrorCode ierr; 2586427ac75SLisandro Dalcin TSEvent event; 2592dc7a7e3SShri Abhyankar PetscReal t; 2602dc7a7e3SShri Abhyankar Vec U; 2612dc7a7e3SShri Abhyankar PetscInt i; 2622dc7a7e3SShri Abhyankar PetscReal dt; 2632dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 264031fbad4SShri Abhyankar PetscBool forwardsolve=PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 2659e12be75SShri Abhyankar PetscInt fvalue_sign,fvalueprev_sign; 2662dc7a7e3SShri Abhyankar 2672dc7a7e3SShri Abhyankar PetscFunctionBegin; 2686427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2696427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 2706427ac75SLisandro Dalcin event = ts->event; 2712dc7a7e3SShri Abhyankar 2722dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 2732dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 2742dc7a7e3SShri Abhyankar 2752dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 2762dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 2772dc7a7e3SShri Abhyankar /* Take initial time step */ 2782dc7a7e3SShri Abhyankar dt = event->initial_timestep; 2792dc7a7e3SShri Abhyankar ts->time_step = dt; 2802dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2812dc7a7e3SShri Abhyankar } 2822dc7a7e3SShri Abhyankar 2832dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 2842dc7a7e3SShri Abhyankar event->tstepend = t; 2852dc7a7e3SShri Abhyankar } 2862dc7a7e3SShri Abhyankar 2872dc7a7e3SShri Abhyankar event->nevents_zero = 0; 2882dc7a7e3SShri Abhyankar 2896427ac75SLisandro Dalcin ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); 2909e12be75SShri Abhyankar 2912dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 292e3005195SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 2932dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 2949e12be75SShri Abhyankar continue; 295006e6a18SShri Abhyankar } 296d94325d3SShri Abhyankar fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 297d94325d3SShri Abhyankar fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 298e3005195SShri Abhyankar if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) { 299d94325d3SShri Abhyankar switch (event->direction[i]) { 300d94325d3SShri Abhyankar case -1: 301d94325d3SShri Abhyankar if (fvalue_sign < 0) { 3022dc7a7e3SShri Abhyankar rollback = 1; 3032dc7a7e3SShri Abhyankar /* Compute linearly interpolated new time step */ 304e105d053SSatish Balay dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3056427ac75SLisandro Dalcin if(event->monitor) { 3066427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent : Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 307006e6a18SShri Abhyankar } 3082dc7a7e3SShri Abhyankar } 309d94325d3SShri Abhyankar break; 310d94325d3SShri Abhyankar case 1: 311d94325d3SShri Abhyankar if (fvalue_sign > 0) { 312d94325d3SShri Abhyankar rollback = 1; 313d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 314d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3156427ac75SLisandro Dalcin if(event->monitor) { 3166427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent : Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 317006e6a18SShri Abhyankar } 3182dc7a7e3SShri Abhyankar } 319d94325d3SShri Abhyankar break; 320d94325d3SShri Abhyankar case 0: 321d94325d3SShri Abhyankar rollback = 1; 322d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 323d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 3246427ac75SLisandro Dalcin if(event->monitor) { 3256427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent : Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 326006e6a18SShri Abhyankar } 327d94325d3SShri Abhyankar break; 328d94325d3SShri Abhyankar } 329d94325d3SShri Abhyankar } 330d94325d3SShri Abhyankar } 331d94325d3SShri Abhyankar if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 332d94325d3SShri Abhyankar 3332dc7a7e3SShri Abhyankar in[0] = event->status; 3342dc7a7e3SShri Abhyankar in[1] = rollback; 335b2566f29SBarry Smith ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 3362dc7a7e3SShri Abhyankar 337128a6747SShri Abhyankar event->status = (TSEventStatus)out[0]; 3382dc7a7e3SShri Abhyankar rollback = out[1]; 3392dc7a7e3SShri Abhyankar if (rollback) { 3402dc7a7e3SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 3412dc7a7e3SShri Abhyankar } 3422dc7a7e3SShri Abhyankar 3439e12be75SShri Abhyankar if (event->status == TSEVENT_ZERO) { 3449e12be75SShri Abhyankar for(i=0; i < event->nevents; i++) { 3459e12be75SShri Abhyankar if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 3469e12be75SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 3476427ac75SLisandro Dalcin if(event->monitor) { 3486427ac75SLisandro Dalcin ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent : Event %D zero crossing at time %g\n",i,(double)t);CHKERRQ(ierr); 3499e12be75SShri Abhyankar } 3509e12be75SShri Abhyankar } 3519e12be75SShri Abhyankar } 3529e12be75SShri Abhyankar 3536427ac75SLisandro Dalcin ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve);CHKERRQ(ierr); 3549e12be75SShri Abhyankar 3559e12be75SShri Abhyankar dt = event->tstepend-t; 3569e12be75SShri Abhyankar if(PetscAbsReal(dt) < PETSC_SMALL) dt += event->initial_timestep; 3579e12be75SShri Abhyankar ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 3589e12be75SShri Abhyankar PetscFunctionReturn(0); 3599e12be75SShri Abhyankar } 3609e12be75SShri Abhyankar 3612dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 3622dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 363c540466cSShri Abhyankar ts->steps--; 364c540466cSShri Abhyankar ts->total_steps--; 365db84a1feSShri Abhyankar ts->reason = TS_CONVERGED_ITERATING; 3662dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 3672dc7a7e3SShri Abhyankar } else { 3682dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 3692dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 3702dc7a7e3SShri Abhyankar } 3712dc7a7e3SShri Abhyankar event->ptime_prev = t; 3722dc7a7e3SShri Abhyankar if (event->status == TSEVENT_PROCESSING) { 3732dc7a7e3SShri Abhyankar dt = event->tstepend - event->ptime_prev; 3742dc7a7e3SShri Abhyankar } 3752dc7a7e3SShri Abhyankar } 3769e12be75SShri Abhyankar 377b2566f29SBarry Smith ierr = MPIU_Allreduce(&dt,&(ts->time_step),1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 3782dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3792dc7a7e3SShri Abhyankar } 3802dc7a7e3SShri Abhyankar 381d0578d90SShri Abhyankar #undef __FUNCT__ 3826427ac75SLisandro Dalcin #define __FUNCT__ "TSAdjointEventHandler" 3836427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts) 384d0578d90SShri Abhyankar { 385d0578d90SShri Abhyankar PetscErrorCode ierr; 3866427ac75SLisandro Dalcin TSEvent event; 387d0578d90SShri Abhyankar PetscReal t; 388d0578d90SShri Abhyankar Vec U; 389d0578d90SShri Abhyankar PetscInt ctr; 390d0578d90SShri Abhyankar PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 391d0578d90SShri Abhyankar 392d0578d90SShri Abhyankar PetscFunctionBegin; 3936427ac75SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3946427ac75SLisandro Dalcin if (!ts->event) PetscFunctionReturn(0); 3956427ac75SLisandro Dalcin event = ts->event; 396d0578d90SShri Abhyankar 397d0578d90SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 398d0578d90SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 399d0578d90SShri Abhyankar 400d0578d90SShri Abhyankar ctr = event->recorder.ctr-1; 401bcbf8bb3SShri Abhyankar if(ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 402d0578d90SShri Abhyankar /* Call the user postevent function */ 403d0578d90SShri Abhyankar if(event->postevent) { 4046427ac75SLisandro Dalcin ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 405d0578d90SShri Abhyankar event->recorder.ctr--; 406d0578d90SShri Abhyankar } 407d0578d90SShri Abhyankar } 408d0578d90SShri Abhyankar 409d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 410d0578d90SShri Abhyankar PetscFunctionReturn(0); 411d0578d90SShri Abhyankar } 412