xref: /petsc/src/ts/event/tsevent.c (revision d569cc17702d2acb6ceff0c04cc0a2e697b140cb)
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);
37e2cdd850SShri Abhyankar   ierr = PetscFree((*event)->fvalue_right);CHKERRQ(ierr);
38e2cdd850SShri Abhyankar   ierr = PetscFree((*event)->zerocrossing);CHKERRQ(ierr);
39e2cdd850SShri 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);
44a4ffd976SShri Abhyankar 
45a4ffd976SShri Abhyankar 
46a4ffd976SShri Abhyankar   for (i=0; i < (*event)->recsize; i++) {
477dbe0728SLisandro Dalcin     ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr);
487dbe0728SLisandro Dalcin   }
49a4ffd976SShri Abhyankar   ierr = PetscFree((*event)->recorder.eventidx);CHKERRQ(ierr);
50a4ffd976SShri Abhyankar   ierr = PetscFree((*event)->recorder.nevents);CHKERRQ(ierr);
51a4ffd976SShri Abhyankar   ierr = PetscFree((*event)->recorder.stepnum);CHKERRQ(ierr);
52a4ffd976SShri Abhyankar   ierr = PetscFree((*event)->recorder.time);CHKERRQ(ierr);
53a4ffd976SShri Abhyankar 
547dbe0728SLisandro Dalcin   ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr);
557dbe0728SLisandro Dalcin   ierr = PetscFree(*event);CHKERRQ(ierr);
567dbe0728SLisandro Dalcin   PetscFunctionReturn(0);
577dbe0728SLisandro Dalcin }
587dbe0728SLisandro Dalcin 
597dbe0728SLisandro Dalcin #undef __FUNCT__
60e3005195SShri Abhyankar #define __FUNCT__ "TSSetEventTolerances"
61e3005195SShri Abhyankar /*@
62e3005195SShri Abhyankar    TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
63e3005195SShri Abhyankar 
64e3005195SShri Abhyankar    Logically Collective
65e3005195SShri Abhyankar 
66e3005195SShri Abhyankar    Input Arguments:
67e3005195SShri Abhyankar +  ts - time integration context
68e3005195SShri Abhyankar .  tol - scalar tolerance, PETSC_DECIDE to leave current value
69e3005195SShri Abhyankar -  vtol - array of tolerances or NULL, used in preference to tol if present
70e3005195SShri Abhyankar 
712589fa24SLisandro Dalcin    Options Database Keys:
722589fa24SLisandro Dalcin .  -ts_event_tol <tol> tolerance for event zero crossing
73e3005195SShri Abhyankar 
74e3005195SShri Abhyankar    Notes:
75f25fe674SLisandro Dalcin    Must call TSSetEventHandler() before setting the tolerances.
76e3005195SShri Abhyankar 
77e3005195SShri Abhyankar    The size of vtol is equal to the number of events.
78e3005195SShri Abhyankar 
79e3005195SShri Abhyankar    Level: beginner
80e3005195SShri Abhyankar 
81f25fe674SLisandro Dalcin .seealso: TS, TSEvent, TSSetEventHandler()
82e3005195SShri Abhyankar @*/
836427ac75SLisandro Dalcin PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
84e3005195SShri Abhyankar {
856427ac75SLisandro Dalcin   TSEvent        event;
86e3005195SShri Abhyankar   PetscInt       i;
87e3005195SShri Abhyankar 
88e3005195SShri Abhyankar   PetscFunctionBegin;
896427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
906427ac75SLisandro Dalcin   if (vtol) PetscValidRealPointer(vtol,3);
91f25fe674SLisandro Dalcin   if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
926427ac75SLisandro Dalcin 
936427ac75SLisandro Dalcin   event = ts->event;
94e3005195SShri Abhyankar   if (vtol) {
95e3005195SShri Abhyankar     for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
96e3005195SShri Abhyankar   } else {
97e3005195SShri Abhyankar     if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
98e3005195SShri Abhyankar       for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
99e3005195SShri Abhyankar     }
100e3005195SShri Abhyankar   }
101e3005195SShri Abhyankar   PetscFunctionReturn(0);
102e3005195SShri Abhyankar }
103e3005195SShri Abhyankar 
104e3005195SShri Abhyankar #undef __FUNCT__
105f25fe674SLisandro Dalcin #define __FUNCT__ "TSSetEventHandler"
1062dc7a7e3SShri Abhyankar /*@C
107f25fe674SLisandro Dalcin    TSSetEventHandler - Sets a monitoring function used for detecting events
1082dc7a7e3SShri Abhyankar 
1092dc7a7e3SShri Abhyankar    Logically Collective on TS
1102dc7a7e3SShri Abhyankar 
1112dc7a7e3SShri Abhyankar    Input Parameters:
1122dc7a7e3SShri Abhyankar +  ts - the TS context obtained from TSCreate()
1132dc7a7e3SShri Abhyankar .  nevents - number of local events
114d94325d3SShri Abhyankar .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
115d94325d3SShri Abhyankar                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
116d94325d3SShri Abhyankar .  terminate - flag to indicate whether time stepping should be terminated after
117d94325d3SShri Abhyankar                event is detected (one for each event)
1186427ac75SLisandro Dalcin .  eventhandler - event monitoring routine
1192dc7a7e3SShri Abhyankar .  postevent - [optional] post-event function
1202589fa24SLisandro Dalcin -  ctx       - [optional] user-defined context for private data for the
1212dc7a7e3SShri Abhyankar                event monitor and post event routine (use NULL if no
1222dc7a7e3SShri Abhyankar                context is desired)
1232dc7a7e3SShri Abhyankar 
1246427ac75SLisandro Dalcin    Calling sequence of eventhandler:
1252589fa24SLisandro Dalcin    PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
1262dc7a7e3SShri Abhyankar 
1272dc7a7e3SShri Abhyankar    Input Parameters:
1282dc7a7e3SShri Abhyankar +  ts  - the TS context
1292dc7a7e3SShri Abhyankar .  t   - current time
1302dc7a7e3SShri Abhyankar .  U   - current iterate
1316427ac75SLisandro Dalcin -  ctx - [optional] context passed with eventhandler
1322dc7a7e3SShri Abhyankar 
1332dc7a7e3SShri Abhyankar    Output parameters:
134d94325d3SShri Abhyankar .  fvalue    - function value of events at time t
1352dc7a7e3SShri Abhyankar 
1362dc7a7e3SShri Abhyankar    Calling sequence of postevent:
1372589fa24SLisandro Dalcin    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero[], PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
1382dc7a7e3SShri Abhyankar 
1392dc7a7e3SShri Abhyankar    Input Parameters:
1402dc7a7e3SShri Abhyankar +  ts - the TS context
1412dc7a7e3SShri Abhyankar .  nevents_zero - number of local events whose event function is zero
1422dc7a7e3SShri Abhyankar .  events_zero  - indices of local events which have reached zero
1432dc7a7e3SShri Abhyankar .  t            - current time
1442dc7a7e3SShri Abhyankar .  U            - current solution
145031fbad4SShri Abhyankar .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
1466427ac75SLisandro Dalcin -  ctx          - the context passed with eventhandler
1472dc7a7e3SShri Abhyankar 
1482dc7a7e3SShri Abhyankar    Level: intermediate
1492dc7a7e3SShri Abhyankar 
1502dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor
1512dc7a7e3SShri Abhyankar 
1522dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
1532dc7a7e3SShri Abhyankar @*/
1542589fa24SLisandro 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)
1552dc7a7e3SShri Abhyankar {
1562dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
1572dc7a7e3SShri Abhyankar   TSEvent        event;
158d94325d3SShri Abhyankar   PetscInt       i;
159006e6a18SShri Abhyankar   PetscBool      flg;
160a6c783d2SShri Abhyankar #if defined PETSC_USE_REAL_SINGLE
161a6c783d2SShri Abhyankar   PetscReal      tol=1e-4;
162a6c783d2SShri Abhyankar #else
163*d569cc17SSatish Balay   PetscReal      tol=1e-6;
164a6c783d2SShri Abhyankar #endif
1652dc7a7e3SShri Abhyankar 
1662dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1676427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1686427ac75SLisandro Dalcin   PetscValidIntPointer(direction,2);
1696427ac75SLisandro Dalcin   PetscValidIntPointer(terminate,3);
1706427ac75SLisandro Dalcin 
1716427ac75SLisandro Dalcin   ierr = PetscNewLog(ts,&event);CHKERRQ(ierr);
172854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr);
173854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr);
174e2cdd850SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->fvalue_right);CHKERRQ(ierr);
175e2cdd850SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->zerocrossing);CHKERRQ(ierr);
176e2cdd850SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->side);CHKERRQ(ierr);
177854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr);
178854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr);
179e3005195SShri Abhyankar   ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr);
180d94325d3SShri Abhyankar   for (i=0; i < nevents; i++) {
181d94325d3SShri Abhyankar     event->direction[i] = direction[i];
182d94325d3SShri Abhyankar     event->terminate[i] = terminate[i];
183e2cdd850SShri Abhyankar     event->zerocrossing[i] = PETSC_FALSE;
184e2cdd850SShri Abhyankar     event->side[i] = 0;
185d94325d3SShri Abhyankar   }
186854ce69bSBarry Smith   ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr);
1872589fa24SLisandro Dalcin   event->nevents = nevents;
1886427ac75SLisandro Dalcin   event->eventhandler = eventhandler;
1892dc7a7e3SShri Abhyankar   event->postevent = postevent;
1906427ac75SLisandro Dalcin   event->ctx = ctx;
1912dc7a7e3SShri Abhyankar 
192a4ffd976SShri Abhyankar   event->recsize = 12;  /* Initial size of the recorder */
193a9514d71SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr);
1942dc7a7e3SShri Abhyankar   {
195e3005195SShri Abhyankar     ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr);
196006e6a18SShri Abhyankar     ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr);
197a4ffd976SShri Abhyankar     ierr = PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);CHKERRQ(ierr);
1982dc7a7e3SShri Abhyankar   }
1999e12be75SShri Abhyankar   PetscOptionsEnd();
200a4ffd976SShri Abhyankar 
201a4ffd976SShri Abhyankar   ierr = PetscMalloc1(event->recsize,&event->recorder.time);CHKERRQ(ierr);
202a4ffd976SShri Abhyankar   ierr = PetscMalloc1(event->recsize,&event->recorder.stepnum);CHKERRQ(ierr);
203a4ffd976SShri Abhyankar   ierr = PetscMalloc1(event->recsize,&event->recorder.nevents);CHKERRQ(ierr);
204a4ffd976SShri Abhyankar   ierr = PetscMalloc1(event->recsize,&event->recorder.eventidx);CHKERRQ(ierr);
205a4ffd976SShri Abhyankar   for (i=0; i < event->recsize; i++) {
206a4ffd976SShri Abhyankar     ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr);
207a4ffd976SShri Abhyankar   }
208a4ffd976SShri Abhyankar 
209e3005195SShri Abhyankar   for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
2106427ac75SLisandro Dalcin   if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);}
2119e12be75SShri Abhyankar 
2126427ac75SLisandro Dalcin   ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr);
213d94325d3SShri Abhyankar   ts->event = event;
2142dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
2152dc7a7e3SShri Abhyankar }
2162dc7a7e3SShri Abhyankar 
217a4ffd976SShri Abhyankar #undef __FUNCT__
218a4ffd976SShri Abhyankar #define __FUNCT__ "TSEventRecorderResize"
219a4ffd976SShri Abhyankar /*
220a4ffd976SShri Abhyankar   TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
221a4ffd976SShri Abhyankar                           is reached.
222a4ffd976SShri Abhyankar */
223a4ffd976SShri Abhyankar PetscErrorCode TSEventRecorderResize(TSEvent event)
224a4ffd976SShri Abhyankar {
225a4ffd976SShri Abhyankar   PetscErrorCode ierr;
226a4ffd976SShri Abhyankar   PetscReal      *time;
227a4ffd976SShri Abhyankar   PetscInt       *stepnum;
228a4ffd976SShri Abhyankar   PetscInt       *nevents;
229a4ffd976SShri Abhyankar   PetscInt       **eventidx;
230a4ffd976SShri Abhyankar   PetscInt       i,fact=2;
231a4ffd976SShri Abhyankar 
232a4ffd976SShri Abhyankar   PetscFunctionBegin;
233a4ffd976SShri Abhyankar 
234a4ffd976SShri Abhyankar   /* Create large arrays */
235a4ffd976SShri Abhyankar   ierr = PetscMalloc1(fact*event->recsize,&time);CHKERRQ(ierr);
236a4ffd976SShri Abhyankar   ierr = PetscMalloc1(fact*event->recsize,&stepnum);CHKERRQ(ierr);
237a4ffd976SShri Abhyankar   ierr = PetscMalloc1(fact*event->recsize,&nevents);CHKERRQ(ierr);
238a4ffd976SShri Abhyankar   ierr = PetscMalloc1(fact*event->recsize,&eventidx);CHKERRQ(ierr);
239a4ffd976SShri Abhyankar   for (i=0; i < fact*event->recsize; i++) {
240a4ffd976SShri Abhyankar     ierr = PetscMalloc1(event->nevents,&eventidx[i]);CHKERRQ(ierr);
241a4ffd976SShri Abhyankar   }
242a4ffd976SShri Abhyankar 
243a4ffd976SShri Abhyankar   /* Copy over data */
244a4ffd976SShri Abhyankar   ierr = PetscMemcpy(time,event->recorder.time,event->recsize*sizeof(PetscReal));CHKERRQ(ierr);
245a4ffd976SShri Abhyankar   ierr = PetscMemcpy(stepnum,event->recorder.stepnum,event->recsize*sizeof(PetscInt));CHKERRQ(ierr);
246a4ffd976SShri Abhyankar   ierr = PetscMemcpy(nevents,event->recorder.nevents,event->recsize*sizeof(PetscInt));CHKERRQ(ierr);
247a4ffd976SShri Abhyankar   for(i=0; i < event->recsize; i++) {
248a4ffd976SShri Abhyankar     ierr = PetscMemcpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]*sizeof(PetscInt));CHKERRQ(ierr);
249a4ffd976SShri Abhyankar   }
250a4ffd976SShri Abhyankar 
251a4ffd976SShri Abhyankar   /* Destroy old arrays */
252a4ffd976SShri Abhyankar   for (i=0; i < event->recsize; i++) {
253a4ffd976SShri Abhyankar     ierr = PetscFree(event->recorder.eventidx[i]);CHKERRQ(ierr);
254a4ffd976SShri Abhyankar   }
255a4ffd976SShri Abhyankar   ierr = PetscFree(event->recorder.eventidx);CHKERRQ(ierr);
256a4ffd976SShri Abhyankar   ierr = PetscFree(event->recorder.nevents);CHKERRQ(ierr);
257a4ffd976SShri Abhyankar   ierr = PetscFree(event->recorder.stepnum);CHKERRQ(ierr);
258a4ffd976SShri Abhyankar   ierr = PetscFree(event->recorder.time);CHKERRQ(ierr);
259a4ffd976SShri Abhyankar 
260a4ffd976SShri Abhyankar   /* Set pointers */
261a4ffd976SShri Abhyankar   event->recorder.time = time;
262a4ffd976SShri Abhyankar   event->recorder.stepnum = stepnum;
263a4ffd976SShri Abhyankar   event->recorder.nevents = nevents;
264a4ffd976SShri Abhyankar   event->recorder.eventidx = eventidx;
265a4ffd976SShri Abhyankar 
266a4ffd976SShri Abhyankar   /* Double size */
267a4ffd976SShri Abhyankar   event->recsize *= fact;
268a4ffd976SShri Abhyankar 
269a4ffd976SShri Abhyankar   PetscFunctionReturn(0);
270a4ffd976SShri Abhyankar }
271a4ffd976SShri Abhyankar 
272031fbad4SShri Abhyankar /*
2734597913aSLisandro Dalcin    Helper rutine to handle user postenvents and recording
274031fbad4SShri Abhyankar */
275031fbad4SShri Abhyankar #undef __FUNCT__
276031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent"
2774597913aSLisandro Dalcin static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
2782dc7a7e3SShri Abhyankar {
2792dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
2802dc7a7e3SShri Abhyankar   TSEvent        event = ts->event;
2812dc7a7e3SShri Abhyankar   PetscBool      terminate = PETSC_FALSE;
282d0578d90SShri Abhyankar   PetscInt       i,ctr,stepnum;
2832dc7a7e3SShri Abhyankar   PetscBool      ts_terminate;
2844597913aSLisandro Dalcin   PetscBool      forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
2852dc7a7e3SShri Abhyankar 
2862dc7a7e3SShri Abhyankar   PetscFunctionBegin;
2872dc7a7e3SShri Abhyankar   if (event->postevent) {
2884597913aSLisandro Dalcin     ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
2892dc7a7e3SShri Abhyankar   }
2904597913aSLisandro Dalcin 
2914597913aSLisandro Dalcin   /* Handle termination events */
2924597913aSLisandro Dalcin   for (i=0; i < event->nevents_zero; i++) {
2934597913aSLisandro Dalcin     terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]);
2942dc7a7e3SShri Abhyankar   }
295b2566f29SBarry Smith   ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr);
2969e12be75SShri Abhyankar   if (ts_terminate) {
2972dc7a7e3SShri Abhyankar     ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);
2982dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
2992dc7a7e3SShri Abhyankar   } else {
3002dc7a7e3SShri Abhyankar     event->status = TSEVENT_RESET_NEXTSTEP;
3012dc7a7e3SShri Abhyankar   }
302f7aea88cSShri Abhyankar 
3034597913aSLisandro Dalcin   event->ptime_prev = t;
3044597913aSLisandro Dalcin   /* Reset event residual functions as states might get changed by the postevent callback */
3054597913aSLisandro Dalcin   if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);}
3064597913aSLisandro Dalcin   /* Cache current event residual functions */
3074597913aSLisandro Dalcin   for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
3084597913aSLisandro Dalcin 
309d0578d90SShri Abhyankar   /* Record the event in the event recorder */
310d0578d90SShri Abhyankar   ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr);
311f7aea88cSShri Abhyankar   ctr = event->recorder.ctr;
312a4ffd976SShri Abhyankar   if (ctr == event->recsize) {
313a4ffd976SShri Abhyankar     ierr = TSEventRecorderResize(event);CHKERRQ(ierr);
314a4ffd976SShri Abhyankar   }
315f7aea88cSShri Abhyankar   event->recorder.time[ctr] = t;
316d0578d90SShri Abhyankar   event->recorder.stepnum[ctr] = stepnum;
3174597913aSLisandro Dalcin   event->recorder.nevents[ctr] = event->nevents_zero;
3184597913aSLisandro Dalcin   for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
319f7aea88cSShri Abhyankar   event->recorder.ctr++;
3202dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
3212dc7a7e3SShri Abhyankar }
3222dc7a7e3SShri Abhyankar 
323e2cdd850SShri Abhyankar /* Uses Anderson-Bjorck variant of regular falsi method */
324a3a8645aSShri Abhyankar PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t, PetscReal tright, PetscScalar fleft, PetscScalar f, PetscScalar fright, PetscInt side,PetscReal dt)
32538bf2713SShri Abhyankar {
326e2cdd850SShri Abhyankar   PetscReal scal=1.0;
327e2cdd850SShri Abhyankar   if(PetscRealPart(fleft)*PetscRealPart(f) < 0) {
328e2cdd850SShri Abhyankar     if(side == 1) {
329a6c783d2SShri Abhyankar       scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
330a6c783d2SShri Abhyankar       if(scal < PETSC_SMALL) scal = 0.5;
331e2cdd850SShri Abhyankar     }
332a3a8645aSShri Abhyankar     dt = PetscMin(dt,(scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft);
333e2cdd850SShri Abhyankar   } else {
334e2cdd850SShri Abhyankar     if(side == -1) {
335a6c783d2SShri Abhyankar       scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
336a6c783d2SShri Abhyankar       if(scal < PETSC_SMALL) scal = 0.5;
337e2cdd850SShri Abhyankar     }
338a3a8645aSShri Abhyankar     dt = PetscMin(dt,(PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t);
339e2cdd850SShri Abhyankar   }
34038bf2713SShri Abhyankar   return dt;
34138bf2713SShri Abhyankar }
342e2cdd850SShri Abhyankar 
34338bf2713SShri Abhyankar 
3442dc7a7e3SShri Abhyankar #undef __FUNCT__
3456427ac75SLisandro Dalcin #define __FUNCT__ "TSEventHandler"
3466427ac75SLisandro Dalcin PetscErrorCode TSEventHandler(TS ts)
3472dc7a7e3SShri Abhyankar {
3482dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
3496427ac75SLisandro Dalcin   TSEvent        event;
3502dc7a7e3SShri Abhyankar   PetscReal      t;
3512dc7a7e3SShri Abhyankar   Vec            U;
3522dc7a7e3SShri Abhyankar   PetscInt       i;
3537dbe0728SLisandro Dalcin   PetscReal      dt,dt_min;
3542dc7a7e3SShri Abhyankar   PetscInt       rollback=0,in[2],out[2];
3559e12be75SShri Abhyankar   PetscInt       fvalue_sign,fvalueprev_sign;
3562dc7a7e3SShri Abhyankar 
3572dc7a7e3SShri Abhyankar   PetscFunctionBegin;
3586427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
3596427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
3606427ac75SLisandro Dalcin   event = ts->event;
3612dc7a7e3SShri Abhyankar 
3622dc7a7e3SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
3637dbe0728SLisandro Dalcin   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
3642dc7a7e3SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
3652dc7a7e3SShri Abhyankar 
3667dbe0728SLisandro Dalcin   if (event->status == TSEVENT_NONE) {
3677dbe0728SLisandro Dalcin     event->timestep_prev = dt;
3687dbe0728SLisandro Dalcin   }
3697dbe0728SLisandro Dalcin 
3702dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_RESET_NEXTSTEP) {
3717dbe0728SLisandro Dalcin     /* Restore previous time step */
3727dbe0728SLisandro Dalcin     dt = event->timestep_prev;
3737dbe0728SLisandro Dalcin     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
3742dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
3752dc7a7e3SShri Abhyankar   }
3762dc7a7e3SShri Abhyankar 
3772dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_NONE) {
3787dbe0728SLisandro Dalcin     event->ptime_end = t;
3792dc7a7e3SShri Abhyankar   }
3802dc7a7e3SShri Abhyankar 
3816427ac75SLisandro Dalcin   ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);
3829e12be75SShri Abhyankar 
3832dc7a7e3SShri Abhyankar   for (i=0; i < event->nevents; i++) {
384e3005195SShri Abhyankar     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
3852dc7a7e3SShri Abhyankar       event->status = TSEVENT_ZERO;
386e2cdd850SShri Abhyankar       event->fvalue_right[i] = event->fvalue[i];
3879e12be75SShri Abhyankar       continue;
388006e6a18SShri Abhyankar     }
389d94325d3SShri Abhyankar     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
390d94325d3SShri Abhyankar     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
391e3005195SShri Abhyankar     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
392d94325d3SShri Abhyankar       switch (event->direction[i]) {
393d94325d3SShri Abhyankar       case -1:
394d94325d3SShri Abhyankar         if (fvalue_sign < 0) {
3952dc7a7e3SShri Abhyankar 	  rollback = 1;
396e2cdd850SShri Abhyankar 
397e2cdd850SShri Abhyankar 	  /* Compute new time step */
398e2cdd850SShri 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);
399e2cdd850SShri Abhyankar 
4006427ac75SLisandro Dalcin 	  if (event->monitor) {
40138bf2713SShri 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);
402006e6a18SShri Abhyankar 	  }
403e2cdd850SShri Abhyankar 	  event->fvalue_right[i] = event->fvalue[i];
404e2cdd850SShri Abhyankar 	  event->side[i] = 1;
405e2cdd850SShri Abhyankar 
406e2cdd850SShri Abhyankar 	  if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
407e2cdd850SShri Abhyankar 	  event->status = TSEVENT_LOCATED_INTERVAL;
4082dc7a7e3SShri Abhyankar         }
409d94325d3SShri Abhyankar         break;
410d94325d3SShri Abhyankar       case 1:
411d94325d3SShri Abhyankar         if (fvalue_sign > 0) {
412d94325d3SShri Abhyankar 	  rollback = 1;
413e2cdd850SShri Abhyankar 
414e2cdd850SShri Abhyankar  	  /* Compute new time step */
415e2cdd850SShri 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);
416e2cdd850SShri Abhyankar 
4176427ac75SLisandro Dalcin 	  if (event->monitor) {
41838bf2713SShri 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);
419006e6a18SShri Abhyankar 	  }
420e2cdd850SShri Abhyankar 	  event->fvalue_right[i] = event->fvalue[i];
421e2cdd850SShri Abhyankar 	  event->side[i] = 1;
422e2cdd850SShri Abhyankar 
423e2cdd850SShri Abhyankar      	  if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
424e2cdd850SShri Abhyankar 	  event->status = TSEVENT_LOCATED_INTERVAL;
4252dc7a7e3SShri Abhyankar         }
426d94325d3SShri Abhyankar         break;
427d94325d3SShri Abhyankar       case 0:
428d94325d3SShri Abhyankar 	rollback = 1;
429e2cdd850SShri Abhyankar 
430e2cdd850SShri Abhyankar 	/* Compute new time step */
431e2cdd850SShri 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);
432e2cdd850SShri Abhyankar 
4336427ac75SLisandro Dalcin 	if (event->monitor) {
43438bf2713SShri 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);
435006e6a18SShri Abhyankar 	}
436e2cdd850SShri Abhyankar 	event->fvalue_right[i] = event->fvalue[i];
437e2cdd850SShri Abhyankar 	event->side[i] = 1;
438e2cdd850SShri Abhyankar 
439e2cdd850SShri Abhyankar 	if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
440e2cdd850SShri Abhyankar 	event->status = TSEVENT_LOCATED_INTERVAL;
441e2cdd850SShri Abhyankar 
442d94325d3SShri Abhyankar         break;
443d94325d3SShri Abhyankar       }
444d94325d3SShri Abhyankar     }
445d94325d3SShri Abhyankar   }
4467dbe0728SLisandro Dalcin 
4472589fa24SLisandro Dalcin   in[0] = event->status; in[1] = rollback;
4487dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
4492589fa24SLisandro Dalcin   event->status = (TSEventStatus)out[0]; rollback = out[1];
4502589fa24SLisandro Dalcin   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
4512dc7a7e3SShri Abhyankar 
4527dbe0728SLisandro Dalcin   event->nevents_zero = 0;
4539e12be75SShri Abhyankar   if (event->status == TSEVENT_ZERO) {
4549e12be75SShri Abhyankar     for (i=0; i < event->nevents; i++) {
4559e12be75SShri Abhyankar       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
4569e12be75SShri Abhyankar         event->events_zero[event->nevents_zero++] = i;
4576427ac75SLisandro Dalcin         if (event->monitor) {
45838bf2713SShri 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);
4599e12be75SShri Abhyankar         }
460e2cdd850SShri Abhyankar 	event->zerocrossing[i] = PETSC_FALSE;
4619e12be75SShri Abhyankar       }
462e2cdd850SShri Abhyankar       event->side[i] = 0;
4639e12be75SShri Abhyankar     }
4644597913aSLisandro Dalcin     ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr);
4659e12be75SShri Abhyankar 
466e2cdd850SShri Abhyankar     dt = event->ptime_end - t;
4677dbe0728SLisandro Dalcin     if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */
4689e12be75SShri Abhyankar     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
46938bf2713SShri Abhyankar     event->iterctr = 0;
4709e12be75SShri Abhyankar     PetscFunctionReturn(0);
4719e12be75SShri Abhyankar   }
4729e12be75SShri Abhyankar 
4732dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_LOCATED_INTERVAL) {
4742dc7a7e3SShri Abhyankar     ierr = TSRollBack(ts);CHKERRQ(ierr);
475c540466cSShri Abhyankar     ts->steps--;
476c540466cSShri Abhyankar     ts->total_steps--;
4772589fa24SLisandro Dalcin     ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr);
4782dc7a7e3SShri Abhyankar     event->status = TSEVENT_PROCESSING;
479e2cdd850SShri Abhyankar 
480e2cdd850SShri Abhyankar     event->ptime_right = t;
4812dc7a7e3SShri Abhyankar   } else {
4822dc7a7e3SShri Abhyankar     for(i=0; i < event->nevents; i++) {
483e2cdd850SShri Abhyankar       if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) {
484e2cdd850SShri Abhyankar 	/* Compute new time step */
485e2cdd850SShri 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);
486e2cdd850SShri Abhyankar 
487e2cdd850SShri Abhyankar 	event->side[i] = -1;
488e2cdd850SShri Abhyankar       }
4892dc7a7e3SShri Abhyankar       event->fvalue_prev[i] = event->fvalue[i];
4902dc7a7e3SShri Abhyankar     }
49138bf2713SShri Abhyankar 
492e2cdd850SShri Abhyankar     if (event->monitor && event->status == TSEVENT_PROCESSING) {
49338bf2713SShri 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);
4942dc7a7e3SShri Abhyankar     }
49538bf2713SShri Abhyankar     event->ptime_prev = t;
49638bf2713SShri Abhyankar   }
49738bf2713SShri Abhyankar 
49838bf2713SShri Abhyankar   if(event->status == TSEVENT_PROCESSING) event->iterctr++;
4999e12be75SShri Abhyankar 
5007dbe0728SLisandro Dalcin   ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
5017dbe0728SLisandro Dalcin   ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr);
5022dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
5032dc7a7e3SShri Abhyankar }
5042dc7a7e3SShri Abhyankar 
505d0578d90SShri Abhyankar #undef __FUNCT__
5066427ac75SLisandro Dalcin #define __FUNCT__ "TSAdjointEventHandler"
5076427ac75SLisandro Dalcin PetscErrorCode TSAdjointEventHandler(TS ts)
508d0578d90SShri Abhyankar {
509d0578d90SShri Abhyankar   PetscErrorCode ierr;
5106427ac75SLisandro Dalcin   TSEvent        event;
511d0578d90SShri Abhyankar   PetscReal      t;
512d0578d90SShri Abhyankar   Vec            U;
513d0578d90SShri Abhyankar   PetscInt       ctr;
514d0578d90SShri Abhyankar   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
515d0578d90SShri Abhyankar 
516d0578d90SShri Abhyankar   PetscFunctionBegin;
5176427ac75SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5186427ac75SLisandro Dalcin   if (!ts->event) PetscFunctionReturn(0);
5196427ac75SLisandro Dalcin   event = ts->event;
520d0578d90SShri Abhyankar 
521d0578d90SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
522d0578d90SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
523d0578d90SShri Abhyankar 
524d0578d90SShri Abhyankar   ctr = event->recorder.ctr-1;
525bcbf8bb3SShri Abhyankar   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
526d0578d90SShri Abhyankar     /* Call the user postevent function */
527d0578d90SShri Abhyankar     if (event->postevent) {
5286427ac75SLisandro Dalcin       ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
529d0578d90SShri Abhyankar       event->recorder.ctr--;
530d0578d90SShri Abhyankar     }
531d0578d90SShri Abhyankar   }
532d0578d90SShri Abhyankar 
533d0578d90SShri Abhyankar   PetscBarrier((PetscObject)ts);
534d0578d90SShri Abhyankar   PetscFunctionReturn(0);
535d0578d90SShri Abhyankar }
536