Actual source code: tsevent.c
petsc-3.14.6 2021-03-30
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 = NULL; 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 used 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 one to set a time-step that is used immediately following an event interval.
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.
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 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: .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
173: @*/
174: 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)
175: {
177: TSEvent event;
178: PetscInt i;
179: PetscBool flg;
180: #if defined PETSC_USE_REAL_SINGLE
181: PetscReal tol=1e-4;
182: #else
183: PetscReal tol=1e-6;
184: #endif
188: if (nevents) {
191: }
193: PetscNewLog(ts,&event);
194: PetscMalloc1(nevents,&event->fvalue);
195: PetscMalloc1(nevents,&event->fvalue_prev);
196: PetscMalloc1(nevents,&event->fvalue_right);
197: PetscMalloc1(nevents,&event->zerocrossing);
198: PetscMalloc1(nevents,&event->side);
199: PetscMalloc1(nevents,&event->direction);
200: PetscMalloc1(nevents,&event->terminate);
201: PetscMalloc1(nevents,&event->vtol);
202: for (i=0; i < nevents; i++) {
203: event->direction[i] = direction[i];
204: event->terminate[i] = terminate[i];
205: event->zerocrossing[i] = PETSC_FALSE;
206: event->side[i] = 0;
207: }
208: PetscMalloc1(nevents,&event->events_zero);
209: event->nevents = nevents;
210: event->eventhandler = eventhandler;
211: event->postevent = postevent;
212: event->ctx = ctx;
213: event->timestep_posteventinterval = ts->time_step;
215: event->recsize = 8; /* Initial size of the recorder */
216: PetscOptionsBegin(((PetscObject)ts)->comm,((PetscObject)ts)->prefix,"TS Event options","TS");
217: {
218: PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL);
219: PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);
220: PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);
221: PetscOptionsReal("-ts_event_post_eventinterval_step","Time step after event interval","",event->timestep_posteventinterval,&event->timestep_posteventinterval,NULL);
222: }
223: PetscOptionsEnd();
225: PetscMalloc1(event->recsize,&event->recorder.time);
226: PetscMalloc1(event->recsize,&event->recorder.stepnum);
227: PetscMalloc1(event->recsize,&event->recorder.nevents);
228: PetscMalloc1(event->recsize,&event->recorder.eventidx);
229: for (i=0; i < event->recsize; i++) {
230: PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);
231: }
232: /* Initialize the event recorder */
233: event->recorder.ctr = 0;
235: for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
236: if (flg) {PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);}
238: TSEventDestroy(&ts->event);
239: ts->event = event;
240: ts->event->refct = 1;
241: return(0);
242: }
244: /*
245: TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
246: is reached.
247: */
248: static PetscErrorCode TSEventRecorderResize(TSEvent event)
249: {
251: PetscReal *time;
252: PetscInt *stepnum;
253: PetscInt *nevents;
254: PetscInt **eventidx;
255: PetscInt i,fact=2;
259: /* Create large arrays */
260: PetscMalloc1(fact*event->recsize,&time);
261: PetscMalloc1(fact*event->recsize,&stepnum);
262: PetscMalloc1(fact*event->recsize,&nevents);
263: PetscMalloc1(fact*event->recsize,&eventidx);
264: for (i=0; i < fact*event->recsize; i++) {
265: PetscMalloc1(event->nevents,&eventidx[i]);
266: }
268: /* Copy over data */
269: PetscArraycpy(time,event->recorder.time,event->recsize);
270: PetscArraycpy(stepnum,event->recorder.stepnum,event->recsize);
271: PetscArraycpy(nevents,event->recorder.nevents,event->recsize);
272: for (i=0; i < event->recsize; i++) {
273: PetscArraycpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]);
274: }
276: /* Destroy old arrays */
277: for (i=0; i < event->recsize; i++) {
278: PetscFree(event->recorder.eventidx[i]);
279: }
280: PetscFree(event->recorder.eventidx);
281: PetscFree(event->recorder.nevents);
282: PetscFree(event->recorder.stepnum);
283: PetscFree(event->recorder.time);
285: /* Set pointers */
286: event->recorder.time = time;
287: event->recorder.stepnum = stepnum;
288: event->recorder.nevents = nevents;
289: event->recorder.eventidx = eventidx;
291: /* Double size */
292: event->recsize *= fact;
294: return(0);
295: }
297: /*
298: Helper routine to handle user postevents and recording
299: */
300: static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
301: {
303: TSEvent event = ts->event;
304: PetscBool terminate = PETSC_FALSE;
305: PetscBool restart = PETSC_FALSE;
306: PetscInt i,ctr,stepnum;
307: PetscBool inflag[2],outflag[2];
308: PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
311: if (event->postevent) {
312: PetscObjectState state_prev,state_post;
313: PetscObjectStateGet((PetscObject)U,&state_prev);
314: (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);
315: PetscObjectStateGet((PetscObject)U,&state_post);
316: if (state_prev != state_post) restart = PETSC_TRUE;
317: }
319: /* Handle termination events and step restart */
320: for (i=0; i<event->nevents_zero; i++) if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
321: inflag[0] = restart; inflag[1] = terminate;
322: MPIU_Allreduce(inflag,outflag,2,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);
323: restart = outflag[0]; terminate = outflag[1];
324: if (restart) {TSRestartStep(ts);}
325: if (terminate) {TSSetConvergedReason(ts,TS_CONVERGED_EVENT);}
326: event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
328: /* Reset event residual functions as states might get changed by the postevent callback */
329: if (event->postevent) {
330: VecLockReadPush(U);
331: (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);
332: VecLockReadPop(U);
333: }
335: /* Cache current time and event residual functions */
336: event->ptime_prev = t;
337: for (i=0; i<event->nevents; i++)
338: event->fvalue_prev[i] = event->fvalue[i];
340: /* Record the event in the event recorder */
341: TSGetStepNumber(ts,&stepnum);
342: ctr = event->recorder.ctr;
343: if (ctr == event->recsize) {
344: TSEventRecorderResize(event);
345: }
346: event->recorder.time[ctr] = t;
347: event->recorder.stepnum[ctr] = stepnum;
348: event->recorder.nevents[ctr] = event->nevents_zero;
349: for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
350: event->recorder.ctr++;
351: return(0);
352: }
354: /* Uses Anderson-Bjorck variant of regula falsi method */
355: PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt)
356: {
357: PetscReal new_dt, scal = 1.0;
358: if (PetscRealPart(fleft)*PetscRealPart(f) < 0) {
359: if (side == 1) {
360: scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
361: if (scal < PETSC_SMALL) scal = 0.5;
362: }
363: new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
364: } else {
365: if (side == -1) {
366: scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
367: if (scal < PETSC_SMALL) scal = 0.5;
368: }
369: new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t;
370: }
371: return PetscMin(dt,new_dt);
372: }
375: PetscErrorCode TSEventHandler(TS ts)
376: {
378: TSEvent event;
379: PetscReal t;
380: Vec U;
381: PetscInt i;
382: PetscReal dt,dt_min;
383: PetscInt rollback=0,in[2],out[2];
384: PetscInt fvalue_sign,fvalueprev_sign;
388: if (!ts->event) return(0);
389: event = ts->event;
391: TSGetTime(ts,&t);
392: TSGetTimeStep(ts,&dt);
393: TSGetSolution(ts,&U);
395: if (event->status == TSEVENT_NONE) {
396: event->timestep_prev = dt;
397: }
399: if (event->status == TSEVENT_RESET_NEXTSTEP) {
400: dt = event->timestep_posteventinterval;
401: if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
402: PetscReal maxdt = ts->max_time-t;
404: dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt);
405: }
406: TSSetTimeStep(ts,dt);
407: event->status = TSEVENT_NONE;
408: }
410: if (event->status == TSEVENT_NONE) {
411: event->ptime_end = t;
412: }
414: VecLockReadPush(U);
415: (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);
416: VecLockReadPop(U);
418: for (i=0; i < event->nevents; i++) {
419: if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
420: event->status = TSEVENT_ZERO;
421: event->fvalue_right[i] = event->fvalue[i];
422: continue;
423: }
424: fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
425: fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
426: if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
427: switch (event->direction[i]) {
428: case -1:
429: if (fvalue_sign < 0) {
430: rollback = 1;
432: /* Compute new time step */
433: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
435: if (event->monitor) {
436: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
437: }
438: event->fvalue_right[i] = event->fvalue[i];
439: event->side[i] = 1;
441: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
442: event->status = TSEVENT_LOCATED_INTERVAL;
443: }
444: break;
445: case 1:
446: if (fvalue_sign > 0) {
447: rollback = 1;
449: /* Compute new time step */
450: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
452: if (event->monitor) {
453: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
454: }
455: event->fvalue_right[i] = event->fvalue[i];
456: event->side[i] = 1;
458: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
459: event->status = TSEVENT_LOCATED_INTERVAL;
460: }
461: break;
462: case 0:
463: rollback = 1;
465: /* Compute new time step */
466: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
468: if (event->monitor) {
469: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
470: }
471: event->fvalue_right[i] = event->fvalue[i];
472: event->side[i] = 1;
474: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
475: event->status = TSEVENT_LOCATED_INTERVAL;
477: break;
478: }
479: }
480: }
482: in[0] = event->status; in[1] = rollback;
483: MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));
484: event->status = (TSEventStatus)out[0]; rollback = out[1];
485: if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
487: event->nevents_zero = 0;
488: if (event->status == TSEVENT_ZERO) {
489: for (i=0; i < event->nevents; i++) {
490: if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
491: event->events_zero[event->nevents_zero++] = i;
492: if (event->monitor) {
493: PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);
494: }
495: event->zerocrossing[i] = PETSC_FALSE;
496: }
497: event->side[i] = 0;
498: }
499: TSPostEvent(ts,t,U);
501: dt = event->ptime_end - t;
502: if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */
503: dt = event->timestep_prev;
504: event->status = TSEVENT_NONE;
505: }
506: if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
507: PetscReal maxdt = ts->max_time-t;
509: dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt);
510: }
511: TSSetTimeStep(ts,dt);
512: event->iterctr = 0;
513: return(0);
514: }
516: if (event->status == TSEVENT_LOCATED_INTERVAL) {
517: TSRollBack(ts);
518: TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);
519: event->status = TSEVENT_PROCESSING;
520: event->ptime_right = t;
521: } else {
522: for (i=0; i < event->nevents; i++) {
523: if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) {
524: /* Compute new time step */
525: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
526: event->side[i] = -1;
527: }
528: event->fvalue_prev[i] = event->fvalue[i];
529: }
530: if (event->monitor && event->status == TSEVENT_PROCESSING) {
531: 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);
532: }
533: event->ptime_prev = t;
534: }
536: if (event->status == TSEVENT_PROCESSING) event->iterctr++;
538: MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
539: TSSetTimeStep(ts,dt_min);
540: return(0);
541: }
543: PetscErrorCode TSAdjointEventHandler(TS ts)
544: {
546: TSEvent event;
547: PetscReal t;
548: Vec U;
549: PetscInt ctr;
550: PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
554: if (!ts->event) return(0);
555: event = ts->event;
557: TSGetTime(ts,&t);
558: TSGetSolution(ts,&U);
560: ctr = event->recorder.ctr-1;
561: if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
562: /* Call the user postevent function */
563: if (event->postevent) {
564: (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);
565: event->recorder.ctr--;
566: }
567: }
569: PetscBarrier((PetscObject)ts);
570: return(0);
571: }