Actual source code: tsevent.c

petsc-3.8.0 2017-09-26
Report Typos and Errors
  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,"","TS Event options","");
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      ts_terminate,ts_restart;
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:   MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);
291:   if (ts_terminate) {TSSetConvergedReason(ts,TS_CONVERGED_EVENT);}
292:   event->status = ts_terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP;
293:   MPIU_Allreduce(&restart,&ts_restart,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);
294:   if (ts_restart) ts->steprestart = PETSC_TRUE;

296:   event->ptime_prev = t;
297:   /* Reset event residual functions as states might get changed by the postevent callback */
298:   if (event->postevent) {(*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);}
299:   /* Cache current event residual functions */
300:   for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i];

302:   /* Record the event in the event recorder */
303:   TSGetStepNumber(ts,&stepnum);
304:   ctr = event->recorder.ctr;
305:   if (ctr == event->recsize) {
306:     TSEventRecorderResize(event);
307:   }
308:   event->recorder.time[ctr] = t;
309:   event->recorder.stepnum[ctr] = stepnum;
310:   event->recorder.nevents[ctr] = event->nevents_zero;
311:   for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i];
312:   event->recorder.ctr++;
313:   return(0);
314: }

316: /* Uses Anderson-Bjorck variant of regula falsi method */
317: PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t,PetscReal tright,PetscScalar fleft,PetscScalar f,PetscScalar fright,PetscInt side,PetscReal dt)
318: {
319:   PetscReal new_dt, scal = 1.0;
320:   if (PetscRealPart(fleft)*PetscRealPart(f) < 0) {
321:     if (side == 1) {
322:       scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright);
323:       if (scal < PETSC_SMALL) scal = 0.5;
324:     }
325:     new_dt = (scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft;
326:   } else {
327:     if (side == -1) {
328:       scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft);
329:       if (scal < PETSC_SMALL) scal = 0.5;
330:     }
331:     new_dt = (PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t;
332:   }
333:   return PetscMin(dt,new_dt);
334: }


337: PetscErrorCode TSEventHandler(TS ts)
338: {
340:   TSEvent        event;
341:   PetscReal      t;
342:   Vec            U;
343:   PetscInt       i;
344:   PetscReal      dt,dt_min;
345:   PetscInt       rollback=0,in[2],out[2];
346:   PetscInt       fvalue_sign,fvalueprev_sign;

350:   if (!ts->event) return(0);
351:   event = ts->event;

353:   TSGetTime(ts,&t);
354:   TSGetTimeStep(ts,&dt);
355:   TSGetSolution(ts,&U);

357:   if (event->status == TSEVENT_NONE) {
358:     if (ts->steps == 1) /* After first successful step */
359:       event->timestep_orig = ts->ptime - ts->ptime_prev;
360:     event->timestep_prev = dt;
361:   }

363:   if (event->status == TSEVENT_RESET_NEXTSTEP) {
364:     /* Restore time step */
365:     dt = PetscMin(event->timestep_orig,event->timestep_prev);
366:     TSSetTimeStep(ts,dt);
367:     event->status = TSEVENT_NONE;
368:   }

370:   if (event->status == TSEVENT_NONE) {
371:     event->ptime_end = t;
372:   }

374:   (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);

376:   for (i=0; i < event->nevents; i++) {
377:     if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
378:       event->status = TSEVENT_ZERO;
379:       event->fvalue_right[i] = event->fvalue[i];
380:       continue;
381:     }
382:     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
383:     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
384:     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) {
385:       switch (event->direction[i]) {
386:       case -1:
387:         if (fvalue_sign < 0) {
388:           rollback = 1;

390:           /* Compute new time step */
391:           dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);

393:           if (event->monitor) {
394:             PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
395:           }
396:           event->fvalue_right[i] = event->fvalue[i];
397:           event->side[i] = 1;

399:           if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
400:           event->status = TSEVENT_LOCATED_INTERVAL;
401:         }
402:         break;
403:       case 1:
404:         if (fvalue_sign > 0) {
405:           rollback = 1;

407:           /* Compute new time step */
408:           dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);

410:           if (event->monitor) {
411:             PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
412:           }
413:           event->fvalue_right[i] = event->fvalue[i];
414:           event->side[i] = 1;

416:           if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
417:           event->status = TSEVENT_LOCATED_INTERVAL;
418:         }
419:         break;
420:       case 0:
421:         rollback = 1;

423:         /* Compute new time step */
424:         dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);

426:         if (event->monitor) {
427:           PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);
428:         }
429:         event->fvalue_right[i] = event->fvalue[i];
430:         event->side[i] = 1;

432:         if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE;
433:         event->status = TSEVENT_LOCATED_INTERVAL;

435:         break;
436:       }
437:     }
438:   }

440:   in[0] = event->status; in[1] = rollback;
441:   MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));
442:   event->status = (TSEventStatus)out[0]; rollback = out[1];
443:   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;

445:   event->nevents_zero = 0;
446:   if (event->status == TSEVENT_ZERO) {
447:     for (i=0; i < event->nevents; i++) {
448:       if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) {
449:         event->events_zero[event->nevents_zero++] = i;
450:         if (event->monitor) {
451:           PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);
452:         }
453:         event->zerocrossing[i] = PETSC_FALSE;
454:       }
455:       event->side[i] = 0;
456:     }
457:     TSPostEvent(ts,t,U);

459:     dt = event->ptime_end - t;
460:     if (PetscAbsReal(dt) < PETSC_SMALL) dt += PetscMin(event->timestep_orig,event->timestep_prev); /* XXX Should be done better */
461:     TSSetTimeStep(ts,dt);
462:     event->iterctr = 0;
463:     return(0);
464:   }

466:   if (event->status == TSEVENT_LOCATED_INTERVAL) {
467:     TSRollBack(ts);
468:     TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);
469:     event->status = TSEVENT_PROCESSING;
470:     event->ptime_right = t;
471:   } else {
472:     for (i=0; i < event->nevents; i++) {
473:       if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) {
474:         /* Compute new time step */
475:         dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt);
476:         event->side[i] = -1;
477:       }
478:       event->fvalue_prev[i] = event->fvalue[i];
479:     }
480:     if (event->monitor && event->status == TSEVENT_PROCESSING) {
481:       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);
482:     }
483:     event->ptime_prev = t;
484:   }

486:   if (event->status == TSEVENT_PROCESSING) event->iterctr++;

488:   MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));
489:   TSSetTimeStep(ts,dt_min);
490:   return(0);
491: }

493: PetscErrorCode TSAdjointEventHandler(TS ts)
494: {
496:   TSEvent        event;
497:   PetscReal      t;
498:   Vec            U;
499:   PetscInt       ctr;
500:   PetscBool      forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */

504:   if (!ts->event) return(0);
505:   event = ts->event;

507:   TSGetTime(ts,&t);
508:   TSGetSolution(ts,&U);

510:   ctr = event->recorder.ctr-1;
511:   if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) {
512:     /* Call the user postevent function */
513:     if (event->postevent) {
514:       (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);
515:       event->recorder.ctr--;
516:     }
517:   }

519:   PetscBarrier((PetscObject)ts);
520:   return(0);
521: }