xref: /petsc/src/ts/event/tsevent.c (revision 4597913a1f019bcfb59a05687b71c0ab204566e7)
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);
177dbe0728SLisandro 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__
277dbe0728SLisandro Dalcin #define __FUNCT__ "TSEventDestroy"
287dbe0728SLisandro Dalcin PetscErrorCode TSEventDestroy(TSEvent *event)
297dbe0728SLisandro Dalcin {
307dbe0728SLisandro Dalcin   PetscErrorCode ierr;
317dbe0728SLisandro Dalcin   PetscInt       i;
327dbe0728SLisandro Dalcin 
337dbe0728SLisandro Dalcin   PetscFunctionBegin;
347dbe0728SLisandro Dalcin   PetscValidPointer(event,1);
357dbe0728SLisandro Dalcin   if (!*event) PetscFunctionReturn(0);
367dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr);
377dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr);
387dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->direction);CHKERRQ(ierr);
397dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->terminate);CHKERRQ(ierr);
407dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr);
417dbe0728SLisandro Dalcin   ierr = PetscFree((*event)->vtol);CHKERRQ(ierr);
427dbe0728SLisandro Dalcin   for (i=0; i < MAXEVENTRECORDERS; i++) {
437dbe0728SLisandro Dalcin     ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr);
447dbe0728SLisandro Dalcin   }
457dbe0728SLisandro Dalcin   ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr);
467dbe0728SLisandro Dalcin   ierr = PetscFree(*event);CHKERRQ(ierr);
477dbe0728SLisandro Dalcin   PetscFunctionReturn(0);
487dbe0728SLisandro Dalcin }
497dbe0728SLisandro Dalcin 
507dbe0728SLisandro 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 
192031fbad4SShri Abhyankar /*
193*4597913aSLisandro Dalcin    Helper rutine to handle user postenvents and recording
194031fbad4SShri Abhyankar */
195031fbad4SShri Abhyankar #undef __FUNCT__
196031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent"
197*4597913aSLisandro Dalcin static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
1982dc7a7e3SShri Abhyankar {
1992dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
2002dc7a7e3SShri Abhyankar   TSEvent        event = ts->event;
2012dc7a7e3SShri Abhyankar   PetscBool      terminate = PETSC_FALSE;
202d0578d90SShri Abhyankar   PetscInt       i,ctr,stepnum;
2032dc7a7e3SShri Abhyankar   PetscBool      ts_terminate;
204*4597913aSLisandro Dalcin   PetscBool      forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
2052dc7a7e3SShri Abhyankar 
2062dc7a7e3SShri Abhyankar   PetscFunctionBegin;
2072dc7a7e3SShri Abhyankar   if (event->postevent) {
208*4597913aSLisandro Dalcin     ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
2092dc7a7e3SShri Abhyankar   }
210*4597913aSLisandro Dalcin 
211*4597913aSLisandro Dalcin   /* Handle termination events */
212*4597913aSLisandro Dalcin   for (i=0; i < event->nevents_zero; i++) {
213*4597913aSLisandro Dalcin     terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]);
2142dc7a7e3SShri Abhyankar   }
215b2566f29SBarry Smith   ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr);
2169e12be75SShri Abhyankar   if (ts_terminate) {
2172dc7a7e3SShri Abhyankar     ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);
2182dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
2192dc7a7e3SShri Abhyankar   } else {
2202dc7a7e3SShri Abhyankar     event->status = TSEVENT_RESET_NEXTSTEP;
2212dc7a7e3SShri Abhyankar   }
222f7aea88cSShri Abhyankar 
223*4597913aSLisandro Dalcin   event->ptime_prev = t;
224*4597913aSLisandro Dalcin   /* Reset event residual functions as states might get changed by the postevent callback */
225*4597913aSLisandro Dalcin   if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);}
226*4597913aSLisandro Dalcin   /* Cache current event residual functions */
227*4597913aSLisandro Dalcin   for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
228*4597913aSLisandro Dalcin 
229d0578d90SShri Abhyankar   /* Record the event in the event recorder */
230d0578d90SShri Abhyankar   ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr);
231f7aea88cSShri Abhyankar   ctr = event->recorder.ctr;
2326c4ed002SBarry Smith   if (ctr == MAXEVENTRECORDERS) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS);
233f7aea88cSShri Abhyankar   event->recorder.time[ctr] = t;
234d0578d90SShri Abhyankar   event->recorder.stepnum[ctr] = stepnum;
235*4597913aSLisandro Dalcin   event->recorder.nevents[ctr] = event->nevents_zero;
236*4597913aSLisandro Dalcin   for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
237f7aea88cSShri Abhyankar   event->recorder.ctr++;
2382dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
2392dc7a7e3SShri Abhyankar }
2402dc7a7e3SShri Abhyankar 
2412dc7a7e3SShri Abhyankar #undef __FUNCT__
2426427ac75SLisandro Dalcin #define __FUNCT__ "TSEventHandler"
2436427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts)
2442dc7a7e3SShri Abhyankar {
2452dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
2466427ac75SLisandro Dalcin   TSEvent        event;
2472dc7a7e3SShri Abhyankar   PetscReal      t;
2482dc7a7e3SShri Abhyankar   Vec            U;
2492dc7a7e3SShri Abhyankar   PetscInt       i;
2507dbe0728SLisandro Dalcin   PetscReal      dt,dt_min;
2512dc7a7e3SShri Abhyankar   PetscInt       rollback=0,in[2],out[2];
2529e12be75SShri Abhyankar   PetscInt       fvalue_sign,fvalueprev_sign;
2532dc7a7e3SShri Abhyankar 
2542dc7a7e3SShri Abhyankar   PetscFunctionBegin;
2556427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2566427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
2576427ac75SLisandro Dalcin   event = ts->event;
2582dc7a7e3SShri Abhyankar 
2592dc7a7e3SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
2607dbe0728SLisandro Dalcin   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
2612dc7a7e3SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
2622dc7a7e3SShri Abhyankar 
2637dbe0728SLisandro Dalcin   if (event->status == TSEVENT_NONE) {
2647dbe0728SLisandro Dalcin     event->timestep_prev = dt;
2657dbe0728SLisandro Dalcin   }
2667dbe0728SLisandro Dalcin 
2672dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_RESET_NEXTSTEP) {
2687dbe0728SLisandro Dalcin     /* Restore previous time step */
2697dbe0728SLisandro Dalcin     dt = event->timestep_prev;
2707dbe0728SLisandro Dalcin     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
2712dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
2722dc7a7e3SShri Abhyankar   }
2732dc7a7e3SShri Abhyankar 
2742dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_NONE) {
2757dbe0728SLisandro Dalcin     event->ptime_end = t;
2762dc7a7e3SShri Abhyankar   }
2772dc7a7e3SShri Abhyankar 
2786427ac75SLisandro Dalcin   ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);
2799e12be75SShri Abhyankar 
2802dc7a7e3SShri Abhyankar   for (i=0; i < event->nevents; i++) {
281e3005195SShri Abhyankar     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
2822dc7a7e3SShri Abhyankar       event->status = TSEVENT_ZERO;
2839e12be75SShri Abhyankar       continue;
284006e6a18SShri Abhyankar     }
285d94325d3SShri Abhyankar     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
286d94325d3SShri Abhyankar     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
287e3005195SShri Abhyankar     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
288d94325d3SShri Abhyankar       switch (event->direction[i]) {
289d94325d3SShri Abhyankar       case -1:
290d94325d3SShri Abhyankar         if (fvalue_sign < 0) {
2912dc7a7e3SShri Abhyankar           rollback = 1;
2922dc7a7e3SShri Abhyankar           /* Compute linearly interpolated new time step */
293e105d053SSatish Balay           dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
2946427ac75SLisandro Dalcin           if (event->monitor) {
2956427ac75SLisandro Dalcin             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
296006e6a18SShri Abhyankar           }
2972dc7a7e3SShri Abhyankar         }
298d94325d3SShri Abhyankar         break;
299d94325d3SShri Abhyankar       case 1:
300d94325d3SShri Abhyankar         if (fvalue_sign > 0) {
301d94325d3SShri Abhyankar           rollback = 1;
302d94325d3SShri Abhyankar           /* Compute linearly interpolated new time step */
303d94325d3SShri Abhyankar           dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
3046427ac75SLisandro Dalcin           if (event->monitor) {
3056427ac75SLisandro Dalcin             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
306006e6a18SShri Abhyankar           }
3072dc7a7e3SShri Abhyankar         }
308d94325d3SShri Abhyankar         break;
309d94325d3SShri Abhyankar       case 0:
310d94325d3SShri Abhyankar         rollback = 1;
311d94325d3SShri Abhyankar         /* Compute linearly interpolated new time step */
312d94325d3SShri Abhyankar         dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
3136427ac75SLisandro Dalcin         if (event->monitor) {
3146427ac75SLisandro Dalcin           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
315006e6a18SShri Abhyankar         }
316d94325d3SShri Abhyankar         break;
317d94325d3SShri Abhyankar       }
318d94325d3SShri Abhyankar     }
319d94325d3SShri Abhyankar   }
3207dbe0728SLisandro Dalcin 
321d94325d3SShri Abhyankar   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
3222589fa24SLisandro Dalcin   in[0] = event->status; in[1] = rollback;
3237dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
3242589fa24SLisandro Dalcin   event->status = (TSEventStatus)out[0]; rollback = out[1];
3252589fa24SLisandro Dalcin   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
3262dc7a7e3SShri Abhyankar 
3277dbe0728SLisandro Dalcin   event->nevents_zero = 0;
3289e12be75SShri Abhyankar   if (event->status == TSEVENT_ZERO) {
3299e12be75SShri Abhyankar     for (i=0; i < event->nevents; i++) {
3309e12be75SShri Abhyankar       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
3319e12be75SShri Abhyankar         event->events_zero[event->nevents_zero++] = i;
3326427ac75SLisandro Dalcin         if (event->monitor) {
3336427ac75SLisandro Dalcin           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g\n",i,(double)t);CHKERRQ(ierr);
3349e12be75SShri Abhyankar         }
3359e12be75SShri Abhyankar       }
3369e12be75SShri Abhyankar     }
337*4597913aSLisandro Dalcin     ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr);
3389e12be75SShri Abhyankar 
3397dbe0728SLisandro Dalcin     dt = event->ptime_end - event->ptime_prev;
3407dbe0728SLisandro Dalcin     if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */
3419e12be75SShri Abhyankar     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
3429e12be75SShri Abhyankar     PetscFunctionReturn(0);
3439e12be75SShri Abhyankar   }
3449e12be75SShri Abhyankar 
3452dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_LOCATED_INTERVAL) {
3462dc7a7e3SShri Abhyankar     ierr = TSRollBack(ts);CHKERRQ(ierr);
347c540466cSShri Abhyankar     ts->steps--;
348c540466cSShri Abhyankar     ts->total_steps--;
3492589fa24SLisandro Dalcin     ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr);
3502dc7a7e3SShri Abhyankar     event->status = TSEVENT_PROCESSING;
3512dc7a7e3SShri Abhyankar   } else {
3522dc7a7e3SShri Abhyankar     for (i=0; i < event->nevents; i++) {
3532dc7a7e3SShri Abhyankar       event->fvalue_prev[i] = event->fvalue[i];
3542dc7a7e3SShri Abhyankar     }
3552dc7a7e3SShri Abhyankar     event->ptime_prev = t;
3562dc7a7e3SShri Abhyankar     if (event->status == TSEVENT_PROCESSING) {
3577dbe0728SLisandro Dalcin       dt = event->ptime_end - event->ptime_prev;
3582dc7a7e3SShri Abhyankar     }
3592dc7a7e3SShri Abhyankar   }
3609e12be75SShri Abhyankar 
3617dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
3627dbe0728SLisandro Dalcin   ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr);
3632dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
3642dc7a7e3SShri Abhyankar }
3652dc7a7e3SShri Abhyankar 
366d0578d90SShri Abhyankar #undef __FUNCT__
3676427ac75SLisandro Dalcin #define __FUNCT__ "TSAdjointEventHandler"
3686427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts)
369d0578d90SShri Abhyankar {
370d0578d90SShri Abhyankar   PetscErrorCode ierr;
3716427ac75SLisandro Dalcin   TSEvent        event;
372d0578d90SShri Abhyankar   PetscReal      t;
373d0578d90SShri Abhyankar   Vec            U;
374d0578d90SShri Abhyankar   PetscInt       ctr;
375d0578d90SShri Abhyankar   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
376d0578d90SShri Abhyankar 
377d0578d90SShri Abhyankar   PetscFunctionBegin;
3786427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3796427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
3806427ac75SLisandro Dalcin   event = ts->event;
381d0578d90SShri Abhyankar 
382d0578d90SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
383d0578d90SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
384d0578d90SShri Abhyankar 
385d0578d90SShri Abhyankar   ctr = event->recorder.ctr-1;
386bcbf8bb3SShri Abhyankar   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
387d0578d90SShri Abhyankar     /* Call the user postevent function */
388d0578d90SShri Abhyankar     if (event->postevent) {
3896427ac75SLisandro Dalcin       ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
390d0578d90SShri Abhyankar       event->recorder.ctr--;
391d0578d90SShri Abhyankar     }
392d0578d90SShri Abhyankar   }
393d0578d90SShri Abhyankar 
394d0578d90SShri Abhyankar   PetscBarrier((PetscObject)ts);
395d0578d90SShri Abhyankar   PetscFunctionReturn(0);
396d0578d90SShri Abhyankar }
397