xref: /petsc/src/ts/event/tsevent.c (revision 7dbe07286e5afe7161ea7127ae4c682ef0a81679)
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);
17*7dbe0728SLisandro Dalcin   ierr = TSGetTimeStep(ts,&event->timestep_prev);CHKERRQ(ierr);
186fea3669SShri Abhyankar 
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__
27*7dbe0728SLisandro Dalcin #define __FUNCT__ "TSEventDestroy"
28*7dbe0728SLisandro Dalcin PetscErrorCode TSEventDestroy(TSEvent *event)
29*7dbe0728SLisandro Dalcin {
30*7dbe0728SLisandro Dalcin   PetscErrorCode ierr;
31*7dbe0728SLisandro Dalcin   PetscInt       i;
32*7dbe0728SLisandro Dalcin 
33*7dbe0728SLisandro Dalcin   PetscFunctionBegin;
34*7dbe0728SLisandro Dalcin   PetscValidPointer(event,1);
35*7dbe0728SLisandro Dalcin   if (!*event) PetscFunctionReturn(0);
36*7dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr);
37*7dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr);
38*7dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->direction);CHKERRQ(ierr);
39*7dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->terminate);CHKERRQ(ierr);
40*7dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr);
41*7dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->vtol);CHKERRQ(ierr);
42*7dbe0728SLisandro Dalcin   for (i=0; i < MAXEVENTRECORDERS; i++) {
43*7dbe0728SLisandro Dalcin     ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr);
44*7dbe0728SLisandro Dalcin   }
45*7dbe0728SLisandro Dalcin   ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr);
46*7dbe0728SLisandro Dalcin   ierr = PetscFree(*event);CHKERRQ(ierr);
47*7dbe0728SLisandro Dalcin   PetscFunctionReturn(0);
48*7dbe0728SLisandro Dalcin }
49*7dbe0728SLisandro Dalcin 
50*7dbe0728SLisandro Dalcin #undef __FUNCT__
51e3005195SShri Abhyankar #define __FUNCT__ "TSSetEventTolerances"
52e3005195SShri Abhyankar /*@
53e3005195SShri Abhyankar    TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
54e3005195SShri Abhyankar 
55e3005195SShri Abhyankar    Logically Collective
56e3005195SShri Abhyankar 
57e3005195SShri Abhyankar    Input Arguments:
58e3005195SShri Abhyankar +  ts - time integration context
59e3005195SShri Abhyankar .  tol - scalar tolerance, PETSC_DECIDE to leave current value
60e3005195SShri Abhyankar -  vtol - array of tolerances or NULL, used in preference to tol if present
61e3005195SShri Abhyankar 
622589fa24SLisandro Dalcin    Options Database Keys:
632589fa24SLisandro Dalcin .  -ts_event_tol <tol> tolerance for event zero crossing
64e3005195SShri Abhyankar 
65e3005195SShri Abhyankar    Notes:
66f25fe674SLisandro Dalcin    Must call TSSetEventHandler() before setting the tolerances.
67e3005195SShri Abhyankar 
68e3005195SShri Abhyankar    The size of vtol is equal to the number of events.
69e3005195SShri Abhyankar 
70e3005195SShri Abhyankar    Level: beginner
71e3005195SShri Abhyankar 
72f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler()
73e3005195SShri Abhyankar @*/
746427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
75e3005195SShri Abhyankar {
766427ac75SLisandro Dalcin   TSEvent        event;
77e3005195SShri Abhyankar   PetscInt       i;
78e3005195SShri Abhyankar 
79e3005195SShri Abhyankar   PetscFunctionBegin;
806427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
816427ac75SLisandro Dalcin   if (vtol) PetscValidRealPointer(vtol,3);
82f25fe674SLisandro Dalcin   if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
836427ac75SLisandro Dalcin 
846427ac75SLisandro Dalcin   event = ts->event;
85e3005195SShri Abhyankar   if (vtol) {
86e3005195SShri Abhyankar     for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
87e3005195SShri Abhyankar   } else {
88e3005195SShri Abhyankar     if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
89e3005195SShri Abhyankar       for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
90e3005195SShri Abhyankar     }
91e3005195SShri Abhyankar   }
92e3005195SShri Abhyankar   PetscFunctionReturn(0);
93e3005195SShri Abhyankar }
94e3005195SShri Abhyankar 
95e3005195SShri Abhyankar #undef __FUNCT__
96f25fe674SLisandro Dalcin #define __FUNCT__ "TSSetEventHandler"
972dc7a7e3SShri Abhyankar /*@C
98f25fe674SLisandro Dalcin    TSSetEventHandler - Sets a monitoring function used for detecting events
992dc7a7e3SShri Abhyankar 
1002dc7a7e3SShri Abhyankar    Logically Collective on TS
1012dc7a7e3SShri Abhyankar 
1022dc7a7e3SShri Abhyankar    Input Parameters:
1032dc7a7e3SShri Abhyankar +  ts - the TS context obtained from TSCreate()
1042dc7a7e3SShri Abhyankar .  nevents - number of local events
105d94325d3SShri Abhyankar .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
106d94325d3SShri Abhyankar                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
107d94325d3SShri Abhyankar .  terminate - flag to indicate whether time stepping should be terminated after
108d94325d3SShri Abhyankar                event is detected (one for each event)
1096427ac75SLisandro Dalcin .  eventhandler - event monitoring routine
1102dc7a7e3SShri Abhyankar .  postevent - [optional] post-event function
1112589fa24SLisandro Dalcin -  ctx       - [optional] user-defined context for private data for the
1122dc7a7e3SShri Abhyankar                event monitor and post event routine (use NULL if no
1132dc7a7e3SShri Abhyankar                context is desired)
1142dc7a7e3SShri Abhyankar 
1156427ac75SLisandro Dalcin    Calling sequence of eventhandler:
1162589fa24SLisandro Dalcin    PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
1172dc7a7e3SShri Abhyankar 
1182dc7a7e3SShri Abhyankar    Input Parameters:
1192dc7a7e3SShri Abhyankar +  ts  - the TS context
1202dc7a7e3SShri Abhyankar .  t   - current time
1212dc7a7e3SShri Abhyankar .  U   - current iterate
1226427ac75SLisandro Dalcin -  ctx - [optional] context passed with eventhandler
1232dc7a7e3SShri Abhyankar 
1242dc7a7e3SShri Abhyankar    Output parameters:
125d94325d3SShri Abhyankar .  fvalue    - function value of events at time t
1262dc7a7e3SShri Abhyankar 
1272dc7a7e3SShri Abhyankar    Calling sequence of postevent:
1282589fa24SLisandro Dalcin    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero[], PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
1292dc7a7e3SShri Abhyankar 
1302dc7a7e3SShri Abhyankar    Input Parameters:
1312dc7a7e3SShri Abhyankar +  ts - the TS context
1322dc7a7e3SShri Abhyankar .  nevents_zero - number of local events whose event function is zero
1332dc7a7e3SShri Abhyankar .  events_zero  - indices of local events which have reached zero
1342dc7a7e3SShri Abhyankar .  t            - current time
1352dc7a7e3SShri Abhyankar .  U            - current solution
136031fbad4SShri Abhyankar .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
1376427ac75SLisandro Dalcin -  ctx          - the context passed with eventhandler
1382dc7a7e3SShri Abhyankar 
1392dc7a7e3SShri Abhyankar    Level: intermediate
1402dc7a7e3SShri Abhyankar 
1412dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor
1422dc7a7e3SShri Abhyankar 
1432dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
1442dc7a7e3SShri Abhyankar @*/
1452589fa24SLisandro 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)
1462dc7a7e3SShri Abhyankar {
1472dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
1482dc7a7e3SShri Abhyankar   TSEvent        event;
149d94325d3SShri Abhyankar   PetscInt       i;
150006e6a18SShri Abhyankar   PetscBool      flg;
151e3005195SShri Abhyankar   PetscReal      tol=1e-6;
1522dc7a7e3SShri Abhyankar 
1532dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1546427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1556427ac75SLisandro Dalcin   PetscValidIntPointer(direction,2);
1566427ac75SLisandro Dalcin   PetscValidIntPointer(terminate,3);
1576427ac75SLisandro Dalcin 
1586427ac75SLisandro Dalcin   ierr = PetscNewLog(ts,&event);CHKERRQ(ierr);
159854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr);
160854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr);
161854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr);
162854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr);
163e3005195SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr);
164d94325d3SShri Abhyankar   for (i=0; i < nevents; i++) {
165d94325d3SShri Abhyankar     event->direction[i] = direction[i];
166d94325d3SShri Abhyankar     event->terminate[i] = terminate[i];
167d94325d3SShri Abhyankar   }
168854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr);
1692589fa24SLisandro Dalcin   event->nevents = nevents;
1706427ac75SLisandro Dalcin   event->eventhandler = eventhandler;
1712dc7a7e3SShri Abhyankar   event->postevent = postevent;
1726427ac75SLisandro Dalcin   event->ctx = ctx;
1732dc7a7e3SShri Abhyankar 
174f7aea88cSShri Abhyankar   for (i=0; i < MAXEVENTRECORDERS; i++) {
1752589fa24SLisandro Dalcin     ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr);
176f7aea88cSShri Abhyankar   }
177f7aea88cSShri Abhyankar 
178a9514d71SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr);
1792dc7a7e3SShri Abhyankar   {
180e3005195SShri Abhyankar     ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr);
181006e6a18SShri Abhyankar     ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr);
1822dc7a7e3SShri Abhyankar   }
1839e12be75SShri Abhyankar   PetscOptionsEnd();
184e3005195SShri Abhyankar   for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
1856427ac75SLisandro Dalcin   if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);}
1869e12be75SShri Abhyankar 
1876427ac75SLisandro Dalcin   ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr);
188d94325d3SShri Abhyankar   ts->event = event;
1892dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
1902dc7a7e3SShri Abhyankar }
1912dc7a7e3SShri Abhyankar 
1922dc7a7e3SShri Abhyankar #undef __FUNCT__
1932dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent"
194031fbad4SShri Abhyankar /*
195031fbad4SShri Abhyankar    TSPostEvent - Does post event processing by calling the user-defined postevent function
196031fbad4SShri Abhyankar 
197031fbad4SShri Abhyankar    Logically Collective on TS
198031fbad4SShri Abhyankar 
199031fbad4SShri Abhyankar    Input Parameters:
200031fbad4SShri Abhyankar +  ts - the TS context
201031fbad4SShri Abhyankar .  nevents_zero - number of local events whose event function is zero
202031fbad4SShri Abhyankar .  events_zero  - indices of local events which have reached zero
203031fbad4SShri Abhyankar .  t            - current time
204031fbad4SShri Abhyankar .  U            - current solution
2056427ac75SLisandro Dalcin -  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
206031fbad4SShri Abhyankar 
207031fbad4SShri Abhyankar    Level: intermediate
208031fbad4SShri Abhyankar 
209031fbad4SShri Abhyankar .keywords: TS, event, set, monitor
210031fbad4SShri Abhyankar 
211f25fe674SLisandro Dalcin .seealso: TSSetEventHandler(), TSEvent
212031fbad4SShri Abhyankar */
213031fbad4SShri Abhyankar #undef __FUNCT__
214031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent"
2156427ac75SLisandro Dalcin PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve)
2162dc7a7e3SShri Abhyankar {
2172dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
2182dc7a7e3SShri Abhyankar   TSEvent        event=ts->event;
2192dc7a7e3SShri Abhyankar   PetscBool      terminate=PETSC_FALSE;
220d0578d90SShri Abhyankar   PetscInt       i,ctr,stepnum;
2212dc7a7e3SShri Abhyankar   PetscBool      ts_terminate;
2222dc7a7e3SShri Abhyankar 
2232dc7a7e3SShri Abhyankar   PetscFunctionBegin;
2242dc7a7e3SShri Abhyankar   if (event->postevent) {
2256427ac75SLisandro Dalcin     ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
2262dc7a7e3SShri Abhyankar   }
2272dc7a7e3SShri Abhyankar   for (i=0; i < nevents_zero; i++) {
228e105d053SSatish Balay     terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]);
2292dc7a7e3SShri Abhyankar   }
230b2566f29SBarry Smith   ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr);
2319e12be75SShri Abhyankar   if (ts_terminate) {
2322dc7a7e3SShri Abhyankar     ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);
2332dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
2342dc7a7e3SShri Abhyankar   } else {
2352dc7a7e3SShri Abhyankar     event->status = TSEVENT_RESET_NEXTSTEP;
2362dc7a7e3SShri Abhyankar   }
237f7aea88cSShri Abhyankar 
238d0578d90SShri Abhyankar   /* Record the event in the event recorder */
239d0578d90SShri Abhyankar   ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr);
240f7aea88cSShri Abhyankar   ctr = event->recorder.ctr;
2416c4ed002SBarry Smith   if (ctr == MAXEVENTRECORDERS) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS);
242f7aea88cSShri Abhyankar   event->recorder.time[ctr] = t;
243d0578d90SShri Abhyankar   event->recorder.stepnum[ctr] = stepnum;
244f7aea88cSShri Abhyankar   event->recorder.nevents[ctr] = nevents_zero;
245f7aea88cSShri Abhyankar   for (i=0; i < nevents_zero; i++) event->recorder.eventidx[ctr][i] = events_zero[i];
246f7aea88cSShri Abhyankar   event->recorder.ctr++;
247f7aea88cSShri Abhyankar 
24873967516SShri Abhyankar   /* Reset the event residual functions as states might get changed by the postevent callback */
24973967516SShri Abhyankar   event->ptime_prev = t;
250*7dbe0728SLisandro Dalcin   ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);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;
263*7dbe0728SLisandro Dalcin   PetscReal      dt,dt_min;
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);
274*7dbe0728SLisandro Dalcin   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
2752dc7a7e3SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
2762dc7a7e3SShri Abhyankar 
277*7dbe0728SLisandro Dalcin   if (event->status == TSEVENT_NONE) {
278*7dbe0728SLisandro Dalcin     event->timestep_prev = dt;
279*7dbe0728SLisandro Dalcin   }
280*7dbe0728SLisandro Dalcin 
2812dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_RESET_NEXTSTEP) {
282*7dbe0728SLisandro Dalcin     /* Restore previous time step */
283*7dbe0728SLisandro Dalcin     dt = event->timestep_prev;
284*7dbe0728SLisandro Dalcin     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
2852dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
2862dc7a7e3SShri Abhyankar   }
2872dc7a7e3SShri Abhyankar 
2882dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_NONE) {
289*7dbe0728SLisandro Dalcin     event->ptime_end = t;
2902dc7a7e3SShri Abhyankar   }
2912dc7a7e3SShri Abhyankar 
2926427ac75SLisandro Dalcin   ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);
2939e12be75SShri Abhyankar 
2942dc7a7e3SShri Abhyankar   for (i=0; i < event->nevents; i++) {
295e3005195SShri Abhyankar     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
2962dc7a7e3SShri Abhyankar       event->status = TSEVENT_ZERO;
2979e12be75SShri Abhyankar       continue;
298006e6a18SShri Abhyankar     }
299d94325d3SShri Abhyankar     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
300d94325d3SShri Abhyankar     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
301e3005195SShri Abhyankar     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
302d94325d3SShri Abhyankar       switch (event->direction[i]) {
303d94325d3SShri Abhyankar       case -1:
304d94325d3SShri Abhyankar         if (fvalue_sign < 0) {
3052dc7a7e3SShri Abhyankar           rollback = 1;
3062dc7a7e3SShri Abhyankar           /* Compute linearly interpolated new time step */
307e105d053SSatish Balay           dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
3086427ac75SLisandro Dalcin           if (event->monitor) {
3096427ac75SLisandro Dalcin             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
310006e6a18SShri Abhyankar           }
3112dc7a7e3SShri Abhyankar         }
312d94325d3SShri Abhyankar         break;
313d94325d3SShri Abhyankar       case 1:
314d94325d3SShri Abhyankar         if (fvalue_sign > 0) {
315d94325d3SShri Abhyankar           rollback = 1;
316d94325d3SShri Abhyankar           /* Compute linearly interpolated new time step */
317d94325d3SShri Abhyankar           dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
3186427ac75SLisandro Dalcin           if (event->monitor) {
3196427ac75SLisandro Dalcin             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
320006e6a18SShri Abhyankar           }
3212dc7a7e3SShri Abhyankar         }
322d94325d3SShri Abhyankar         break;
323d94325d3SShri Abhyankar       case 0:
324d94325d3SShri Abhyankar         rollback = 1;
325d94325d3SShri Abhyankar         /* Compute linearly interpolated new time step */
326d94325d3SShri Abhyankar         dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
3276427ac75SLisandro Dalcin         if (event->monitor) {
3286427ac75SLisandro Dalcin           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
329006e6a18SShri Abhyankar         }
330d94325d3SShri Abhyankar         break;
331d94325d3SShri Abhyankar       }
332d94325d3SShri Abhyankar     }
333d94325d3SShri Abhyankar   }
334*7dbe0728SLisandro Dalcin 
335d94325d3SShri Abhyankar   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
3362589fa24SLisandro Dalcin   in[0] = event->status; in[1] = rollback;
337*7dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
3382589fa24SLisandro Dalcin   event->status = (TSEventStatus)out[0]; rollback = out[1];
3392589fa24SLisandro Dalcin   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
3402dc7a7e3SShri Abhyankar 
341*7dbe0728SLisandro Dalcin   event->nevents_zero = 0;
3429e12be75SShri Abhyankar   if (event->status == TSEVENT_ZERO) {
3439e12be75SShri Abhyankar     for (i=0; i < event->nevents; i++) {
3449e12be75SShri Abhyankar       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
3459e12be75SShri Abhyankar         event->events_zero[event->nevents_zero++] = i;
3466427ac75SLisandro Dalcin         if (event->monitor) {
3476427ac75SLisandro Dalcin           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g\n",i,(double)t);CHKERRQ(ierr);
3489e12be75SShri Abhyankar         }
3499e12be75SShri Abhyankar       }
3509e12be75SShri Abhyankar     }
3516427ac75SLisandro Dalcin     ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve);CHKERRQ(ierr);
3529e12be75SShri Abhyankar 
353*7dbe0728SLisandro Dalcin     event->ptime_prev = t;
354*7dbe0728SLisandro Dalcin     dt = event->ptime_end - event->ptime_prev;
355*7dbe0728SLisandro Dalcin     if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */
3569e12be75SShri Abhyankar     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
3579e12be75SShri Abhyankar     PetscFunctionReturn(0);
3589e12be75SShri Abhyankar   }
3599e12be75SShri Abhyankar 
3602dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_LOCATED_INTERVAL) {
3612dc7a7e3SShri Abhyankar     ierr = TSRollBack(ts);CHKERRQ(ierr);
362c540466cSShri Abhyankar     ts->steps--;
363c540466cSShri Abhyankar     ts->total_steps--;
3642589fa24SLisandro Dalcin     ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr);
3652dc7a7e3SShri Abhyankar     event->status = TSEVENT_PROCESSING;
3662dc7a7e3SShri Abhyankar   } else {
3672dc7a7e3SShri Abhyankar     for (i=0; i < event->nevents; i++) {
3682dc7a7e3SShri Abhyankar       event->fvalue_prev[i] = event->fvalue[i];
3692dc7a7e3SShri Abhyankar     }
3702dc7a7e3SShri Abhyankar     event->ptime_prev = t;
3712dc7a7e3SShri Abhyankar     if (event->status == TSEVENT_PROCESSING) {
372*7dbe0728SLisandro Dalcin       dt = event->ptime_end - event->ptime_prev;
3732dc7a7e3SShri Abhyankar     }
3742dc7a7e3SShri Abhyankar   }
3759e12be75SShri Abhyankar 
376*7dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
377*7dbe0728SLisandro Dalcin   ierr = TSSetTimeStep(ts,dt_min);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