12dc7a7e3SShri Abhyankar 2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 32dc7a7e3SShri Abhyankar 42dc7a7e3SShri Abhyankar #undef __FUNCT__ 56fea3669SShri Abhyankar #define __FUNCT__ "TSEventMonitorInitialize" 66fea3669SShri Abhyankar /* 76fea3669SShri Abhyankar TSEventMonitorInitialize - Initializes TSEvent for TSSolve 86fea3669SShri Abhyankar */ 96fea3669SShri Abhyankar PetscErrorCode TSEventMonitorInitialize(TS ts) 106fea3669SShri Abhyankar { 116fea3669SShri Abhyankar PetscErrorCode ierr; 126fea3669SShri Abhyankar PetscReal t; 136fea3669SShri Abhyankar Vec U; 146fea3669SShri Abhyankar TSEvent event=ts->event; 156fea3669SShri Abhyankar 166fea3669SShri Abhyankar PetscFunctionBegin; 176fea3669SShri Abhyankar 186fea3669SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 196fea3669SShri Abhyankar ierr = TSGetTimeStep(ts,&event->initial_timestep);CHKERRQ(ierr); 206fea3669SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 216fea3669SShri Abhyankar event->ptime_prev = t; 226fea3669SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,event->monitorcontext);CHKERRQ(ierr); 236fea3669SShri Abhyankar 246fea3669SShri Abhyankar /* Initialize the event recorder */ 256fea3669SShri Abhyankar event->recorder.ctr = 0; 266fea3669SShri Abhyankar 276fea3669SShri Abhyankar PetscFunctionReturn(0); 286fea3669SShri Abhyankar } 296fea3669SShri Abhyankar 306fea3669SShri Abhyankar #undef __FUNCT__ 312dc7a7e3SShri Abhyankar #define __FUNCT__ "TSSetEventMonitor" 322dc7a7e3SShri Abhyankar /*@C 332dc7a7e3SShri Abhyankar TSSetEventMonitor - Sets a monitoring function used for detecting events 342dc7a7e3SShri Abhyankar 352dc7a7e3SShri Abhyankar Logically Collective on TS 362dc7a7e3SShri Abhyankar 372dc7a7e3SShri Abhyankar Input Parameters: 382dc7a7e3SShri Abhyankar + ts - the TS context obtained from TSCreate() 392dc7a7e3SShri Abhyankar . nevents - number of local events 40d94325d3SShri Abhyankar . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 41d94325d3SShri Abhyankar +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 42d94325d3SShri Abhyankar . terminate - flag to indicate whether time stepping should be terminated after 43d94325d3SShri Abhyankar event is detected (one for each event) 442dc7a7e3SShri Abhyankar . eventmonitor - event monitoring routine 452dc7a7e3SShri Abhyankar . postevent - [optional] post-event function 462dc7a7e3SShri Abhyankar - mectx - [optional] user-defined context for private data for the 472dc7a7e3SShri Abhyankar event monitor and post event routine (use NULL if no 482dc7a7e3SShri Abhyankar context is desired) 492dc7a7e3SShri Abhyankar 502dc7a7e3SShri Abhyankar Calling sequence of eventmonitor: 51d94325d3SShri Abhyankar PetscErrorCode EventMonitor(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void* mectx) 522dc7a7e3SShri Abhyankar 532dc7a7e3SShri Abhyankar Input Parameters: 542dc7a7e3SShri Abhyankar + ts - the TS context 552dc7a7e3SShri Abhyankar . t - current time 562dc7a7e3SShri Abhyankar . U - current iterate 572dc7a7e3SShri Abhyankar - ctx - [optional] context passed with eventmonitor 582dc7a7e3SShri Abhyankar 592dc7a7e3SShri Abhyankar Output parameters: 60d94325d3SShri Abhyankar . fvalue - function value of events at time t 612dc7a7e3SShri Abhyankar 622dc7a7e3SShri Abhyankar Calling sequence of postevent: 63031fbad4SShri Abhyankar PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero, PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 642dc7a7e3SShri Abhyankar 652dc7a7e3SShri Abhyankar Input Parameters: 662dc7a7e3SShri Abhyankar + ts - the TS context 672dc7a7e3SShri Abhyankar . nevents_zero - number of local events whose event function is zero 682dc7a7e3SShri Abhyankar . events_zero - indices of local events which have reached zero 692dc7a7e3SShri Abhyankar . t - current time 702dc7a7e3SShri Abhyankar . U - current solution 71031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 722dc7a7e3SShri Abhyankar - ctx - the context passed with eventmonitor 732dc7a7e3SShri Abhyankar 742dc7a7e3SShri Abhyankar Level: intermediate 752dc7a7e3SShri Abhyankar 762dc7a7e3SShri Abhyankar .keywords: TS, event, set, monitor 772dc7a7e3SShri Abhyankar 782dc7a7e3SShri Abhyankar .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 792dc7a7e3SShri Abhyankar @*/ 80031fbad4SShri 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,PetscBool,void*),void *mectx) 812dc7a7e3SShri Abhyankar { 822dc7a7e3SShri Abhyankar PetscErrorCode ierr; 832dc7a7e3SShri Abhyankar TSEvent event; 84d94325d3SShri Abhyankar PetscInt i; 85*006e6a18SShri Abhyankar PetscBool flg; 862dc7a7e3SShri Abhyankar 872dc7a7e3SShri Abhyankar PetscFunctionBegin; 8842ea6711SShri ierr = PetscNew(&event);CHKERRQ(ierr); 89854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 90854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 91854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 92854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 93d94325d3SShri Abhyankar for (i=0; i < nevents; i++) { 94d94325d3SShri Abhyankar event->direction[i] = direction[i]; 95d94325d3SShri Abhyankar event->terminate[i] = terminate[i]; 96d94325d3SShri Abhyankar } 97854ce69bSBarry Smith ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 982dc7a7e3SShri Abhyankar event->monitor = eventmonitor; 992dc7a7e3SShri Abhyankar event->postevent = postevent; 1002dc7a7e3SShri Abhyankar event->monitorcontext = (void*)mectx; 1012dc7a7e3SShri Abhyankar event->nevents = nevents; 1022dc7a7e3SShri Abhyankar 103f7aea88cSShri Abhyankar for(i=0; i < MAXEVENTRECORDERS; i++) { 104f7aea88cSShri Abhyankar ierr = PetscMalloc1(nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 105f7aea88cSShri Abhyankar } 106f7aea88cSShri Abhyankar 1072dc7a7e3SShri Abhyankar ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"TS Event options","");CHKERRQ(ierr); 1082dc7a7e3SShri Abhyankar { 1092dc7a7e3SShri Abhyankar event->tol = 1.0e-6; 1102dc7a7e3SShri Abhyankar ierr = PetscOptionsReal("-ts_event_tol","","",event->tol,&event->tol,NULL);CHKERRQ(ierr); 111*006e6a18SShri Abhyankar ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr); 1122dc7a7e3SShri Abhyankar } 1132dc7a7e3SShri Abhyankar ierr = PetscOptionsEnd();CHKERRQ(ierr); 114*006e6a18SShri Abhyankar 115*006e6a18SShri Abhyankar if(flg) { 116*006e6a18SShri Abhyankar ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)ts),"stdout",&event->mon);CHKERRQ(ierr); 117*006e6a18SShri Abhyankar } 118d94325d3SShri Abhyankar ts->event = event; 1192dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1202dc7a7e3SShri Abhyankar } 1212dc7a7e3SShri Abhyankar 1222dc7a7e3SShri Abhyankar #undef __FUNCT__ 1232dc7a7e3SShri Abhyankar #define __FUNCT__ "TSPostEvent" 124031fbad4SShri Abhyankar /* 125031fbad4SShri Abhyankar TSPostEvent - Does post event processing by calling the user-defined postevent function 126031fbad4SShri Abhyankar 127031fbad4SShri Abhyankar Logically Collective on TS 128031fbad4SShri Abhyankar 129031fbad4SShri Abhyankar Input Parameters: 130031fbad4SShri Abhyankar + ts - the TS context 131031fbad4SShri Abhyankar . nevents_zero - number of local events whose event function is zero 132031fbad4SShri Abhyankar . events_zero - indices of local events which have reached zero 133031fbad4SShri Abhyankar . t - current time 134031fbad4SShri Abhyankar . U - current solution 135031fbad4SShri Abhyankar . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 136031fbad4SShri Abhyankar - ctx - the context passed with eventmonitor 137031fbad4SShri Abhyankar 138031fbad4SShri Abhyankar Level: intermediate 139031fbad4SShri Abhyankar 140031fbad4SShri Abhyankar .keywords: TS, event, set, monitor 141031fbad4SShri Abhyankar 142031fbad4SShri Abhyankar .seealso: TSSetEventMonitor(),TSEvent 143031fbad4SShri Abhyankar */ 144031fbad4SShri Abhyankar #undef __FUNCT__ 145031fbad4SShri Abhyankar #define __FUNCT__ "TSPostEvent" 146031fbad4SShri Abhyankar PetscErrorCode TSPostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void *ctx) 1472dc7a7e3SShri Abhyankar { 1482dc7a7e3SShri Abhyankar PetscErrorCode ierr; 1492dc7a7e3SShri Abhyankar TSEvent event=ts->event; 1502dc7a7e3SShri Abhyankar PetscBool terminate=PETSC_FALSE; 151d0578d90SShri Abhyankar PetscInt i,ctr,stepnum; 1522dc7a7e3SShri Abhyankar PetscBool ts_terminate; 1532dc7a7e3SShri Abhyankar 1542dc7a7e3SShri Abhyankar PetscFunctionBegin; 1552dc7a7e3SShri Abhyankar if (event->postevent) { 156031fbad4SShri Abhyankar ierr = (*event->postevent)(ts,nevents_zero,events_zero,t,U,forwardsolve,ctx);CHKERRQ(ierr); 1572dc7a7e3SShri Abhyankar } 1582dc7a7e3SShri Abhyankar for(i = 0; i < nevents_zero;i++) { 159e105d053SSatish Balay terminate = (PetscBool)(terminate || event->terminate[events_zero[i]]); 1602dc7a7e3SShri Abhyankar } 161b4549700SJed Brown ierr = MPI_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 1622dc7a7e3SShri Abhyankar if (terminate) { 1632dc7a7e3SShri Abhyankar ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 1642dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 1652dc7a7e3SShri Abhyankar } else { 1662dc7a7e3SShri Abhyankar event->status = TSEVENT_RESET_NEXTSTEP; 1672dc7a7e3SShri Abhyankar } 168f7aea88cSShri Abhyankar 169d0578d90SShri Abhyankar /* Record the event in the event recorder */ 170d0578d90SShri Abhyankar ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr); 171f7aea88cSShri Abhyankar ctr = event->recorder.ctr; 172f7aea88cSShri Abhyankar if (ctr == MAXEVENTRECORDERS) { 173f7aea88cSShri Abhyankar SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Exceeded limit (=%d) for number of events recorded",MAXEVENTRECORDERS); 174f7aea88cSShri Abhyankar } 175f7aea88cSShri Abhyankar event->recorder.time[ctr] = t; 176d0578d90SShri Abhyankar event->recorder.stepnum[ctr] = stepnum; 177f7aea88cSShri Abhyankar event->recorder.nevents[ctr] = nevents_zero; 178f7aea88cSShri Abhyankar for(i=0; i < nevents_zero; i++) event->recorder.eventidx[ctr][i] = events_zero[i]; 179f7aea88cSShri Abhyankar event->recorder.ctr++; 180f7aea88cSShri Abhyankar 18173967516SShri Abhyankar /* Reset the event residual functions as states might get changed by the postevent callback */ 18273967516SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue_prev,event->monitorcontext);CHKERRQ(ierr); 18373967516SShri Abhyankar event->ptime_prev = t; 18473967516SShri Abhyankar 1852dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 1862dc7a7e3SShri Abhyankar } 1872dc7a7e3SShri Abhyankar 1882dc7a7e3SShri Abhyankar #undef __FUNCT__ 1892dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitorDestroy" 1902dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitorDestroy(TSEvent *event) 1912dc7a7e3SShri Abhyankar { 1922dc7a7e3SShri Abhyankar PetscErrorCode ierr; 193f7aea88cSShri Abhyankar PetscInt i; 1942dc7a7e3SShri Abhyankar 1952dc7a7e3SShri Abhyankar PetscFunctionBegin; 1962dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 1972dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 1982dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 199d94325d3SShri Abhyankar ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 2002dc7a7e3SShri Abhyankar ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 201f7aea88cSShri Abhyankar for(i=0; i < MAXEVENTRECORDERS; i++) { 202302440fdSBarry Smith ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr); 203f7aea88cSShri Abhyankar } 204*006e6a18SShri Abhyankar ierr = PetscViewerDestroy(&(*event)->mon);CHKERRQ(ierr); 2052dc7a7e3SShri Abhyankar ierr = PetscFree(*event);CHKERRQ(ierr); 2062dc7a7e3SShri Abhyankar *event = NULL; 2072dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2082dc7a7e3SShri Abhyankar } 2092dc7a7e3SShri Abhyankar 2102dc7a7e3SShri Abhyankar #undef __FUNCT__ 2112dc7a7e3SShri Abhyankar #define __FUNCT__ "TSEventMonitor" 2122dc7a7e3SShri Abhyankar PetscErrorCode TSEventMonitor(TS ts) 2132dc7a7e3SShri Abhyankar { 2142dc7a7e3SShri Abhyankar PetscErrorCode ierr; 2152dc7a7e3SShri Abhyankar TSEvent event=ts->event; 2162dc7a7e3SShri Abhyankar PetscReal t; 2172dc7a7e3SShri Abhyankar Vec U; 2182dc7a7e3SShri Abhyankar PetscInt i; 2192dc7a7e3SShri Abhyankar PetscReal dt; 220b4549700SJed Brown TSEventStatus status = event->status; 2212dc7a7e3SShri Abhyankar PetscInt rollback=0,in[2],out[2]; 222031fbad4SShri Abhyankar PetscBool forwardsolve=PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 2232dc7a7e3SShri Abhyankar 2242dc7a7e3SShri Abhyankar PetscFunctionBegin; 2252dc7a7e3SShri Abhyankar 2262dc7a7e3SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 2272dc7a7e3SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 2282dc7a7e3SShri Abhyankar 2292dc7a7e3SShri Abhyankar ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 2302dc7a7e3SShri Abhyankar if (event->status == TSEVENT_RESET_NEXTSTEP) { 2312dc7a7e3SShri Abhyankar /* Take initial time step */ 2322dc7a7e3SShri Abhyankar dt = event->initial_timestep; 2332dc7a7e3SShri Abhyankar ts->time_step = dt; 2342dc7a7e3SShri Abhyankar event->status = TSEVENT_NONE; 2352dc7a7e3SShri Abhyankar } 2362dc7a7e3SShri Abhyankar 2372dc7a7e3SShri Abhyankar if (event->status == TSEVENT_NONE) { 2382dc7a7e3SShri Abhyankar event->tstepend = t; 2392dc7a7e3SShri Abhyankar } 2402dc7a7e3SShri Abhyankar 2412dc7a7e3SShri Abhyankar event->nevents_zero = 0; 2422dc7a7e3SShri Abhyankar 243d94325d3SShri Abhyankar ierr = (*event->monitor)(ts,t,U,event->fvalue,event->monitorcontext);CHKERRQ(ierr); 2442dc7a7e3SShri Abhyankar for (i=0; i < event->nevents; i++) { 245e105d053SSatish Balay if (PetscAbsScalar(event->fvalue[i]) < event->tol) { 2462dc7a7e3SShri Abhyankar event->status = TSEVENT_ZERO; 2472dc7a7e3SShri Abhyankar event->events_zero[event->nevents_zero++] = i; 248*006e6a18SShri Abhyankar if(event->mon) { 249*006e6a18SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->mon,"TSEvent : Event %D zero crossing at time %g\n",i,(double)t);CHKERRQ(ierr); 250*006e6a18SShri Abhyankar } 2512dc7a7e3SShri Abhyankar } 2522dc7a7e3SShri Abhyankar } 2532dc7a7e3SShri Abhyankar 2542dc7a7e3SShri Abhyankar status = event->status; 255b4549700SJed Brown ierr = MPI_Allreduce((PetscEnum*)&status,(PetscEnum*)&event->status,1,MPIU_ENUM,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 2562dc7a7e3SShri Abhyankar 2572dc7a7e3SShri Abhyankar if (event->status == TSEVENT_ZERO) { 258031fbad4SShri Abhyankar ierr = TSPostEvent(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->monitorcontext);CHKERRQ(ierr); 259c540466cSShri Abhyankar dt = event->tstepend-t; 260bcbf8bb3SShri Abhyankar if(PetscAbsReal(dt) < PETSC_SMALL) dt += event->initial_timestep; 26193fbeba1SShri Abhyankar ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 2622dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 2632dc7a7e3SShri Abhyankar } 2642dc7a7e3SShri Abhyankar 2652dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 266d94325d3SShri Abhyankar PetscInt fvalue_sign,fvalueprev_sign; 267d94325d3SShri Abhyankar fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 268d94325d3SShri Abhyankar fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 26993fbeba1SShri Abhyankar if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->tol)) { 270d94325d3SShri Abhyankar switch (event->direction[i]) { 271d94325d3SShri Abhyankar case -1: 272d94325d3SShri Abhyankar if (fvalue_sign < 0) { 2732dc7a7e3SShri Abhyankar rollback = 1; 2742dc7a7e3SShri Abhyankar /* Compute linearly interpolated new time step */ 275e105d053SSatish Balay dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 276*006e6a18SShri Abhyankar if(event->mon) { 277*006e6a18SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->mon,"TSEvent : Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 278*006e6a18SShri Abhyankar } 2792dc7a7e3SShri Abhyankar } 280d94325d3SShri Abhyankar break; 281d94325d3SShri Abhyankar case 1: 282d94325d3SShri Abhyankar if (fvalue_sign > 0) { 283d94325d3SShri Abhyankar rollback = 1; 284d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 285d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 286*006e6a18SShri Abhyankar if(event->mon) { 287*006e6a18SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->mon,"TSEvent : Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 288*006e6a18SShri Abhyankar } 2892dc7a7e3SShri Abhyankar } 290d94325d3SShri Abhyankar break; 291d94325d3SShri Abhyankar case 0: 292d94325d3SShri Abhyankar rollback = 1; 293d94325d3SShri Abhyankar /* Compute linearly interpolated new time step */ 294d94325d3SShri Abhyankar dt = PetscMin(dt,PetscRealPart(-event->fvalue_prev[i]*(t - event->ptime_prev)/(event->fvalue[i] - event->fvalue_prev[i]))); 295*006e6a18SShri Abhyankar if(event->mon) { 296*006e6a18SShri Abhyankar ierr = PetscViewerASCIIPrintf(event->mon,"TSEvent : Event %D interval located [%g - %g]\n",i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 297*006e6a18SShri Abhyankar } 298d94325d3SShri Abhyankar break; 299d94325d3SShri Abhyankar } 300d94325d3SShri Abhyankar } 301d94325d3SShri Abhyankar } 302d94325d3SShri Abhyankar if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 303d94325d3SShri Abhyankar 3042dc7a7e3SShri Abhyankar in[0] = event->status; 3052dc7a7e3SShri Abhyankar in[1] = rollback; 3062dc7a7e3SShri Abhyankar ierr = MPI_Allreduce(in,out,2,MPIU_INT,MPI_MAX,((PetscObject)ts)->comm);CHKERRQ(ierr); 3072dc7a7e3SShri Abhyankar 3082dc7a7e3SShri Abhyankar rollback = out[1]; 3092dc7a7e3SShri Abhyankar if (rollback) { 3102dc7a7e3SShri Abhyankar event->status = TSEVENT_LOCATED_INTERVAL; 3112dc7a7e3SShri Abhyankar } 3122dc7a7e3SShri Abhyankar 3132dc7a7e3SShri Abhyankar if (event->status == TSEVENT_LOCATED_INTERVAL) { 3142dc7a7e3SShri Abhyankar ierr = TSRollBack(ts);CHKERRQ(ierr); 315c540466cSShri Abhyankar ts->steps--; 316c540466cSShri Abhyankar ts->total_steps--; 317db84a1feSShri Abhyankar ts->reason = TS_CONVERGED_ITERATING; 3182dc7a7e3SShri Abhyankar event->status = TSEVENT_PROCESSING; 3192dc7a7e3SShri Abhyankar } else { 3202dc7a7e3SShri Abhyankar for (i = 0; i < event->nevents; i++) { 3212dc7a7e3SShri Abhyankar event->fvalue_prev[i] = event->fvalue[i]; 3222dc7a7e3SShri Abhyankar } 3232dc7a7e3SShri Abhyankar event->ptime_prev = t; 3242dc7a7e3SShri Abhyankar if (event->status == TSEVENT_PROCESSING) { 3252dc7a7e3SShri Abhyankar dt = event->tstepend - event->ptime_prev; 3262dc7a7e3SShri Abhyankar } 3272dc7a7e3SShri Abhyankar } 328a9b180a6SBarry Smith ierr = MPI_Allreduce(&dt,&(ts->time_step),1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 3292dc7a7e3SShri Abhyankar PetscFunctionReturn(0); 3302dc7a7e3SShri Abhyankar } 3312dc7a7e3SShri Abhyankar 332d0578d90SShri Abhyankar #undef __FUNCT__ 333d0578d90SShri Abhyankar #define __FUNCT__ "TSAdjointEventMonitor" 334d0578d90SShri Abhyankar PetscErrorCode TSAdjointEventMonitor(TS ts) 335d0578d90SShri Abhyankar { 336d0578d90SShri Abhyankar PetscErrorCode ierr; 337d0578d90SShri Abhyankar TSEvent event=ts->event; 338d0578d90SShri Abhyankar PetscReal t; 339d0578d90SShri Abhyankar Vec U; 340d0578d90SShri Abhyankar PetscInt ctr; 341d0578d90SShri Abhyankar PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 342d0578d90SShri Abhyankar 343d0578d90SShri Abhyankar PetscFunctionBegin; 344d0578d90SShri Abhyankar 345d0578d90SShri Abhyankar ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 346d0578d90SShri Abhyankar ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 347d0578d90SShri Abhyankar 348d0578d90SShri Abhyankar ctr = event->recorder.ctr-1; 349bcbf8bb3SShri Abhyankar if(ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 350d0578d90SShri Abhyankar /* Call the user postevent function */ 351d0578d90SShri Abhyankar if(event->postevent) { 352d0578d90SShri Abhyankar ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->monitorcontext);CHKERRQ(ierr); 353d0578d90SShri Abhyankar event->recorder.ctr--; 354d0578d90SShri Abhyankar } 355d0578d90SShri Abhyankar } 356d0578d90SShri Abhyankar 357d0578d90SShri Abhyankar PetscBarrier((PetscObject)ts); 358d0578d90SShri Abhyankar PetscFunctionReturn(0); 359d0578d90SShri Abhyankar } 360d0578d90SShri Abhyankar 361d0578d90SShri Abhyankar 362d0578d90SShri Abhyankar 363