xref: /petsc/src/ts/event/tsevent.c (revision 1ebf93c6b7d760d592de6ebe6cdc0debaa3caf75)
1 #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "TSEventInitialize"
5 /*
6   TSEventInitialize - Initializes TSEvent for TSSolve
7 */
8 PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U)
9 {
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   if (!event) PetscFunctionReturn(0);
14   PetscValidPointer(event,1);
15   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
16   PetscValidHeaderSpecific(U,VEC_CLASSID,4);
17   ierr = TSGetTimeStep(ts,&event->timestep_prev);CHKERRQ(ierr);
18 
19   event->ptime_prev = t;
20   ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr);
21   /* Initialize the event recorder */
22   event->recorder.ctr = 0;
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "TSEventDestroy"
28 PetscErrorCode TSEventDestroy(TSEvent *event)
29 {
30   PetscErrorCode ierr;
31   PetscInt       i;
32 
33   PetscFunctionBegin;
34   PetscValidPointer(event,1);
35   if (!*event) PetscFunctionReturn(0);
36   ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr);
37   ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr);
38   ierr = PetscFree((*event)->direction);CHKERRQ(ierr);
39   ierr = PetscFree((*event)->terminate);CHKERRQ(ierr);
40   ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr);
41   ierr = PetscFree((*event)->vtol);CHKERRQ(ierr);
42   for (i=0; i < MAXEVENTRECORDERS; i++) {
43     ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr);
44   }
45   ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr);
46   ierr = PetscFree(*event);CHKERRQ(ierr);
47   PetscFunctionReturn(0);
48 }
49 
50 #undef __FUNCT__
51 #define __FUNCT__ "TSSetEventTolerances"
52 /*@
53    TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
54 
55    Logically Collective
56 
57    Input Arguments:
58 +  ts - time integration context
59 .  tol - scalar tolerance, PETSC_DECIDE to leave current value
60 -  vtol - array of tolerances or NULL, used in preference to tol if present
61 
62    Options Database Keys:
63 .  -ts_event_tol <tol> tolerance for event zero crossing
64 
65    Notes:
66    Must call TSSetEventHandler() before setting the tolerances.
67 
68    The size of vtol is equal to the number of events.
69 
70    Level: beginner
71 
72 .seealso: TS, TSEvent, TSSetEventHandler()
73 @*/
74 PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
75 {
76   TSEvent        event;
77   PetscInt       i;
78 
79   PetscFunctionBegin;
80   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
81   if (vtol) PetscValidRealPointer(vtol,3);
82   if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
83 
84   event = ts->event;
85   if (vtol) {
86     for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
87   } else {
88     if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
89       for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
90     }
91   }
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "TSSetEventHandler"
97 /*@C
98    TSSetEventHandler - Sets a monitoring function used for detecting events
99 
100    Logically Collective on TS
101 
102    Input Parameters:
103 +  ts - the TS context obtained from TSCreate()
104 .  nevents - number of local events
105 .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
106                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
107 .  terminate - flag to indicate whether time stepping should be terminated after
108                event is detected (one for each event)
109 .  eventhandler - event monitoring routine
110 .  postevent - [optional] post-event function
111 -  ctx       - [optional] user-defined context for private data for the
112                event monitor and post event routine (use NULL if no
113                context is desired)
114 
115    Calling sequence of eventhandler:
116    PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
117 
118    Input Parameters:
119 +  ts  - the TS context
120 .  t   - current time
121 .  U   - current iterate
122 -  ctx - [optional] context passed with eventhandler
123 
124    Output parameters:
125 .  fvalue    - function value of events at time t
126 
127    Calling sequence of postevent:
128    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero[], PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
129 
130    Input Parameters:
131 +  ts - the TS context
132 .  nevents_zero - number of local events whose event function is zero
133 .  events_zero  - indices of local events which have reached zero
134 .  t            - current time
135 .  U            - current solution
136 .  forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
137 -  ctx          - the context passed with eventhandler
138 
139    Level: intermediate
140 
141 .keywords: TS, event, set, monitor
142 
143 .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
144 @*/
145 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)
146 {
147   PetscErrorCode ierr;
148   TSEvent        event;
149   PetscInt       i;
150   PetscBool      flg;
151   PetscReal      tol=1e-6;
152 
153   PetscFunctionBegin;
154   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
155   PetscValidIntPointer(direction,2);
156   PetscValidIntPointer(terminate,3);
157 
158   ierr = PetscNewLog(ts,&event);CHKERRQ(ierr);
159   ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr);
160   ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr);
161   ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr);
162   ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr);
163   ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr);
164   for (i=0; i < nevents; i++) {
165     event->direction[i] = direction[i];
166     event->terminate[i] = terminate[i];
167   }
168   ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr);
169   event->nevents = nevents;
170   event->eventhandler = eventhandler;
171   event->postevent = postevent;
172   event->ctx = ctx;
173 
174   for (i=0; i < MAXEVENTRECORDERS; i++) {
175     ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr);
176   }
177 
178   ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr);
179   {
180     ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr);
181     ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr);
182   }
183   PetscOptionsEnd();
184   for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
185   if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);}
186 
187   ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr);
188   ts->event = event;
189   PetscFunctionReturn(0);
190 }
191 
192 /*
193    Helper rutine to handle user postenvents and recording
194 */
195 #undef __FUNCT__
196 #define __FUNCT__ "TSPostEvent"
197 static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
198 {
199   PetscErrorCode ierr;
200   TSEvent        event = ts->event;
201   PetscBool      terminate = PETSC_FALSE;
202   PetscInt       i,ctr,stepnum;
203   PetscBool      ts_terminate;
204   PetscBool      forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
205 
206   PetscFunctionBegin;
207   if (event->postevent) {
208     ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
209   }
210 
211   /* Handle termination events */
212   for (i=0; i < event->nevents_zero; i++) {
213     terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]);
214   }
215   ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr);
216   if (ts_terminate) {
217     ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);
218     event->status = TSEVENT_NONE;
219   } else {
220     event->status = TSEVENT_RESET_NEXTSTEP;
221   }
222 
223   event->ptime_prev = t;
224   /* Reset event residual functions as states might get changed by the postevent callback */
225   if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);}
226   /* Cache current event residual functions */
227   for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];
228 
229   /* Record the event in the event recorder */
230   ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr);
231   ctr = event->recorder.ctr;
232   if (ctr == MAXEVENTRECORDERS) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS);
233   event->recorder.time[ctr] = t;
234   event->recorder.stepnum[ctr] = stepnum;
235   event->recorder.nevents[ctr] = event->nevents_zero;
236   for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
237   event->recorder.ctr++;
238   PetscFunctionReturn(0);
239 }
240 
241 #undef __FUNCT__
242 #define __FUNCT__ "TSEventHandler"
243 PetscErrorCode TSEventHandler(TS ts)
244 {
245   PetscErrorCode ierr;
246   TSEvent        event;
247   PetscReal      t;
248   Vec            U;
249   PetscInt       i;
250   PetscReal      dt,dt_min;
251   PetscInt       rollback=0,in[2],out[2];
252   PetscInt       fvalue_sign,fvalueprev_sign;
253 
254   PetscFunctionBegin;
255   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
256   if (!ts->event) PetscFunctionReturn(0);
257   event = ts->event;
258 
259   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
260   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
261   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
262 
263   if (event->status == TSEVENT_NONE) {
264     event->timestep_prev = dt;
265   }
266 
267   if (event->status == TSEVENT_RESET_NEXTSTEP) {
268     /* Restore previous time step */
269     dt = event->timestep_prev;
270     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
271     event->status = TSEVENT_NONE;
272   }
273 
274   if (event->status == TSEVENT_NONE) {
275     event->ptime_end = t;
276   }
277 
278   ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);
279 
280   for (i=0; i < event->nevents; i++) {
281     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
282       event->status = TSEVENT_ZERO;
283       continue;
284     }
285     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
286     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
287     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
288       switch (event->direction[i]) {
289       case -1:
290         if (fvalue_sign < 0) {
291           rollback = 1;
292           /* Compute linearly interpolated new time step */
293           dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
294           if (event->monitor) {
295             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
296           }
297         }
298         break;
299       case 1:
300         if (fvalue_sign > 0) {
301           rollback = 1;
302           /* Compute linearly interpolated new time step */
303           dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
304           if (event->monitor) {
305             ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
306           }
307         }
308         break;
309       case 0:
310         rollback = 1;
311         /* Compute linearly interpolated new time step */
312         dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
313         if (event->monitor) {
314           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr);
315         }
316         break;
317       }
318     }
319   }
320 
321   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
322   in[0] = event->status; in[1] = rollback;
323   ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
324   event->status = (TSEventStatus)out[0]; rollback = out[1];
325   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
326 
327   event->nevents_zero = 0;
328   if (event->status == TSEVENT_ZERO) {
329     for (i=0; i < event->nevents; i++) {
330       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
331         event->events_zero[event->nevents_zero++] = i;
332         if (event->monitor) {
333           ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g\n",i,(double)t);CHKERRQ(ierr);
334         }
335       }
336     }
337     ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr);
338 
339     dt = event->ptime_end - event->ptime_prev;
340     if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */
341     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
342     PetscFunctionReturn(0);
343   }
344 
345   if (event->status == TSEVENT_LOCATED_INTERVAL) {
346     ierr = TSRollBack(ts);CHKERRQ(ierr);
347     ts->steps--;
348     ts->total_steps--;
349     ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr);
350     event->status = TSEVENT_PROCESSING;
351   } else {
352     for (i=0; i < event->nevents; i++) {
353       event->fvalue_prev[i] = event->fvalue[i];
354     }
355     event->ptime_prev = t;
356     if (event->status == TSEVENT_PROCESSING) {
357       dt = event->ptime_end - event->ptime_prev;
358     }
359   }
360 
361   ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr);
362   ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr);
363   PetscFunctionReturn(0);
364 }
365 
366 #undef __FUNCT__
367 #define __FUNCT__ "TSAdjointEventHandler"
368 PetscErrorCode TSAdjointEventHandler(TS ts)
369 {
370   PetscErrorCode ierr;
371   TSEvent        event;
372   PetscReal      t;
373   Vec            U;
374   PetscInt       ctr;
375   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
376 
377   PetscFunctionBegin;
378   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
379   if (!ts->event) PetscFunctionReturn(0);
380   event = ts->event;
381 
382   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
383   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
384 
385   ctr = event->recorder.ctr-1;
386   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
387     /* Call the user postevent function */
388     if (event->postevent) {
389       ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr);
390       event->recorder.ctr--;
391     }
392   }
393 
394   PetscBarrier((PetscObject)ts);
395   PetscFunctionReturn(0);
396 }
397