Actual source code: tsevent.c
petsc-3.10.5 2019-03-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: TSSetEventTolerances - Set tolerances for event zero crossings when using event handler
57: Logically Collective
59: Input Arguments:
60: + ts - time integration context
61: . tol - scalar tolerance, PETSC_DECIDE to leave current value
62: - vtol - array of tolerances or NULL, used in preference to tol if present
64: Options Database Keys:
65: . -ts_event_tol <tol> tolerance for event zero crossing
67: Notes:
68: Must call TSSetEventHandler() before setting the tolerances.
70: The size of vtol is equal to the number of events.
72: Level: beginner
74: .seealso: TS, TSEvent, TSSetEventHandler()
75: @*/
76: PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[])
77: {
78: TSEvent event;
79: PetscInt i;
84: if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()");
86: event = ts->event;
87: if (vtol) {
88: for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i];
89: } else {
90: if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) {
91: for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
92: }
93: }
94: return(0);
95: }
97: /*@C
98: TSSetEventHandler - Sets a monitoring function used for detecting events
100: Logically Collective on TS
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)
115: Calling sequence of eventhandler:
116: PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx)
118: Input Parameters:
119: + ts - the TS context
120: . t - current time
121: . U - current iterate
122: - ctx - [optional] context passed with eventhandler
124: Output parameters:
125: . fvalue - function value of events at time t
127: Calling sequence of postevent:
128: PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
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
139: Level: intermediate
141: .keywords: TS, event, set
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: {
148: TSEvent event;
149: PetscInt i;
150: PetscBool flg;
151: #if defined PETSC_USE_REAL_SINGLE
152: PetscReal tol=1e-4;
153: #else
154: PetscReal tol=1e-6;
155: #endif
159: if (nevents) {
162: }
164: PetscNewLog(ts,&event);
165: PetscMalloc1(nevents,&event->fvalue);
166: PetscMalloc1(nevents,&event->fvalue_prev);
167: PetscMalloc1(nevents,&event->fvalue_right);
168: PetscMalloc1(nevents,&event->zerocrossing);
169: PetscMalloc1(nevents,&event->side);
170: PetscMalloc1(nevents,&event->direction);
171: PetscMalloc1(nevents,&event->terminate);
172: PetscMalloc1(nevents,&event->vtol);
173: for (i=0; i < nevents; i++) {
174: event->direction[i] = direction[i];
175: event->terminate[i] = terminate[i];
176: event->zerocrossing[i] = PETSC_FALSE;
177: event->side[i] = 0;
178: }
179: PetscMalloc1(nevents,&event->events_zero);
180: event->nevents = nevents;
181: event->eventhandler = eventhandler;
182: event->postevent = postevent;
183: event->ctx = ctx;
185: event->recsize = 8; /* Initial size of the recorder */
186: PetscOptionsBegin(((PetscObject)ts)->comm,((PetscObject)ts)->prefix,"TS Event options","TS");
187: {
188: PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL);
189: PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);
190: PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);
191: }
192: PetscOptionsEnd();
194: PetscMalloc1(event->recsize,&event->recorder.time);
195: PetscMalloc1(event->recsize,&event->recorder.stepnum);
196: PetscMalloc1(event->recsize,&event->recorder.nevents);
197: PetscMalloc1(event->recsize,&event->recorder.eventidx);
198: for (i=0; i < event->recsize; i++) {
199: PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);
200: }
201: /* Initialize the event recorder */
202: event->recorder.ctr = 0;
204: for (i=0; i < event->nevents; i++) event->vtol[i] = tol;
205: if (flg) {PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);}
207: TSEventDestroy(&ts->event);
208: ts->event = event;
209: ts->event->refct = 1;
210: return(0);
211: }
213: /*
214: TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize)
215: is reached.
216: */
217: static PetscErrorCode TSEventRecorderResize(TSEvent event)
218: {
220: PetscReal *time;
221: PetscInt *stepnum;
222: PetscInt *nevents;
223: PetscInt **eventidx;
224: PetscInt i,fact=2;
228: /* Create large arrays */
229: PetscMalloc1(fact*event->recsize,&time);
230: PetscMalloc1(fact*event->recsize,&stepnum);
231: PetscMalloc1(fact*event->recsize,&nevents);
232: PetscMalloc1(fact*event->recsize,&eventidx);
233: for (i=0; i < fact*event->recsize; i++) {
234: PetscMalloc1(event->nevents,&eventidx[i]);
235: }
237: /* Copy over data */
238: PetscMemcpy(time,event->recorder.time,event->recsize*sizeof(PetscReal));
239: PetscMemcpy(stepnum,event->recorder.stepnum,event->recsize*sizeof(PetscInt));
240: PetscMemcpy(nevents,event->recorder.nevents,event->recsize*sizeof(PetscInt));
241: for (i=0; i < event->recsize; i++) {
242: PetscMemcpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]*sizeof(PetscInt));
243: }
245: /* Destroy old arrays */
246: for (i=0; i < event->recsize; i++) {
247: PetscFree(event->recorder.eventidx[i]);
248: }
249: PetscFree(event->recorder.eventidx);
250: PetscFree(event->recorder.nevents);
251: PetscFree(event->recorder.stepnum);
252: PetscFree(event->recorder.time);
254: /* Set pointers */
255: event->recorder.time = time;
256: event->recorder.stepnum = stepnum;
257: event->recorder.nevents = nevents;
258: event->recorder.eventidx = eventidx;
260: /* Double size */
261: event->recsize *= fact;
263: return(0);
264: }
266: /*
267: Helper rutine to handle user postenvents and recording
268: */
269: static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U)
270: {
272: TSEvent event = ts->event;
273: PetscBool terminate = PETSC_FALSE;
274: PetscBool restart = PETSC_FALSE;
275: PetscInt i,ctr,stepnum;
276: PetscBool inflag[2],outflag[2];
277: PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */
280: if (event->postevent) {
281: PetscObjectState state_prev,state_post;
282: PetscObjectStateGet((PetscObject)U,&state_prev);
283: (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);
284: PetscObjectStateGet((PetscObject)U,&state_post);
285: if (state_prev != state_post) restart = PETSC_TRUE;
286: }
288: /* Handle termination events and step restart */
289: for (i=0; i<event->nevents_zero; i++) if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE;
290: inflag[0] = restart; inflag[1] = terminate;
291: MPIU_Allreduce(inflag,outflag,2,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);
292: restart = outflag[0]; terminate = outflag[1];
293: if (restart) {TSRestartStep(ts);}
294: if (terminate) {TSSetConvergedReason(ts,TS_CONVERGED_EVENT);}
295: event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
297: /* Reset event residual functions as states might get changed by the postevent callback */
298: if (event->postevent) {
299: VecLockPush(U);
300: (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);
301: VecLockPop(U);
302: }
304: /* Cache current time and event residual functions */
305: event->ptime_prev = t;
306: for (i=0; i<event->nevents; i++)
307: event->fvalue_prev[i] = event->fvalue[i];
309: /* Record the event in the event recorder */
310: TSGetStepNumber(ts,&stepnum);
311: ctr = event->recorder.ctr;
312: if (ctr == event->recsize) {
313: TSEventRecorderResize(event);
314: }
315: event->recorder.time[ctr] = t;
316: event->recorder.stepnum[ctr] = stepnum;
317: event->recorder.nevents[ctr] = event->nevents_zero;
318: for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
319: event->recorder.ctr++;
320: return(0);
321: }
323: /* Uses Anderson-Bjorck variant of regula falsi method */
324: PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt)
325: {
326: PetscReal new_dt, scal = 1.0;
327: if (PetscRealPart(fleft)*PetscRealPart(f) < 0) {
328: if (side == 1) {
329: scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
330: if (scal < PETSC_SMALL) scal = 0.5;
331: }
332: new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
333: } else {
334: if (side == -1) {
335: scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
336: if (scal < PETSC_SMALL) scal = 0.5;
337: }
338: new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t;
339: }
340: return PetscMin(dt,new_dt);
341: }
344: PetscErrorCode TSEventHandler(TS ts)
345: {
347: TSEvent event;
348: PetscReal t;
349: Vec U;
350: PetscInt i;
351: PetscReal dt,dt_min;
352: PetscInt rollback=0,in[2],out[2];
353: PetscInt fvalue_sign,fvalueprev_sign;
357: if (!ts->event) return(0);
358: event = ts->event;
360: TSGetTime(ts,&t);
361: TSGetTimeStep(ts,&dt);
362: TSGetSolution(ts,&U);
364: if (event->status == TSEVENT_NONE) {
365: if (ts->steps == 1) /* After first successful step */
366: event->timestep_orig = ts->ptime - ts->ptime_prev;
367: event->timestep_prev = dt;
368: }
370: if (event->status == TSEVENT_RESET_NEXTSTEP) {
371: /* Restore time step */
372: dt = PetscMin(event->timestep_orig,event->timestep_prev);
373: TSSetTimeStep(ts,dt);
374: event->status = TSEVENT_NONE;
375: }
377: if (event->status == TSEVENT_NONE) {
378: event->ptime_end = t;
379: }
381: VecLockPush(U);
382: (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);
383: VecLockPop(U);
385: for (i=0; i < event->nevents; i++) {
386: if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
387: event->status = TSEVENT_ZERO;
388: event->fvalue_right[i] = event->fvalue[i];
389: continue;
390: }
391: fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
392: fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
393: if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
394: switch (event->direction[i]) {
395: case -1:
396: if (fvalue_sign < 0) {
397: rollback = 1;
399: /* Compute new time step */
400: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
402: if (event->monitor) {
403: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
404: }
405: event->fvalue_right[i] = event->fvalue[i];
406: event->side[i] = 1;
408: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
409: event->status = TSEVENT_LOCATED_INTERVAL;
410: }
411: break;
412: case 1:
413: if (fvalue_sign > 0) {
414: rollback = 1;
416: /* Compute new time step */
417: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
419: if (event->monitor) {
420: PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
421: }
422: event->fvalue_right[i] = event->fvalue[i];
423: event->side[i] = 1;
425: if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
426: event->status = TSEVENT_LOCATED_INTERVAL;
427: }
428: break;
429: case 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;
444: break;
445: }
446: }
447: }
449: in[0] = event->status; in[1] = rollback;
450: MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));
451: event->status = (TSEventStatus)out[0]; rollback = out[1];
452: if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
454: event->nevents_zero = 0;
455: if (event->status == TSEVENT_ZERO) {
456: for (i=0; i < event->nevents; i++) {
457: if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
458: event->events_zero[event->nevents_zero++] = i;
459: if (event->monitor) {
460: PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);
461: }
462: event->zerocrossing[i] = PETSC_FALSE;
463: }
464: event->side[i] = 0;
465: }
466: TSPostEvent(ts,t,U);
468: dt = event->ptime_end - t;
469: if (PetscAbsReal(dt) < PETSC_SMALL) dt += PetscMin(event->timestep_orig,event->timestep_prev); /* XXX Should be done better */
470: TSSetTimeStep(ts,dt);
471: event->iterctr = 0;
472: return(0);
473: }
475: if (event->status == TSEVENT_LOCATED_INTERVAL) {
476: TSRollBack(ts);
477: TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);
478: event->status = TSEVENT_PROCESSING;
479: event->ptime_right = t;
480: } else {
481: for (i=0; i < event->nevents; i++) {
482: if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) {
483: /* Compute new time step */
484: dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
485: event->side[i] = -1;
486: }
487: event->fvalue_prev[i] = event->fvalue[i];
488: }
489: if (event->monitor && event->status == TSEVENT_PROCESSING) {
490: 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);
491: }
492: event->ptime_prev = t;
493: }
495: if (event->status == TSEVENT_PROCESSING) event->iterctr++;
497: MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
498: TSSetTimeStep(ts,dt_min);
499: return(0);
500: }
502: PetscErrorCode TSAdjointEventHandler(TS ts)
503: {
505: TSEvent event;
506: PetscReal t;
507: Vec U;
508: PetscInt ctr;
509: PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */
513: if (!ts->event) return(0);
514: event = ts->event;
516: TSGetTime(ts,&t);
517: TSGetSolution(ts,&U);
519: ctr = event->recorder.ctr-1;
520: if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
521: /* Call the user postevent function */
522: if (event->postevent) {
523: (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);
524: event->recorder.ctr--;
525: }
526: }
528: PetscBarrier((PetscObject)ts);
529: return(0);
530: }