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