Actual source code: tsevent.c
petsc-3.11.4 2019-09-28
1: #include <petsc/private/tsimpl.h>
3: /*
4: TSEventInitialize - Initializes TSEvent for TSSolve
5: */
6: PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U)
7: {
11: if (!event) return(0);
15: event->ptime_prev = t;
16: event->iterctr = 0;
17: (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);
18: return(0);
19: }
21: PetscErrorCode TSEventDestroy(TSEvent *event)
22: {
24: PetscInt i;
28: if (!*event) return(0);
29: if (--(*event)->refct > 0) {*event = 0; return(0);}
31: PetscFree((*event)->fvalue);
32: PetscFree((*event)->fvalue_prev);
33: PetscFree((*event)->fvalue_right);
34: PetscFree((*event)->zerocrossing);
35: PetscFree((*event)->side);
36: PetscFree((*event)->direction);
37: PetscFree((*event)->terminate);
38: PetscFree((*event)->events_zero);
39: PetscFree((*event)->vtol);
41: for (i=0; i < (*event)->recsize; i++) {
42: PetscFree((*event)->recorder.eventidx[i]);
43: }
44: PetscFree((*event)->recorder.eventidx);
45: PetscFree((*event)->recorder.nevents);
46: PetscFree((*event)->recorder.stepnum);
47: PetscFree((*event)->recorder.time);
49: PetscViewerDestroy(&(*event)->monitor);
50: PetscFree(*event);
51: return(0);
52: }
54: /*@
55: TSSetPostEventIntervalStep - Set the time-step immediately following the event interval
57: Logically Collective
59: Input Arguments:
60: + ts - time integration context
61: - dt - post event interval step
63: Options Database Keys:
64: . -ts_event_post_eventinterval_step <dt> time-step after event interval
66: Notes:
67: TSSetPostEventIntervalStep allows to set a time-step immediately following an event interval.
68:
69: This function should be called from the postevent function set with TSSetEventHandler().
71: The post event interval time-step should be selected based on the dynamics following the event.
72: If the dynamics are stiff, a conservative (small) step should be used.
73: If not, then a larger time-step can be used.
74:
75: Level: Advanced
76: .seealso: TS, TSEvent, TSSetEventHandler()
77: @*/
78: PetscErrorCode TSSetPostEventIntervalStep(TS ts,PetscReal dt)
79: {
81: ts->event->timestep_posteventinterval = dt;
82: return(0);
83: }
85: /*@
86: TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
88: Logically Collective
90: Input Arguments:
91: + ts - time integration context
92: . tol - scalar tolerance, PETSC_DECIDE to leave current value
93: - vtol - array of tolerances or NULL, used in preference to tol if present
95: Options Database Keys:
96: . -ts_event_tol <tol> tolerance for event zero crossing
98: Notes:
99: Must call TSSetEventHandler() before setting the tolerances.
101: The size of vtol is equal to the number of events.
103: Level: beginner
105: .seealso: TS, TSEvent, TSSetEventHandler()
106: @*/
107: PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
108: {
109: TSEvent event;
110: PetscInt i;
115: if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
117: event = ts->event;
118: if (vtol) {
119: for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
120: } else {
121: if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
122: for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
123: }
124: }
125: return(0);
126: }
128: /*@C
129: TSSetEventHandler - Sets a monitoring function used for detecting events
131: Logically Collective on TS
133: Input Parameters:
134: + ts - the TS context obtained from TSCreate()
135: . nevents - number of local events
136: . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
137: +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
138: . terminate - flag to indicate whether time stepping should be terminated after
139: event is detected (one for each event)
140: . eventhandler - event monitoring routine
141: . postevent - [optional] post-event function
142: - ctx - [optional] user-defined context for private data for the
143: event monitor and post event routine (use NULL if no
144: context is desired)
146: Calling sequence of eventhandler:
147: PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
149: Input Parameters:
150: + ts - the TS context
151: . t - current time
152: . U - current iterate
153: - ctx - [optional] context passed with eventhandler
155: Output parameters:
156: . fvalue - function value of events at time t
158: Calling sequence of postevent:
159: PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
161: Input Parameters:
162: + ts - the TS context
163: . nevents_zero - number of local events whose event function is zero
164: . events_zero - indices of local events which have reached zero
165: . t - current time
166: . U - current solution
167: . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0)
168: - ctx - the context passed with eventhandler
170: Level: intermediate
172: .keywords: TS, event, set
174: .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
175: @*/
176: 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)
177: {
179: TSEvent event;
180: PetscInt i;
181: PetscBool flg;
182: #if defined PETSC_USE_REAL_SINGLE
183: PetscReal tol=1e-4;
184: #else
185: PetscReal tol=1e-6;
186: #endif
190: if (nevents) {
193: }
195: PetscNewLog(ts,&event);
196: PetscMalloc1(nevents,&event->fvalue);
197: PetscMalloc1(nevents,&event->fvalue_prev);
198: PetscMalloc1(nevents,&event->fvalue_right);
199: PetscMalloc1(nevents,&event->zerocrossing);
200: PetscMalloc1(nevents,&event->side);
201: PetscMalloc1(nevents,&event->direction);
202: PetscMalloc1(nevents,&event->terminate);
203: PetscMalloc1(nevents,&event->vtol);
204: for (i=0; i < nevents; i++) {
205: event->direction[i] = direction[i];
206: event->terminate[i] = terminate[i];
207: event->zerocrossing[i] = PETSC_FALSE;
208: event->side[i] = 0;
209: }
210: PetscMalloc1(nevents,&event->events_zero);
211: event->nevents = nevents;
212: event->eventhandler = eventhandler;
213: event->postevent = postevent;
214: event->ctx = ctx;
215: event->timestep_posteventinterval = ts->time_step;
217: event->recsize = 8; /* Initial size of the recorder */
218: PetscOptionsBegin(((PetscObject)ts)->comm,((PetscObject)ts)->prefix,"TS Event options","TS");
219: {
220: PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL);
221: PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);
222: PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);
223: PetscOptionsReal("-ts_event_post_eventinterval_step","Time step after event interval","",event->timestep_posteventinterval,&event->timestep_posteventinterval,NULL);
224: }
225: PetscOptionsEnd();
227: PetscMalloc1(event->recsize,&event->recorder.time);
228: PetscMalloc1(event->recsize,&event->recorder.stepnum);
229: PetscMalloc1(event->recsize,&event->recorder.nevents);
230: PetscMalloc1(event->recsize,&event->recorder.eventidx);
231: for (i=0; i < event->recsize; i++) {
232: PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);
233: }
234: /* Initialize the event recorder */
235: event->recorder.ctr = 0;
237: for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
238: if (flg) {PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);}
240: TSEventDestroy(&ts->event);
241: ts->event = event;
242: ts->event->refct = 1;
243: return(0);
244: }
246: /*
247: TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
248: is reached.
249: */
250: static PetscErrorCode TSEventRecorderResize(TSEvent event)
251: {
253: PetscReal *time;
254: PetscInt *stepnum;
255: PetscInt *nevents;
256: PetscInt **eventidx;
257: PetscInt i,fact=2;
261: /* Create large arrays */
262: PetscMalloc1(fact*event->recsize,&time);
263: PetscMalloc1(fact*event->recsize,&stepnum);
264: PetscMalloc1(fact*event->recsize,&nevents);
265: PetscMalloc1(fact*event->recsize,&eventidx);
266: for (i=0; i < fact*event->recsize; i++) {
267: PetscMalloc1(event->nevents,&eventidx[i]);
268: }
270: /* Copy over data */
271: PetscMemcpy(time,event->recorder.time,event->recsize*sizeof(PetscReal));
272: PetscMemcpy(stepnum,event->recorder.stepnum,event->recsize*sizeof(PetscInt));
273: PetscMemcpy(nevents,event->recorder.nevents,event->recsize*sizeof(PetscInt));
274: for (i=0; i < event->recsize; i++) {
275: PetscMemcpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]*sizeof(PetscInt));
276: }
278: /* Destroy old arrays */
279: for (i=0; i < event->recsize; i++) {
280: PetscFree(event->recorder.eventidx[i]);
281: }
282: PetscFree(event->recorder.eventidx);
283: PetscFree(event->recorder.nevents);
284: PetscFree(event->recorder.stepnum);
285: PetscFree(event->recorder.time);
287: /* Set pointers */
288: event->recorder.time = time;
289: event->recorder.stepnum = stepnum;
290: event->recorder.nevents = nevents;
291: event->recorder.eventidx = eventidx;
293: /* Double size */
294: event->recsize *= fact;
296: return(0);
297: }
299: /*
300: Helper rutine to handle user postenvents and recording
301: */
302: static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
303: {
305: TSEvent event = ts->event;
306: PetscBool terminate = PETSC_FALSE;
307: PetscBool restart = PETSC_FALSE;
308: PetscInt i,ctr,stepnum;
309: PetscBool inflag[2],outflag[2];
310: PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
313: if (event->postevent) {
314: PetscObjectState state_prev,state_post;
315: PetscObjectStateGet((PetscObject)U,&state_prev);
316: (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);
317: PetscObjectStateGet((PetscObject)U,&state_post);
318: if (state_prev != state_post) restart = PETSC_TRUE;
319: }
321: /* Handle termination events and step restart */
322: for (i=0; i<event->nevents_zero; i++) if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
323: inflag[0] = restart; inflag[1] = terminate;
324: MPIU_Allreduce(inflag,outflag,2,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);
325: restart = outflag[0]; terminate = outflag[1];
326: if (restart) {TSRestartStep(ts);}
327: if (terminate) {TSSetConvergedReason(ts,TS_CONVERGED_EVENT);}
328: event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
330: /* Reset event residual functions as states might get changed by the postevent callback */
331: if (event->postevent) {
332: VecLockReadPush(U);
333: (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);
334: VecLockReadPop(U);
335: }
337: /* Cache current time and event residual functions */
338: event->ptime_prev = t;
339: for (i=0; i<event->nevents; i++)
340: event->fvalue_prev[i] = event->fvalue[i];
342: /* Record the event in the event recorder */
343: TSGetStepNumber(ts,&stepnum);
344: ctr = event->recorder.ctr;
345: if (ctr == event->recsize) {
346: TSEventRecorderResize(event);
347: }
348: event->recorder.time[ctr] = t;
349: event->recorder.stepnum[ctr] = stepnum;
350: event->recorder.nevents[ctr] = event->nevents_zero;
351: for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
352: event->recorder.ctr++;
353: return(0);
354: }
356: /* Uses Anderson-Bjorck variant of regula falsi method */
357: PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt)
358: {
359: PetscReal new_dt, scal = 1.0;
360: if (PetscRealPart(fleft)*PetscRealPart(f) < 0) {
361: if (side == 1) {
362: scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
363: if (scal < PETSC_SMALL) scal = 0.5;
364: }
365: new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
366: } else {
367: if (side == -1) {
368: scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
369: if (scal < PETSC_SMALL) scal = 0.5;
370: }
371: new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t;
372: }
373: return PetscMin(dt,new_dt);
374: }
377: PetscErrorCode TSEventHandler(TS ts)
378: {
380: TSEvent event;
381: PetscReal t;
382: Vec U;
383: PetscInt i;
384: PetscReal dt,dt_min;
385: PetscInt rollback=0,in[2],out[2];
386: PetscInt fvalue_sign,fvalueprev_sign;
390: if (!ts->event) return(0);
391: event = ts->event;
393: TSGetTime(ts,&t);
394: TSGetTimeStep(ts,&dt);
395: TSGetSolution(ts,&U);
397: if (event->status == TSEVENT_NONE) {
398: event->timestep_prev = dt;
399: }
401: if (event->status == TSEVENT_RESET_NEXTSTEP) {
402: dt = event->timestep_posteventinterval;
403: TSSetTimeStep(ts,dt);
404: event->status = TSEVENT_NONE;
405: }
407: if (event->status == TSEVENT_NONE) {
408: event->ptime_end = t;
409: }
411: VecLockReadPush(U);
412: (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);
413: VecLockReadPop(U);
415: for (i=0; i < event->nevents; i++) {
416: if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
417: event->status = TSEVENT_ZERO;
418: event->fvalue_right[i] = event->fvalue[i];
419: continue;
420: }
421: fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
422: fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
423: if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
424: switch (event->direction[i]) {
425: case -1:
426: if (fvalue_sign < 0) {
427: rollback = 1;
429: /* Compute new time step */
430: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
432: if (event->monitor) {
433: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
434: }
435: event->fvalue_right[i] = event->fvalue[i];
436: event->side[i] = 1;
438: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
439: event->status = TSEVENT_LOCATED_INTERVAL;
440: }
441: break;
442: case 1:
443: if (fvalue_sign > 0) {
444: rollback = 1;
446: /* Compute new time step */
447: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
449: if (event->monitor) {
450: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
451: }
452: event->fvalue_right[i] = event->fvalue[i];
453: event->side[i] = 1;
455: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
456: event->status = TSEVENT_LOCATED_INTERVAL;
457: }
458: break;
459: case 0:
460: rollback = 1;
462: /* Compute new time step */
463: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
465: if (event->monitor) {
466: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
467: }
468: event->fvalue_right[i] = event->fvalue[i];
469: event->side[i] = 1;
471: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
472: event->status = TSEVENT_LOCATED_INTERVAL;
474: break;
475: }
476: }
477: }
479: in[0] = event->status; in[1] = rollback;
480: MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));
481: event->status = (TSEventStatus)out[0]; rollback = out[1];
482: if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
484: event->nevents_zero = 0;
485: if (event->status == TSEVENT_ZERO) {
486: for (i=0; i < event->nevents; i++) {
487: if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
488: event->events_zero[event->nevents_zero++] = i;
489: if (event->monitor) {
490: PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);
491: }
492: event->zerocrossing[i] = PETSC_FALSE;
493: }
494: event->side[i] = 0;
495: }
496: TSPostEvent(ts,t,U);
498: dt = event->ptime_end - t;
499: if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */
500: dt = event->timestep_prev;
501: event->status = TSEVENT_NONE;
502: }
503: TSSetTimeStep(ts,dt);
504: event->iterctr = 0;
505: return(0);
506: }
508: if (event->status == TSEVENT_LOCATED_INTERVAL) {
509: TSRollBack(ts);
510: TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);
511: event->status = TSEVENT_PROCESSING;
512: event->ptime_right = t;
513: } else {
514: for (i=0; i < event->nevents; i++) {
515: if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) {
516: /* Compute new time step */
517: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
518: event->side[i] = -1;
519: }
520: event->fvalue_prev[i] = event->fvalue[i];
521: }
522: if (event->monitor && event->status == TSEVENT_PROCESSING) {
523: 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);
524: }
525: event->ptime_prev = t;
526: }
528: if (event->status == TSEVENT_PROCESSING) event->iterctr++;
530: MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
531: TSSetTimeStep(ts,dt_min);
532: return(0);
533: }
535: PetscErrorCode TSAdjointEventHandler(TS ts)
536: {
538: TSEvent event;
539: PetscReal t;
540: Vec U;
541: PetscInt ctr;
542: PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
546: if (!ts->event) return(0);
547: event = ts->event;
549: TSGetTime(ts,&t);
550: TSGetSolution(ts,&U);
552: ctr = event->recorder.ctr-1;
553: if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
554: /* Call the user postevent function */
555: if (event->postevent) {
556: (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);
557: event->recorder.ctr--;
558: }
559: }
561: PetscBarrier((PetscObject)ts);
562: return(0);
563: }