1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 3 /* 4 TSEventInitialize - Initializes TSEvent for TSSolve 5 */ 6 PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U) 7 { 8 PetscErrorCode ierr; 9 10 PetscFunctionBegin; 11 if (!event) PetscFunctionReturn(0); 12 PetscValidPointer(event,1); 13 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 14 PetscValidHeaderSpecific(U,VEC_CLASSID,4); 15 event->ptime_prev = t; 16 event->iterctr = 0; 17 ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr); 18 PetscFunctionReturn(0); 19 } 20 21 PetscErrorCode TSEventDestroy(TSEvent *event) 22 { 23 PetscErrorCode ierr; 24 PetscInt i; 25 26 PetscFunctionBegin; 27 PetscValidPointer(event,1); 28 if (!*event) PetscFunctionReturn(0); 29 if (--(*event)->refct > 0) {*event = NULL; PetscFunctionReturn(0);} 30 31 ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 32 ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 33 ierr = PetscFree((*event)->fvalue_right);CHKERRQ(ierr); 34 ierr = PetscFree((*event)->zerocrossing);CHKERRQ(ierr); 35 ierr = PetscFree((*event)->side);CHKERRQ(ierr); 36 ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 37 ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 38 ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 39 ierr = PetscFree((*event)->vtol);CHKERRQ(ierr); 40 41 for (i=0; i < (*event)->recsize; i++) { 42 ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr); 43 } 44 ierr = PetscFree((*event)->recorder.eventidx);CHKERRQ(ierr); 45 ierr = PetscFree((*event)->recorder.nevents);CHKERRQ(ierr); 46 ierr = PetscFree((*event)->recorder.stepnum);CHKERRQ(ierr); 47 ierr = PetscFree((*event)->recorder.time);CHKERRQ(ierr); 48 49 ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr); 50 ierr = PetscFree(*event);CHKERRQ(ierr); 51 PetscFunctionReturn(0); 52 } 53 54 /*@ 55 TSSetPostEventIntervalStep - Set the time-step used immediately following the event interval 56 57 Logically Collective 58 59 Input Arguments: 60 + ts - time integration context 61 - dt - post event interval step 62 63 Options Database Keys: 64 . -ts_event_post_eventinterval_step <dt> time-step after event interval 65 66 Notes: 67 TSSetPostEventIntervalStep allows one to set a time-step that is used immediately following an event interval. 68 69 This function should be called from the postevent function set with TSSetEventHandler(). 70 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. 74 75 Level: Advanced 76 .seealso: TS, TSEvent, TSSetEventHandler() 77 @*/ 78 PetscErrorCode TSSetPostEventIntervalStep(TS ts,PetscReal dt) 79 { 80 PetscFunctionBegin; 81 ts->event->timestep_posteventinterval = dt; 82 PetscFunctionReturn(0); 83 } 84 85 /*@ 86 TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 87 88 Logically Collective 89 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 94 95 Options Database Keys: 96 . -ts_event_tol <tol> tolerance for event zero crossing 97 98 Notes: 99 Must call TSSetEventHandler() before setting the tolerances. 100 101 The size of vtol is equal to the number of events. 102 103 Level: beginner 104 105 .seealso: TS, TSEvent, TSSetEventHandler() 106 @*/ 107 PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[]) 108 { 109 TSEvent event; 110 PetscInt i; 111 112 PetscFunctionBegin; 113 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 114 if (vtol) PetscValidRealPointer(vtol,3); 115 if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()"); 116 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 PetscFunctionReturn(0); 126 } 127 128 /*@C 129 TSSetEventHandler - Sets a function used for detecting events 130 131 Logically Collective on TS 132 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) 145 146 Calling sequence of eventhandler: 147 PetscErrorCode PetscEventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 148 149 Input Parameters: 150 + ts - the TS context 151 . t - current time 152 . U - current iterate 153 - ctx - [optional] context passed with eventhandler 154 155 Output parameters: 156 . fvalue - function value of events at time t 157 158 Calling sequence of postevent: 159 PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero,PetscInt events_zero[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 160 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 169 170 Level: intermediate 171 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 { 176 PetscErrorCode ierr; 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 185 186 PetscFunctionBegin; 187 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 188 if (nevents) { 189 PetscValidIntPointer(direction,2); 190 PetscValidIntPointer(terminate,3); 191 } 192 193 ierr = PetscNewLog(ts,&event);CHKERRQ(ierr); 194 ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 195 ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 196 ierr = PetscMalloc1(nevents,&event->fvalue_right);CHKERRQ(ierr); 197 ierr = PetscMalloc1(nevents,&event->zerocrossing);CHKERRQ(ierr); 198 ierr = PetscMalloc1(nevents,&event->side);CHKERRQ(ierr); 199 ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 200 ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 201 ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr); 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 ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 209 event->nevents = nevents; 210 event->eventhandler = eventhandler; 211 event->postevent = postevent; 212 event->ctx = ctx; 213 event->timestep_posteventinterval = ts->time_step; 214 215 event->recsize = 8; /* Initial size of the recorder */ 216 ierr = PetscOptionsBegin(((PetscObject)ts)->comm,((PetscObject)ts)->prefix,"TS Event options","TS");CHKERRQ(ierr); 217 { 218 ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","TSSetEventTolerances",tol,&tol,NULL);CHKERRQ(ierr); 219 ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr); 220 ierr = PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);CHKERRQ(ierr); 221 ierr = PetscOptionsReal("-ts_event_post_eventinterval_step","Time step after event interval","",event->timestep_posteventinterval,&event->timestep_posteventinterval,NULL);CHKERRQ(ierr); 222 } 223 ierr = PetscOptionsEnd();CHKERRQ(ierr); 224 225 ierr = PetscMalloc1(event->recsize,&event->recorder.time);CHKERRQ(ierr); 226 ierr = PetscMalloc1(event->recsize,&event->recorder.stepnum);CHKERRQ(ierr); 227 ierr = PetscMalloc1(event->recsize,&event->recorder.nevents);CHKERRQ(ierr); 228 ierr = PetscMalloc1(event->recsize,&event->recorder.eventidx);CHKERRQ(ierr); 229 for (i=0; i < event->recsize; i++) { 230 ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 231 } 232 /* Initialize the event recorder */ 233 event->recorder.ctr = 0; 234 235 for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 236 if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);} 237 238 ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr); 239 ts->event = event; 240 ts->event->refct = 1; 241 PetscFunctionReturn(0); 242 } 243 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 { 250 PetscErrorCode ierr; 251 PetscReal *time; 252 PetscInt *stepnum; 253 PetscInt *nevents; 254 PetscInt **eventidx; 255 PetscInt i,fact=2; 256 257 PetscFunctionBegin; 258 259 /* Create large arrays */ 260 ierr = PetscMalloc1(fact*event->recsize,&time);CHKERRQ(ierr); 261 ierr = PetscMalloc1(fact*event->recsize,&stepnum);CHKERRQ(ierr); 262 ierr = PetscMalloc1(fact*event->recsize,&nevents);CHKERRQ(ierr); 263 ierr = PetscMalloc1(fact*event->recsize,&eventidx);CHKERRQ(ierr); 264 for (i=0; i < fact*event->recsize; i++) { 265 ierr = PetscMalloc1(event->nevents,&eventidx[i]);CHKERRQ(ierr); 266 } 267 268 /* Copy over data */ 269 ierr = PetscArraycpy(time,event->recorder.time,event->recsize);CHKERRQ(ierr); 270 ierr = PetscArraycpy(stepnum,event->recorder.stepnum,event->recsize);CHKERRQ(ierr); 271 ierr = PetscArraycpy(nevents,event->recorder.nevents,event->recsize);CHKERRQ(ierr); 272 for (i=0; i < event->recsize; i++) { 273 ierr = PetscArraycpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]);CHKERRQ(ierr); 274 } 275 276 /* Destroy old arrays */ 277 for (i=0; i < event->recsize; i++) { 278 ierr = PetscFree(event->recorder.eventidx[i]);CHKERRQ(ierr); 279 } 280 ierr = PetscFree(event->recorder.eventidx);CHKERRQ(ierr); 281 ierr = PetscFree(event->recorder.nevents);CHKERRQ(ierr); 282 ierr = PetscFree(event->recorder.stepnum);CHKERRQ(ierr); 283 ierr = PetscFree(event->recorder.time);CHKERRQ(ierr); 284 285 /* Set pointers */ 286 event->recorder.time = time; 287 event->recorder.stepnum = stepnum; 288 event->recorder.nevents = nevents; 289 event->recorder.eventidx = eventidx; 290 291 /* Double size */ 292 event->recsize *= fact; 293 294 PetscFunctionReturn(0); 295 } 296 297 /* 298 Helper routine to handle user postevents and recording 299 */ 300 static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U) 301 { 302 PetscErrorCode ierr; 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 */ 309 310 PetscFunctionBegin; 311 if (event->postevent) { 312 PetscObjectState state_prev,state_post; 313 ierr = PetscObjectStateGet((PetscObject)U,&state_prev);CHKERRQ(ierr); 314 ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 315 ierr = PetscObjectStateGet((PetscObject)U,&state_post);CHKERRQ(ierr); 316 if (state_prev != state_post) restart = PETSC_TRUE; 317 } 318 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 ierr = MPIU_Allreduce(inflag,outflag,2,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 323 restart = outflag[0]; terminate = outflag[1]; 324 if (restart) {ierr = TSRestartStep(ts);CHKERRQ(ierr);} 325 if (terminate) {ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);} 326 event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP; 327 328 /* Reset event residual functions as states might get changed by the postevent callback */ 329 if (event->postevent) { 330 ierr = VecLockReadPush(U);CHKERRQ(ierr); 331 ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); 332 ierr = VecLockReadPop(U);CHKERRQ(ierr); 333 } 334 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]; 339 340 /* Record the event in the event recorder */ 341 ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr); 342 ctr = event->recorder.ctr; 343 if (ctr == event->recsize) { 344 ierr = TSEventRecorderResize(event);CHKERRQ(ierr); 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 PetscFunctionReturn(0); 352 } 353 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 } 373 374 375 PetscErrorCode TSEventHandler(TS ts) 376 { 377 PetscErrorCode ierr; 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; 385 386 PetscFunctionBegin; 387 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 388 if (!ts->event) PetscFunctionReturn(0); 389 event = ts->event; 390 391 ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 392 ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 393 ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 394 395 if (event->status == TSEVENT_NONE) { 396 event->timestep_prev = dt; 397 } 398 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; 403 404 dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt); 405 } 406 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 407 event->status = TSEVENT_NONE; 408 } 409 410 if (event->status == TSEVENT_NONE) { 411 event->ptime_end = t; 412 } 413 414 ierr = VecLockReadPush(U);CHKERRQ(ierr); 415 ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); 416 ierr = VecLockReadPop(U);CHKERRQ(ierr); 417 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; 431 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); 434 435 if (event->monitor) { 436 ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 437 } 438 event->fvalue_right[i] = event->fvalue[i]; 439 event->side[i] = 1; 440 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; 448 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); 451 452 if (event->monitor) { 453 ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 454 } 455 event->fvalue_right[i] = event->fvalue[i]; 456 event->side[i] = 1; 457 458 if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 459 event->status = TSEVENT_LOCATED_INTERVAL; 460 } 461 break; 462 case 0: 463 rollback = 1; 464 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); 467 468 if (event->monitor) { 469 ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); 470 } 471 event->fvalue_right[i] = event->fvalue[i]; 472 event->side[i] = 1; 473 474 if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 475 event->status = TSEVENT_LOCATED_INTERVAL; 476 477 break; 478 } 479 } 480 } 481 482 in[0] = event->status; in[1] = rollback; 483 ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 484 event->status = (TSEventStatus)out[0]; rollback = out[1]; 485 if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 486 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 ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);CHKERRQ(ierr); 494 } 495 event->zerocrossing[i] = PETSC_FALSE; 496 } 497 event->side[i] = 0; 498 } 499 ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr); 500 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; 508 509 dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt,maxdt,10*PETSC_MACHINE_EPSILON,0) ? maxdt : dt); 510 } 511 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 512 event->iterctr = 0; 513 PetscFunctionReturn(0); 514 } 515 516 if (event->status == TSEVENT_LOCATED_INTERVAL) { 517 ierr = TSRollBack(ts);CHKERRQ(ierr); 518 ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr); 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 ierr = 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);CHKERRQ(ierr); 532 } 533 event->ptime_prev = t; 534 } 535 536 if (event->status == TSEVENT_PROCESSING) event->iterctr++; 537 538 ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 539 ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr); 540 PetscFunctionReturn(0); 541 } 542 543 PetscErrorCode TSAdjointEventHandler(TS ts) 544 { 545 PetscErrorCode ierr; 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 */ 551 552 PetscFunctionBegin; 553 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 554 if (!ts->event) PetscFunctionReturn(0); 555 event = ts->event; 556 557 ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 558 ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 559 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 ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 565 event->recorder.ctr--; 566 } 567 } 568 569 PetscBarrier((PetscObject)ts); 570 PetscFunctionReturn(0); 571 } 572 573 /*@ 574 TSGetNumEvents - Get the numbers of events set 575 576 Logically Collective 577 578 Input Argument: 579 + ts - the TS context 580 581 Output Argument: 582 . nevents - number of events 583 584 Level: intermediate 585 586 .seealso: TSSetEventHandler() 587 588 @*/ 589 PetscErrorCode TSGetNumEvents(TS ts,PetscInt * nevents) 590 { 591 PetscFunctionBegin; 592 *nevents = ts->event->nevents; 593 PetscFunctionReturn(0); 594 } 595