xref: /petsc/src/ts/event/tsevent.c (revision 7a728e9f8bd0bc88bbf3140fc3ed678ac237f1b5)
12dc7a7e3SShri Abhyankar 
22dc7a7e3SShri Abhyankar #include <petsc-private/tsimpl.h> /*I  "petscts.h" I*/
32dc7a7e3SShri Abhyankar 
42dc7a7e3SShri Abhyankar #undef __FUNCT__
52dc7a7e3SShri Abhyankar #define __FUNCT__ "TSSetEventMonitor"
62dc7a7e3SShri Abhyankar /*@C
72dc7a7e3SShri Abhyankar    TSSetEventMonitor - Sets a monitoring function used for detecting events
82dc7a7e3SShri Abhyankar 
92dc7a7e3SShri Abhyankar    Logically Collective on TS
102dc7a7e3SShri Abhyankar 
112dc7a7e3SShri Abhyankar    Input Parameters:
122dc7a7e3SShri Abhyankar +  ts - the TS context obtained from TSCreate()
132dc7a7e3SShri Abhyankar .  nevents - number of local events
14d94325d3SShri Abhyankar .  direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction,
15d94325d3SShri Abhyankar                +1 => Zero crossing in positive direction, 0 => both ways (one for each event)
16d94325d3SShri Abhyankar .  terminate - flag to indicate whether time stepping should be terminated after
17d94325d3SShri Abhyankar                event is detected (one for each event)
182dc7a7e3SShri Abhyankar .  eventmonitor - event monitoring routine
192dc7a7e3SShri Abhyankar .  postevent - [optional] post-event function
202dc7a7e3SShri Abhyankar -  mectx - [optional] user-defined context for private data for the
212dc7a7e3SShri Abhyankar               event monitor and post event routine (use NULL if no
222dc7a7e3SShri Abhyankar               context is desired)
232dc7a7e3SShri Abhyankar 
242dc7a7e3SShri Abhyankar    Calling sequence of eventmonitor:
25d94325d3SShri Abhyankar    PetscErrorCode EventMonitor(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void* mectx)
262dc7a7e3SShri Abhyankar 
272dc7a7e3SShri Abhyankar    Input Parameters:
282dc7a7e3SShri Abhyankar +  ts  - the TS context
292dc7a7e3SShri Abhyankar .  t   - current time
302dc7a7e3SShri Abhyankar .  U   - current iterate
312dc7a7e3SShri Abhyankar -  ctx - [optional] context passed with eventmonitor
322dc7a7e3SShri Abhyankar 
332dc7a7e3SShri Abhyankar    Output parameters:
34d94325d3SShri Abhyankar .  fvalue    - function value of events at time t
352dc7a7e3SShri Abhyankar 
362dc7a7e3SShri Abhyankar    Calling sequence of postevent:
372dc7a7e3SShri Abhyankar    PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero, PetscReal t,Vec U,void* ctx)
382dc7a7e3SShri Abhyankar 
392dc7a7e3SShri Abhyankar    Input Parameters:
402dc7a7e3SShri Abhyankar +  ts - the TS context
412dc7a7e3SShri Abhyankar .  nevents_zero - number of local events whose event function is zero
422dc7a7e3SShri Abhyankar .  events_zero  - indices of local events which have reached zero
432dc7a7e3SShri Abhyankar .  t            - current time
442dc7a7e3SShri Abhyankar .  U            - current solution
452dc7a7e3SShri Abhyankar -  ctx          - the context passed with eventmonitor
462dc7a7e3SShri Abhyankar 
472dc7a7e3SShri Abhyankar    Level: intermediate
482dc7a7e3SShri Abhyankar 
492dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor
502dc7a7e3SShri Abhyankar 
512dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason()
522dc7a7e3SShri Abhyankar @*/
53d94325d3SShri Abhyankar PetscErrorCode TSSetEventMonitor(TS ts,PetscInt nevents,PetscInt *direction,PetscBool *terminate,PetscErrorCode (*eventmonitor)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,void*),void *mectx)
542dc7a7e3SShri Abhyankar {
552dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
562dc7a7e3SShri Abhyankar   PetscReal      t;
572dc7a7e3SShri Abhyankar   Vec            U;
582dc7a7e3SShri Abhyankar   TSEvent        event;
59d94325d3SShri Abhyankar   PetscInt       i;
602dc7a7e3SShri Abhyankar 
612dc7a7e3SShri Abhyankar   PetscFunctionBegin;
6242ea6711SShri   ierr = PetscNew(&event);CHKERRQ(ierr);
632dc7a7e3SShri Abhyankar   ierr = PetscMalloc(nevents*sizeof(PetscScalar),&event->fvalue);CHKERRQ(ierr);
642dc7a7e3SShri Abhyankar   ierr = PetscMalloc(nevents*sizeof(PetscScalar),&event->fvalue_prev);CHKERRQ(ierr);
652dc7a7e3SShri Abhyankar   ierr = PetscMalloc(nevents*sizeof(PetscInt),&event->direction);CHKERRQ(ierr);
66d94325d3SShri Abhyankar   ierr = PetscMalloc(nevents*sizeof(PetscInt),&event->terminate);CHKERRQ(ierr);
67d94325d3SShri Abhyankar   for (i=0; i < nevents; i++) {
68d94325d3SShri Abhyankar     event->direction[i] = direction[i];
69d94325d3SShri Abhyankar     event->terminate[i] = terminate[i];
70d94325d3SShri Abhyankar   }
712dc7a7e3SShri Abhyankar   ierr = PetscMalloc(nevents*sizeof(PetscInt),&event->events_zero);CHKERRQ(ierr);
722dc7a7e3SShri Abhyankar   event->monitor = eventmonitor;
732dc7a7e3SShri Abhyankar   event->postevent = postevent;
742dc7a7e3SShri Abhyankar   event->monitorcontext = (void*)mectx;
752dc7a7e3SShri Abhyankar   event->nevents = nevents;
762dc7a7e3SShri Abhyankar 
772dc7a7e3SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
782dc7a7e3SShri Abhyankar   ierr = TSGetTimeStep(ts,&event->initial_timestep);CHKERRQ(ierr);
792dc7a7e3SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
802dc7a7e3SShri Abhyankar   event->ptime_prev = t;
81d94325d3SShri Abhyankar   ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,mectx);CHKERRQ(ierr);
822dc7a7e3SShri Abhyankar   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"TS Event options","");CHKERRQ(ierr);
832dc7a7e3SShri Abhyankar   {
842dc7a7e3SShri Abhyankar     event->tol = 1.0e-6;
852dc7a7e3SShri Abhyankar     ierr = PetscOptionsReal("-ts_event_tol","","",event->tol,&event->tol,NULL);CHKERRQ(ierr);
862dc7a7e3SShri Abhyankar   }
872dc7a7e3SShri Abhyankar   ierr = PetscOptionsEnd();CHKERRQ(ierr);
88d94325d3SShri Abhyankar   ts->event = event;
892dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
902dc7a7e3SShri Abhyankar }
912dc7a7e3SShri Abhyankar 
922dc7a7e3SShri Abhyankar #undef __FUNCT__
932dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent"
942dc7a7e3SShri Abhyankar PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,void *ctx)
952dc7a7e3SShri Abhyankar {
962dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
972dc7a7e3SShri Abhyankar   TSEvent        event=ts->event;
982dc7a7e3SShri Abhyankar   PetscBool      terminate=PETSC_FALSE;
992dc7a7e3SShri Abhyankar   PetscInt       i;
1002dc7a7e3SShri Abhyankar   PetscBool      ts_terminate;
1012dc7a7e3SShri Abhyankar 
1022dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1032dc7a7e3SShri Abhyankar   if (event->postevent) {
1042dc7a7e3SShri Abhyankar     ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,ctx);CHKERRQ(ierr);
1052dc7a7e3SShri Abhyankar   }
1062dc7a7e3SShri Abhyankar   for(i = 0; i < nevents_zero;i++) {
107e105d053SSatish Balay     terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]);
1082dc7a7e3SShri Abhyankar   }
109b4549700SJed Brown   ierr = MPI_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr);
1102dc7a7e3SShri Abhyankar   if (terminate) {
1112dc7a7e3SShri Abhyankar     ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);
1122dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
1132dc7a7e3SShri Abhyankar   } else {
1142dc7a7e3SShri Abhyankar     event->status = TSEVENT_RESET_NEXTSTEP;
1152dc7a7e3SShri Abhyankar   }
1162dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
1172dc7a7e3SShri Abhyankar }
1182dc7a7e3SShri Abhyankar 
1192dc7a7e3SShri Abhyankar #undef __FUNCT__
1202dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitorDestroy"
1212dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitorDestroy(TSEvent *event)
1222dc7a7e3SShri Abhyankar {
1232dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
1242dc7a7e3SShri Abhyankar 
1252dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1262dc7a7e3SShri Abhyankar   ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr);
1272dc7a7e3SShri Abhyankar   ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr);
1282dc7a7e3SShri Abhyankar   ierr = PetscFree((*event)->direction);CHKERRQ(ierr);
129d94325d3SShri Abhyankar   ierr = PetscFree((*event)->terminate);CHKERRQ(ierr);
1302dc7a7e3SShri Abhyankar   ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr);
1312dc7a7e3SShri Abhyankar   ierr = PetscFree(*event);CHKERRQ(ierr);
1322dc7a7e3SShri Abhyankar   *event = NULL;
1332dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
1342dc7a7e3SShri Abhyankar }
1352dc7a7e3SShri Abhyankar 
1362dc7a7e3SShri Abhyankar #undef __FUNCT__
1372dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitor"
1382dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitor(TS ts)
1392dc7a7e3SShri Abhyankar {
1402dc7a7e3SShri Abhyankar   PetscErrorCode ierr;
1412dc7a7e3SShri Abhyankar   TSEvent        event=ts->event;
1422dc7a7e3SShri Abhyankar   PetscReal      t;
1432dc7a7e3SShri Abhyankar   Vec            U;
1442dc7a7e3SShri Abhyankar   PetscInt       i;
1452dc7a7e3SShri Abhyankar   PetscReal      dt;
146b4549700SJed Brown   TSEventStatus  status = event->status;
1472dc7a7e3SShri Abhyankar   PetscInt       rollback=0,in[2],out[2];
1482dc7a7e3SShri Abhyankar 
1492dc7a7e3SShri Abhyankar   PetscFunctionBegin;
1502dc7a7e3SShri Abhyankar 
1512dc7a7e3SShri Abhyankar   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
1522dc7a7e3SShri Abhyankar   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
1532dc7a7e3SShri Abhyankar 
1542dc7a7e3SShri Abhyankar   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
1552dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_RESET_NEXTSTEP) {
1562dc7a7e3SShri Abhyankar     /* Take initial time step */
1572dc7a7e3SShri Abhyankar     dt = event->initial_timestep;
1582dc7a7e3SShri Abhyankar     ts->time_step = dt;
1592dc7a7e3SShri Abhyankar     event->status = TSEVENT_NONE;
1602dc7a7e3SShri Abhyankar   }
1612dc7a7e3SShri Abhyankar 
1622dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_NONE) {
1632dc7a7e3SShri Abhyankar     event->tstepend   = t;
1642dc7a7e3SShri Abhyankar   }
1652dc7a7e3SShri Abhyankar 
1662dc7a7e3SShri Abhyankar   event->nevents_zero = 0;
1672dc7a7e3SShri Abhyankar 
168d94325d3SShri Abhyankar   ierr = (*event->monitor)(ts,t,U,event->fvalue,event->monitorcontext);CHKERRQ(ierr);
1692dc7a7e3SShri Abhyankar   if (event->status != TSEVENT_NONE) {
1702dc7a7e3SShri Abhyankar     for (i=0; i < event->nevents; i++) {
171e105d053SSatish Balay       if (PetscAbsScalar(event->fvalue[i]) < event->tol) {
1722dc7a7e3SShri Abhyankar 	event->status = TSEVENT_ZERO;
1732dc7a7e3SShri Abhyankar 	event->events_zero[event->nevents_zero++] = i;
1742dc7a7e3SShri Abhyankar       }
1752dc7a7e3SShri Abhyankar     }
1762dc7a7e3SShri Abhyankar   }
1772dc7a7e3SShri Abhyankar 
1782dc7a7e3SShri Abhyankar   status = event->status;
179b4549700SJed Brown   ierr = MPI_Allreduce((PetscEnum*)&status,(PetscEnum*)&event->status,1,MPIU_ENUM,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr);
1802dc7a7e3SShri Abhyankar 
1812dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_ZERO) {
1822dc7a7e3SShri Abhyankar     dt = event->tstepend-t;
1832dc7a7e3SShri Abhyankar     ts->time_step = dt;
1842dc7a7e3SShri Abhyankar     ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,event->monitorcontext);CHKERRQ(ierr);
1852dc7a7e3SShri Abhyankar     for (i = 0; i < event->nevents; i++) {
1862dc7a7e3SShri Abhyankar       event->fvalue_prev[i] = event->fvalue[i];
1872dc7a7e3SShri Abhyankar     }
1882dc7a7e3SShri Abhyankar     event->ptime_prev  = t;
1892dc7a7e3SShri Abhyankar     PetscFunctionReturn(0);
1902dc7a7e3SShri Abhyankar   }
1912dc7a7e3SShri Abhyankar 
1922dc7a7e3SShri Abhyankar   for (i = 0; i < event->nevents; i++) {
193d94325d3SShri Abhyankar     PetscInt fvalue_sign,fvalueprev_sign;
194d94325d3SShri Abhyankar     fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i]));
195d94325d3SShri Abhyankar     fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i]));
196*7a728e9fSShri Abhyankar     if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) {
197d94325d3SShri Abhyankar       switch (event->direction[i]) {
198d94325d3SShri Abhyankar       case -1:
199d94325d3SShri Abhyankar 	if (fvalue_sign < 0) {
2002dc7a7e3SShri Abhyankar 	  rollback = 1;
2012dc7a7e3SShri Abhyankar 	  /* Compute linearly interpolated new time step */
202e105d053SSatish Balay 	  dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
2032dc7a7e3SShri Abhyankar 	}
204d94325d3SShri Abhyankar 	break;
205d94325d3SShri Abhyankar       case 1:
206d94325d3SShri Abhyankar 	if (fvalue_sign > 0) {
207d94325d3SShri Abhyankar 	  rollback = 1;
208d94325d3SShri Abhyankar 	  /* Compute linearly interpolated new time step */
209d94325d3SShri Abhyankar 	  dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
2102dc7a7e3SShri Abhyankar 	}
211d94325d3SShri Abhyankar 	break;
212d94325d3SShri Abhyankar       case 0:
213d94325d3SShri Abhyankar 	rollback = 1;
214d94325d3SShri Abhyankar 	/* Compute linearly interpolated new time step */
215d94325d3SShri Abhyankar 	dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i])));
216d94325d3SShri Abhyankar 	break;
217d94325d3SShri Abhyankar       }
218d94325d3SShri Abhyankar     }
219d94325d3SShri Abhyankar   }
220d94325d3SShri Abhyankar   if (rollback) event->status = TSEVENT_LOCATED_INTERVAL;
221d94325d3SShri Abhyankar 
2222dc7a7e3SShri Abhyankar   in[0] = event->status;
2232dc7a7e3SShri Abhyankar   in[1] = rollback;
2242dc7a7e3SShri Abhyankar   ierr = MPI_Allreduce(in,out,2,MPIU_INT,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr);
2252dc7a7e3SShri Abhyankar 
2262dc7a7e3SShri Abhyankar   rollback = out[1];
2272dc7a7e3SShri Abhyankar   if (rollback) {
2282dc7a7e3SShri Abhyankar     event->status = TSEVENT_LOCATED_INTERVAL;
2292dc7a7e3SShri Abhyankar   }
2302dc7a7e3SShri Abhyankar 
2312dc7a7e3SShri Abhyankar   if (event->status == TSEVENT_LOCATED_INTERVAL) {
2322dc7a7e3SShri Abhyankar     ierr = TSRollBack(ts);CHKERRQ(ierr);
2332dc7a7e3SShri Abhyankar     event->status = TSEVENT_PROCESSING;
2342dc7a7e3SShri Abhyankar   } else {
2352dc7a7e3SShri Abhyankar     for (i = 0; i < event->nevents; i++) {
2362dc7a7e3SShri Abhyankar       event->fvalue_prev[i] = event->fvalue[i];
2372dc7a7e3SShri Abhyankar     }
2382dc7a7e3SShri Abhyankar     event->ptime_prev  = t;
2392dc7a7e3SShri Abhyankar     if (event->status == TSEVENT_PROCESSING) {
2402dc7a7e3SShri Abhyankar       dt = event->tstepend - event->ptime_prev;
2412dc7a7e3SShri Abhyankar     }
2422dc7a7e3SShri Abhyankar   }
243e6b5be14SSatish Balay   ierr = MPI_Allreduce(&dt,&(ts->time_step),1,MPIU_REAL,MPI_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr);
2442dc7a7e3SShri Abhyankar   PetscFunctionReturn(0);
2452dc7a7e3SShri Abhyankar }
2462dc7a7e3SShri Abhyankar 
247