xref: /petsc/src/ts/event/tsevent.c (revision e2cdd85071a133aa8d689c3e944872f1de3c6ded)
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;
1838bf2713SShri 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);
37*e2cdd850SShri Abhyankar   ierr = PetscFree((*event)->fvalue_right);CHKERRQ(ierr);
38*e2cdd850SShri Abhyankar   ierr = PetscFree((*event)->zerocrossing);CHKERRQ(ierr);
39*e2cdd850SShri Abhyankar   ierr = PetscFree((*event)->side);CHKERRQ(ierr);
407dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->direction);CHKERRQ(ierr);
417dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->terminate);CHKERRQ(ierr);
427dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr);
437dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->vtol);CHKERRQ(ierr);
447dbe0728SLisandro Dalcin   for (i=0; i < MAXEVENTRECORDERS; i++) {
457dbe0728SLisandro Dalcin     ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr);
467dbe0728SLisandro Dalcin   }
477dbe0728SLisandro Dalcin   ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr);
487dbe0728SLisandro Dalcin   ierr = PetscFree(*event);CHKERRQ(ierr);
497dbe0728SLisandro Dalcin   PetscFunctionReturn(0);
507dbe0728SLisandro Dalcin }
517dbe0728SLisandro Dalcin 
527dbe0728SLisandro Dalcin #undef __FUNCT__
53e3005195SShri Abhyankar #define __FUNCT__ "TSSetEventTolerances"
54e3005195SShri Abhyankar /*@
55e3005195SShri Abhyankar    TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
56e3005195SShri Abhyankar 
57e3005195SShri Abhyankar    Logically Collective
58e3005195SShri Abhyankar 
59e3005195SShri Abhyankar    Input Arguments:
60e3005195SShri Abhyankar +  ts - time integration context
61e3005195SShri Abhyankar .  tol - scalar tolerance, PETSC_DECIDE to leave current value
62e3005195SShri Abhyankar -  vtol - array of tolerances or NULL, used in preference to tol if present
63e3005195SShri Abhyankar 
642589fa24SLisandro Dalcin    Options Database Keys:
652589fa24SLisandro Dalcin .  -ts_event_tol <tol> tolerance for event zero crossing
66e3005195SShri Abhyankar 
67e3005195SShri Abhyankar    Notes:
68f25fe674SLisandro Dalcin    Must call TSSetEventHandler() before setting the tolerances.
69e3005195SShri Abhyankar 
70e3005195SShri Abhyankar    The size of vtol is equal to the number of events.
71e3005195SShri Abhyankar 
72e3005195SShri Abhyankar    Level: beginner
73e3005195SShri Abhyankar 
74f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler()
75e3005195SShri Abhyankar @*/
766427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
77e3005195SShri Abhyankar {
786427ac75SLisandro Dalcin   TSEvent        event;
79e3005195SShri Abhyankar   PetscInt       i;
80e3005195SShri Abhyankar 
81e3005195SShri Abhyankar   PetscFunctionBegin;
826427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
836427ac75SLisandro Dalcin   if (vtol) PetscValidRealPointer(vtol,3);
84f25fe674SLisandro Dalcin   if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
856427ac75SLisandro Dalcin 
866427ac75SLisandro Dalcin   event = ts->event;
87e3005195SShri Abhyankar   if (vtol) {
88e3005195SShri Abhyankar     for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
89e3005195SShri Abhyankar   } else {
90e3005195SShri Abhyankar     if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
91e3005195SShri Abhyankar       for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
92e3005195SShri Abhyankar     }
93e3005195SShri Abhyankar   }
94e3005195SShri Abhyankar   PetscFunctionReturn(0);
95e3005195SShri Abhyankar }
96e3005195SShri Abhyankar 
97e3005195SShri Abhyankar #undef __FUNCT__
98f25fe674SLisandro Dalcin #define __FUNCT__ "TSSetEventHandler"
992dc7a7e3SShri Abhyankar /*@C
100f25fe674SLisandro Dalcin    TSSetEventHandler - Sets a monitoring function used for detecting events
1012dc7a7e3SShri Abhyankar 
1022dc7a7e3SShri Abhyankar    Logically Collective on TS
1032dc7a7e3SShri Abhyankar 
1042dc7a7e3SShri Abhyankar    Input Parameters:
1052dc7a7e3SShri Abhyankar +  ts - the TS context obtained from TSCreate()
1062dc7a7e3SShri Abhyankar .  nevents - number of local events
107d94325d3SShri Abhyankar .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
108d94325d3SShri Abhyankar                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
109d94325d3SShri Abhyankar .  terminate - flag to indicate whether time stepping should be terminated after
110d94325d3SShri Abhyankar                event is detected (one for each event)
1116427ac75SLisandro Dalcin .  eventhandler - event monitoring routine
1122dc7a7e3SShri Abhyankar .  postevent - [optional] post-event function
1132589fa24SLisandro Dalcin -  ctx       - [optional] user-defined context for private data for the
1142dc7a7e3SShri Abhyankar                event monitor and post event routine (use NULL if no
1152dc7a7e3SShri Abhyankar                context is desired)
1162dc7a7e3SShri Abhyankar 
1176427ac75SLisandro Dalcin    Calling sequence of eventhandler:
1182589fa24SLisandro Dalcin    PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
1192dc7a7e3SShri Abhyankar 
1202dc7a7e3SShri Abhyankar    Input Parameters:
1212dc7a7e3SShri Abhyankar +  ts  - the TS context
1222dc7a7e3SShri Abhyankar .  t   - current time
1232dc7a7e3SShri Abhyankar .  U   - current iterate
1246427ac75SLisandro Dalcin -  ctx - [optional] context passed with eventhandler
1252dc7a7e3SShri Abhyankar 
1262dc7a7e3SShri Abhyankar    Output parameters:
127d94325d3SShri Abhyankar .  fvalue    - function value of events at time t
1282dc7a7e3SShri Abhyankar 
1292dc7a7e3SShri Abhyankar    Calling sequence of postevent:
1302589fa24SLisandro Dalcin    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero[], PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
1312dc7a7e3SShri Abhyankar 
1322dc7a7e3SShri Abhyankar    Input Parameters:
1332dc7a7e3SShri Abhyankar +  ts - the TS context
1342dc7a7e3SShri Abhyankar .  nevents_zero - number of local events whose event function is zero
1352dc7a7e3SShri Abhyankar .  events_zero  - indices of local events which have reached zero
1362dc7a7e3SShri Abhyankar .  t            - current time
1372dc7a7e3SShri Abhyankar .  U            - current solution
138031fbad4SShri Abhyankar .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
1396427ac75SLisandro Dalcin -  ctx          - the context passed with eventhandler
1402dc7a7e3SShri Abhyankar 
1412dc7a7e3SShri Abhyankar    Level: intermediate
1422dc7a7e3SShri Abhyankar 
1432dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor
1442dc7a7e3SShri Abhyankar 
1452dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
1462dc7a7e3SShri Abhyankar @*/
1472589fa24SLisandro 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)
1482dc7a7e3SShri Abhyankar {
1492dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
1502dc7a7e3SShri Abhyankar   TSEvent        event;
151d94325d3SShri Abhyankar   PetscInt       i;
152006e6a18SShri Abhyankar   PetscBool      flg;
153e3005195SShri Abhyankar   PetscReal      tol=1e-6;
1542dc7a7e3SShri Abhyankar 
1552dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1566427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1576427ac75SLisandro Dalcin   PetscValidIntPointer(direction,2);
1586427ac75SLisandro Dalcin   PetscValidIntPointer(terminate,3);
1596427ac75SLisandro Dalcin 
1606427ac75SLisandro Dalcin   ierr = PetscNewLog(ts,&event);CHKERRQ(ierr);
161854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr);
162854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr);
163*e2cdd850SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->fvalue_right);CHKERRQ(ierr);
164*e2cdd850SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->zerocrossing);CHKERRQ(ierr);
165*e2cdd850SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->side);CHKERRQ(ierr);
166854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr);
167854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr);
168e3005195SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr);
169d94325d3SShri Abhyankar   for (i=0; i < nevents; i++) {
170d94325d3SShri Abhyankar     event->direction[i] = direction[i];
171d94325d3SShri Abhyankar     event->terminate[i] = terminate[i];
172*e2cdd850SShri Abhyankar     event->zerocrossing[i] = PETSC_FALSE;
173*e2cdd850SShri Abhyankar     event->side[i] = 0;
174d94325d3SShri Abhyankar   }
175854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr);
1762589fa24SLisandro Dalcin   event->nevents = nevents;
1776427ac75SLisandro Dalcin   event->eventhandler = eventhandler;
1782dc7a7e3SShri Abhyankar   event->postevent = postevent;
1796427ac75SLisandro Dalcin   event->ctx = ctx;
1802dc7a7e3SShri Abhyankar 
181f7aea88cSShri Abhyankar   for (i=0; i < MAXEVENTRECORDERS; i++) {
1822589fa24SLisandro Dalcin     ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr);
183f7aea88cSShri Abhyankar   }
184f7aea88cSShri Abhyankar 
185a9514d71SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr);
1862dc7a7e3SShri Abhyankar   {
187e3005195SShri Abhyankar     ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr);
188006e6a18SShri Abhyankar     ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr);
1892dc7a7e3SShri Abhyankar   }
1909e12be75SShri Abhyankar   PetscOptionsEnd();
191e3005195SShri Abhyankar   for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
1926427ac75SLisandro Dalcin   if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);}
1939e12be75SShri Abhyankar 
1946427ac75SLisandro Dalcin   ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr);
195d94325d3SShri Abhyankar   ts->event = event;
1962dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
1972dc7a7e3SShri Abhyankar }
1982dc7a7e3SShri Abhyankar 
199031fbad4SShri Abhyankar /*
2004597913aSLisandro Dalcin    Helper rutine to handle user postenvents and recording
201031fbad4SShri Abhyankar */
202031fbad4SShri Abhyankar #undef __FUNCT__
203031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent"
2044597913aSLisandro Dalcin static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
2052dc7a7e3SShri Abhyankar {
2062dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
2072dc7a7e3SShri Abhyankar   TSEvent        event = ts->event;
2082dc7a7e3SShri Abhyankar   PetscBool      terminate = PETSC_FALSE;
209d0578d90SShri Abhyankar   PetscInt       i,ctr,stepnum;
2102dc7a7e3SShri Abhyankar   PetscBool      ts_terminate;
2114597913aSLisandro Dalcin   PetscBool      forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
2122dc7a7e3SShri Abhyankar 
2132dc7a7e3SShri Abhyankar   PetscFunctionBegin;
2142dc7a7e3SShri Abhyankar   if (event->postevent) {
2154597913aSLisandro Dalcin     ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
2162dc7a7e3SShri Abhyankar   }
2174597913aSLisandro Dalcin 
2184597913aSLisandro Dalcin   /* Handle termination events */
2194597913aSLisandro Dalcin   for (i=0; i < event->nevents_zero; i++) {
2204597913aSLisandro Dalcin     terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]);
2212dc7a7e3SShri Abhyankar   }
222b2566f29SBarry Smith   ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr);
2239e12be75SShri Abhyankar   if (ts_terminate) {
2242dc7a7e3SShri Abhyankar     ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);
2252dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
2262dc7a7e3SShri Abhyankar   } else {
2272dc7a7e3SShri Abhyankar     event->status = TSEVENT_RESET_NEXTSTEP;
2282dc7a7e3SShri Abhyankar   }
229f7aea88cSShri Abhyankar 
2304597913aSLisandro Dalcin   event->ptime_prev = t;
2314597913aSLisandro Dalcin   /* Reset event residual functions as states might get changed by the postevent callback */
2324597913aSLisandro Dalcin   if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);}
2334597913aSLisandro Dalcin   /* Cache current event residual functions */
2344597913aSLisandro Dalcin   for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
2354597913aSLisandro Dalcin 
236d0578d90SShri Abhyankar   /* Record the event in the event recorder */
237d0578d90SShri Abhyankar   ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr);
238f7aea88cSShri Abhyankar   ctr = event->recorder.ctr;
2396c4ed002SBarry Smith   if (ctr == MAXEVENTRECORDERS) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS);
240f7aea88cSShri Abhyankar   event->recorder.time[ctr] = t;
241d0578d90SShri Abhyankar   event->recorder.stepnum[ctr] = stepnum;
2424597913aSLisandro Dalcin   event->recorder.nevents[ctr] = event->nevents_zero;
2434597913aSLisandro Dalcin   for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
244f7aea88cSShri Abhyankar   event->recorder.ctr++;
2452dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
2462dc7a7e3SShri Abhyankar }
2472dc7a7e3SShri Abhyankar 
248*e2cdd850SShri Abhyankar /* Uses Anderson-Bjorck variant of regular falsi method */
249*e2cdd850SShri Abhyankar PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t, PetscReal tright, PetscReal fleft, PetscScalar f, PetscScalar fright, PetscInt side,PetscReal dt)
25038bf2713SShri Abhyankar {
251*e2cdd850SShri Abhyankar   PetscReal scal=1.0;
252*e2cdd850SShri Abhyankar   if(PetscRealPart(fleft)*PetscRealPart(f) < 0) {
253*e2cdd850SShri Abhyankar     if(side == 1) {
254*e2cdd850SShri Abhyankar       scal = 1.0 - PetscRealPart(f)/PetscRealPart(fright);
255*e2cdd850SShri Abhyankar       if(scal < 0) scal = 0.5;
256*e2cdd850SShri Abhyankar     }
257*e2cdd850SShri Abhyankar     dt = PetscMin(dt,(scal*fleft*t - f*tleft)/(scal*fleft - f) - tleft);
258*e2cdd850SShri Abhyankar   } else {
259*e2cdd850SShri Abhyankar     if(side == -1) {
260*e2cdd850SShri Abhyankar       scal = 1.0 - PetscRealPart(f)/PetscRealPart(fleft);
261*e2cdd850SShri Abhyankar       if(scal < 0) scal = 0.5;
262*e2cdd850SShri Abhyankar     }
263*e2cdd850SShri Abhyankar     dt = PetscMin(dt,(f*tright - scal*fright*t)/(f - scal*fright) - t);
264*e2cdd850SShri Abhyankar   }
26538bf2713SShri Abhyankar   return dt;
26638bf2713SShri Abhyankar }
267*e2cdd850SShri Abhyankar 
26838bf2713SShri Abhyankar 
2692dc7a7e3SShri Abhyankar #undef __FUNCT__
2706427ac75SLisandro Dalcin #define __FUNCT__ "TSEventHandler"
2716427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts)
2722dc7a7e3SShri Abhyankar {
2732dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
2746427ac75SLisandro Dalcin   TSEvent        event;
2752dc7a7e3SShri Abhyankar   PetscReal      t;
2762dc7a7e3SShri Abhyankar   Vec            U;
2772dc7a7e3SShri Abhyankar   PetscInt       i;
2787dbe0728SLisandro Dalcin   PetscReal      dt,dt_min;
2792dc7a7e3SShri Abhyankar   PetscInt       rollback=0,in[2],out[2];
2809e12be75SShri Abhyankar   PetscInt       fvalue_sign,fvalueprev_sign;
2812dc7a7e3SShri Abhyankar 
2822dc7a7e3SShri Abhyankar   PetscFunctionBegin;
2836427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2846427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
2856427ac75SLisandro Dalcin   event = ts->event;
2862dc7a7e3SShri Abhyankar 
2872dc7a7e3SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
2887dbe0728SLisandro Dalcin   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
2892dc7a7e3SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
2902dc7a7e3SShri Abhyankar 
2917dbe0728SLisandro Dalcin   if (event->status == TSEVENT_NONE) {
2927dbe0728SLisandro Dalcin     event->timestep_prev = dt;
2937dbe0728SLisandro Dalcin   }
2947dbe0728SLisandro Dalcin 
2952dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_RESET_NEXTSTEP) {
2967dbe0728SLisandro Dalcin     /* Restore previous time step */
2977dbe0728SLisandro Dalcin     dt = event->timestep_prev;
2987dbe0728SLisandro Dalcin     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
2992dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
3002dc7a7e3SShri Abhyankar   }
3012dc7a7e3SShri Abhyankar 
3022dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_NONE) {
3037dbe0728SLisandro Dalcin     event->ptime_end = t;
3042dc7a7e3SShri Abhyankar   }
3052dc7a7e3SShri Abhyankar 
3066427ac75SLisandro Dalcin   ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);
3079e12be75SShri Abhyankar 
3082dc7a7e3SShri Abhyankar   for (i=0; i < event->nevents; i++) {
309e3005195SShri Abhyankar     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
3102dc7a7e3SShri Abhyankar       event->status = TSEVENT_ZERO;
311*e2cdd850SShri Abhyankar       event->fvalue_right[i] = event->fvalue[i];
3129e12be75SShri Abhyankar       continue;
313006e6a18SShri Abhyankar     }
314d94325d3SShri Abhyankar     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
315d94325d3SShri Abhyankar     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
316e3005195SShri Abhyankar     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
317d94325d3SShri Abhyankar       switch (event->direction[i]) {
318d94325d3SShri Abhyankar       case -1:
319d94325d3SShri Abhyankar         if (fvalue_sign < 0) {
3202dc7a7e3SShri Abhyankar 	  rollback = 1;
321*e2cdd850SShri Abhyankar 
322*e2cdd850SShri Abhyankar 	  /* Compute new time step */
323*e2cdd850SShri Abhyankar 	  dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
324*e2cdd850SShri Abhyankar 
3256427ac75SLisandro Dalcin 	  if (event->monitor) {
32638bf2713SShri 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);
327006e6a18SShri Abhyankar 	  }
328*e2cdd850SShri Abhyankar 	  event->fvalue_right[i] = event->fvalue[i];
329*e2cdd850SShri Abhyankar 	  event->side[i] = 1;
330*e2cdd850SShri Abhyankar 
331*e2cdd850SShri Abhyankar 	  if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
332*e2cdd850SShri Abhyankar 	  event->status = TSEVENT_LOCATED_INTERVAL;
3332dc7a7e3SShri Abhyankar         }
334d94325d3SShri Abhyankar         break;
335d94325d3SShri Abhyankar       case 1:
336d94325d3SShri Abhyankar         if (fvalue_sign > 0) {
337d94325d3SShri Abhyankar 	  rollback = 1;
338*e2cdd850SShri Abhyankar 
339*e2cdd850SShri Abhyankar  	  /* Compute new time step */
340*e2cdd850SShri Abhyankar 	  dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
341*e2cdd850SShri Abhyankar 
3426427ac75SLisandro Dalcin 	  if (event->monitor) {
34338bf2713SShri 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);
344006e6a18SShri Abhyankar 	  }
345*e2cdd850SShri Abhyankar 	  event->fvalue_right[i] = event->fvalue[i];
346*e2cdd850SShri Abhyankar 	  event->side[i] = 1;
347*e2cdd850SShri Abhyankar 
348*e2cdd850SShri Abhyankar      	  if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
349*e2cdd850SShri Abhyankar 	  event->status = TSEVENT_LOCATED_INTERVAL;
3502dc7a7e3SShri Abhyankar         }
351d94325d3SShri Abhyankar         break;
352d94325d3SShri Abhyankar       case 0:
353d94325d3SShri Abhyankar 	rollback = 1;
354*e2cdd850SShri Abhyankar 
355*e2cdd850SShri Abhyankar 	/* Compute new time step */
356*e2cdd850SShri Abhyankar 	dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
357*e2cdd850SShri Abhyankar 
3586427ac75SLisandro Dalcin 	if (event->monitor) {
35938bf2713SShri 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);
360006e6a18SShri Abhyankar 	}
361*e2cdd850SShri Abhyankar 	event->fvalue_right[i] = event->fvalue[i];
362*e2cdd850SShri Abhyankar 	event->side[i] = 1;
363*e2cdd850SShri Abhyankar 
364*e2cdd850SShri Abhyankar 	if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
365*e2cdd850SShri Abhyankar 	event->status = TSEVENT_LOCATED_INTERVAL;
366*e2cdd850SShri Abhyankar 
367d94325d3SShri Abhyankar         break;
368d94325d3SShri Abhyankar       }
369d94325d3SShri Abhyankar     }
370d94325d3SShri Abhyankar   }
3717dbe0728SLisandro Dalcin 
3722589fa24SLisandro Dalcin   in[0] = event->status; in[1] = rollback;
3737dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
3742589fa24SLisandro Dalcin   event->status = (TSEventStatus)out[0]; rollback = out[1];
3752589fa24SLisandro Dalcin   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
3762dc7a7e3SShri Abhyankar 
3777dbe0728SLisandro Dalcin   event->nevents_zero = 0;
3789e12be75SShri Abhyankar   if (event->status == TSEVENT_ZERO) {
3799e12be75SShri Abhyankar     for (i=0; i < event->nevents; i++) {
3809e12be75SShri Abhyankar       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
3819e12be75SShri Abhyankar         event->events_zero[event->nevents_zero++] = i;
3826427ac75SLisandro Dalcin         if (event->monitor) {
38338bf2713SShri 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);
3849e12be75SShri Abhyankar         }
385*e2cdd850SShri Abhyankar 	event->zerocrossing[i] = PETSC_FALSE;
3869e12be75SShri Abhyankar       }
387*e2cdd850SShri Abhyankar       event->side[i] = 0;
3889e12be75SShri Abhyankar     }
3894597913aSLisandro Dalcin     ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr);
3909e12be75SShri Abhyankar 
391*e2cdd850SShri Abhyankar     dt = event->ptime_end - t;
3927dbe0728SLisandro Dalcin     if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */
3939e12be75SShri Abhyankar     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
39438bf2713SShri Abhyankar     event->iterctr = 0;
3959e12be75SShri Abhyankar     PetscFunctionReturn(0);
3969e12be75SShri Abhyankar   }
3979e12be75SShri Abhyankar 
3982dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_LOCATED_INTERVAL) {
3992dc7a7e3SShri Abhyankar     ierr = TSRollBack(ts);CHKERRQ(ierr);
400c540466cSShri Abhyankar     ts->steps--;
401c540466cSShri Abhyankar     ts->total_steps--;
4022589fa24SLisandro Dalcin     ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr);
4032dc7a7e3SShri Abhyankar     event->status = TSEVENT_PROCESSING;
404*e2cdd850SShri Abhyankar 
405*e2cdd850SShri Abhyankar     event->ptime_right = t;
4062dc7a7e3SShri Abhyankar   } else {
4072dc7a7e3SShri Abhyankar     for(i=0; i < event->nevents; i++) {
408*e2cdd850SShri Abhyankar       if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) {
409*e2cdd850SShri Abhyankar 	/* Compute new time step */
410*e2cdd850SShri Abhyankar 	dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
411*e2cdd850SShri Abhyankar 
412*e2cdd850SShri Abhyankar 	event->side[i] = -1;
413*e2cdd850SShri Abhyankar       }
4142dc7a7e3SShri Abhyankar       event->fvalue_prev[i] = event->fvalue[i];
4152dc7a7e3SShri Abhyankar     }
41638bf2713SShri Abhyankar 
417*e2cdd850SShri Abhyankar     if (event->monitor && event->status == TSEVENT_PROCESSING) {
41838bf2713SShri 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);
4192dc7a7e3SShri Abhyankar     }
42038bf2713SShri Abhyankar     event->ptime_prev = t;
42138bf2713SShri Abhyankar   }
42238bf2713SShri Abhyankar 
42338bf2713SShri Abhyankar   if(event->status == TSEVENT_PROCESSING) event->iterctr++;
4249e12be75SShri Abhyankar 
4257dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
4267dbe0728SLisandro Dalcin   ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr);
4272dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
4282dc7a7e3SShri Abhyankar }
4292dc7a7e3SShri Abhyankar 
430d0578d90SShri Abhyankar #undef __FUNCT__
4316427ac75SLisandro Dalcin #define __FUNCT__ "TSAdjointEventHandler"
4326427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts)
433d0578d90SShri Abhyankar {
434d0578d90SShri Abhyankar   PetscErrorCode ierr;
4356427ac75SLisandro Dalcin   TSEvent        event;
436d0578d90SShri Abhyankar   PetscReal      t;
437d0578d90SShri Abhyankar   Vec            U;
438d0578d90SShri Abhyankar   PetscInt       ctr;
439d0578d90SShri Abhyankar   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
440d0578d90SShri Abhyankar 
441d0578d90SShri Abhyankar   PetscFunctionBegin;
4426427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4436427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
4446427ac75SLisandro Dalcin   event = ts->event;
445d0578d90SShri Abhyankar 
446d0578d90SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
447d0578d90SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
448d0578d90SShri Abhyankar 
449d0578d90SShri Abhyankar   ctr = event->recorder.ctr-1;
450bcbf8bb3SShri Abhyankar   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
451d0578d90SShri Abhyankar     /* Call the user postevent function */
452d0578d90SShri Abhyankar     if (event->postevent) {
4536427ac75SLisandro Dalcin       ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
454d0578d90SShri Abhyankar       event->recorder.ctr--;
455d0578d90SShri Abhyankar     }
456d0578d90SShri Abhyankar   }
457d0578d90SShri Abhyankar 
458d0578d90SShri Abhyankar   PetscBarrier((PetscObject)ts);
459d0578d90SShri Abhyankar   PetscFunctionReturn(0);
460d0578d90SShri Abhyankar }
461