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 PetscFunctionBegin; 9 if (!event) PetscFunctionReturn(PETSC_SUCCESS); 10 PetscAssertPointer(event, 1); 11 PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 12 PetscValidHeaderSpecific(U, VEC_CLASSID, 4); 13 event->ptime_prev = t; 14 event->iterctr = 0; 15 PetscCall((*event->indicator)(ts, t, U, event->fvalue_prev, event->ctx)); 16 PetscFunctionReturn(PETSC_SUCCESS); 17 } 18 19 PetscErrorCode TSEventDestroy(TSEvent *event) 20 { 21 PetscInt i; 22 23 PetscFunctionBegin; 24 PetscAssertPointer(event, 1); 25 if (!*event) PetscFunctionReturn(PETSC_SUCCESS); 26 if (--(*event)->refct > 0) { 27 *event = NULL; 28 PetscFunctionReturn(PETSC_SUCCESS); 29 } 30 31 PetscCall(PetscFree((*event)->fvalue)); 32 PetscCall(PetscFree((*event)->fvalue_prev)); 33 PetscCall(PetscFree((*event)->fvalue_right)); 34 PetscCall(PetscFree((*event)->zerocrossing)); 35 PetscCall(PetscFree((*event)->side)); 36 PetscCall(PetscFree((*event)->direction)); 37 PetscCall(PetscFree((*event)->terminate)); 38 PetscCall(PetscFree((*event)->events_zero)); 39 PetscCall(PetscFree((*event)->vtol)); 40 41 for (i = 0; i < (*event)->recsize; i++) PetscCall(PetscFree((*event)->recorder.eventidx[i])); 42 PetscCall(PetscFree((*event)->recorder.eventidx)); 43 PetscCall(PetscFree((*event)->recorder.nevents)); 44 PetscCall(PetscFree((*event)->recorder.stepnum)); 45 PetscCall(PetscFree((*event)->recorder.time)); 46 47 PetscCall(PetscViewerDestroy(&(*event)->monitor)); 48 PetscCall(PetscFree(*event)); 49 PetscFunctionReturn(PETSC_SUCCESS); 50 } 51 52 /*@ 53 TSSetPostEventIntervalStep - Set the time-step used immediately following an event interval 54 55 Logically Collective 56 57 Input Parameters: 58 + ts - time integration context 59 - dt - post-event interval step 60 61 Options Database Key: 62 . -ts_event_post_eventinterval_step <dt> - time-step after event interval 63 64 Level: advanced 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 post-event 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 .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()` 76 @*/ 77 PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt) 78 { 79 PetscFunctionBegin; 80 ts->event->timestep_posteventinterval = dt; 81 PetscFunctionReturn(PETSC_SUCCESS); 82 } 83 84 /*@ 85 TSSetEventTolerances - Set tolerances for event zero crossings 86 87 Logically Collective 88 89 Input Parameters: 90 + ts - time integration context 91 . tol - scalar tolerance, `PETSC_DECIDE` to leave current value 92 - vtol - array of tolerances or `NULL`, used in preference to `tol` if present 93 94 Options Database Key: 95 . -ts_event_tol <tol> - tolerance for event zero crossing 96 97 Level: beginner 98 99 Notes: 100 Must call `TSSetEventHandler(`) before setting the tolerances. 101 102 The size of `vtol` is equal to the number of events. 103 104 The tolerance is some measure of how close the event function is to zero for the event detector to stop 105 and declare the time of the event has been detected. 106 107 .seealso: [](ch_ts), `TS`, `TSEvent`, `TSSetEventHandler()` 108 @*/ 109 PetscErrorCode TSSetEventTolerances(TS ts, PetscReal tol, PetscReal vtol[]) 110 { 111 TSEvent event; 112 PetscInt i; 113 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 116 if (vtol) PetscAssertPointer(vtol, 3); 117 PetscCheck(ts->event, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must set the events first by calling TSSetEventHandler()"); 118 119 event = ts->event; 120 if (vtol) { 121 for (i = 0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 122 } else { 123 if ((tol != (PetscReal)PETSC_DECIDE) || (tol != (PetscReal)PETSC_DEFAULT)) { 124 for (i = 0; i < event->nevents; i++) event->vtol[i] = tol; 125 } 126 } 127 PetscFunctionReturn(PETSC_SUCCESS); 128 } 129 130 /*@C 131 TSSetEventHandler - Sets functions used for indicating event locations and handling them 132 133 Logically Collective 134 135 Input Parameters: 136 + ts - the `TS` context obtained from `TSCreate()` 137 . nevents - number of local events 138 . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 139 0.1 => Zero crossing in positive direction, 0 => both ways (one for each event) 140 . terminate - flag to indicate whether time stepping should be terminated after an event is detected (one for each event) 141 . indicator - a change in sign of this function (see `direction`) is used to determine an event has occurred 142 . postevent - [optional] post-event function, this function can change properties of the solution, ODE etc at the time of the event 143 - ctx - [optional] user-defined context for private data for the indicator and post-event routine (use `NULL` if no context is desired) 144 145 Calling sequence of `indicator`: 146 + ts - the `TS` context 147 . t - current time 148 . U - current iterate 149 . fvalue - function value of events at time `t` 150 - ctx - the context passed as the final argument to `TSSetEventHandler()` 151 152 Calling sequence of `postevent`: 153 + ts - the `TS` context 154 . nevents_zero - number of local events whose event function has been marked as crossing 0 155 . events_zero - indices of local events which have been marked as crossing 0 156 . t - current time 157 . U - current solution 158 . forwardsolve - Flag to indicate whether `TS` is doing a forward solve (1) or adjoint solve (0) 159 - ctx - the context passed as the final argument to `TSSetEventHandler()` 160 161 Level: intermediate 162 163 .seealso: [](ch_ts), `TSEvent`, `TSCreate()`, `TSSetTimeStep()`, `TSSetConvergedReason()` 164 @*/ 165 PetscErrorCode TSSetEventHandler(TS ts, PetscInt nevents, PetscInt direction[], PetscBool terminate[], PetscErrorCode (*indicator)(TS ts, PetscReal t, Vec U, PetscScalar fvalue[], void *ctx), PetscErrorCode (*postevent)(TS ts, PetscInt nevents_zero, PetscInt events_zero[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx), void *ctx) 166 { 167 TSAdapt adapt; 168 PetscReal hmin; 169 TSEvent event; 170 PetscInt i; 171 PetscBool flg; 172 #if defined PETSC_USE_REAL_SINGLE 173 PetscReal tol = 1e-4; 174 #else 175 PetscReal tol = 1e-6; 176 #endif 177 178 PetscFunctionBegin; 179 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 180 if (nevents) { 181 PetscAssertPointer(direction, 3); 182 PetscAssertPointer(terminate, 4); 183 } 184 185 PetscCall(PetscNew(&event)); 186 PetscCall(PetscMalloc1(nevents, &event->fvalue)); 187 PetscCall(PetscMalloc1(nevents, &event->fvalue_prev)); 188 PetscCall(PetscMalloc1(nevents, &event->fvalue_right)); 189 PetscCall(PetscMalloc1(nevents, &event->zerocrossing)); 190 PetscCall(PetscMalloc1(nevents, &event->side)); 191 PetscCall(PetscMalloc1(nevents, &event->direction)); 192 PetscCall(PetscMalloc1(nevents, &event->terminate)); 193 PetscCall(PetscMalloc1(nevents, &event->vtol)); 194 for (i = 0; i < nevents; i++) { 195 event->direction[i] = direction[i]; 196 event->terminate[i] = terminate[i]; 197 event->zerocrossing[i] = PETSC_FALSE; 198 event->side[i] = 0; 199 } 200 PetscCall(PetscMalloc1(nevents, &event->events_zero)); 201 event->nevents = nevents; 202 event->indicator = indicator; 203 event->postevent = postevent; 204 event->ctx = ctx; 205 event->timestep_posteventinterval = ts->time_step; 206 PetscCall(TSGetAdapt(ts, &adapt)); 207 PetscCall(TSAdaptGetStepLimits(adapt, &hmin, NULL)); 208 event->timestep_min = hmin; 209 210 event->recsize = 8; /* Initial size of the recorder */ 211 PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "TS Event options", "TS"); 212 { 213 PetscCall(PetscOptionsReal("-ts_event_tol", "Scalar event tolerance for zero crossing check", "TSSetEventTolerances", tol, &tol, NULL)); 214 PetscCall(PetscOptionsName("-ts_event_monitor", "Print choices made by event handler", "", &flg)); 215 PetscCall(PetscOptionsInt("-ts_event_recorder_initial_size", "Initial size of event recorder", "", event->recsize, &event->recsize, NULL)); 216 PetscCall(PetscOptionsReal("-ts_event_post_eventinterval_step", "Time step after event interval", "", event->timestep_posteventinterval, &event->timestep_posteventinterval, NULL)); 217 PetscCall(PetscOptionsReal("-ts_event_post_event_step", "Time step after event", "", event->timestep_postevent, &event->timestep_postevent, NULL)); 218 PetscCall(PetscOptionsReal("-ts_event_dt_min", "Minimum time step considered for TSEvent", "", event->timestep_min, &event->timestep_min, NULL)); 219 } 220 PetscOptionsEnd(); 221 222 PetscCall(PetscMalloc1(event->recsize, &event->recorder.time)); 223 PetscCall(PetscMalloc1(event->recsize, &event->recorder.stepnum)); 224 PetscCall(PetscMalloc1(event->recsize, &event->recorder.nevents)); 225 PetscCall(PetscMalloc1(event->recsize, &event->recorder.eventidx)); 226 for (i = 0; i < event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &event->recorder.eventidx[i])); 227 /* Initialize the event recorder */ 228 event->recorder.ctr = 0; 229 230 for (i = 0; i < event->nevents; i++) event->vtol[i] = tol; 231 if (flg) PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF, "stdout", &event->monitor)); 232 233 PetscCall(TSEventDestroy(&ts->event)); 234 ts->event = event; 235 ts->event->refct = 1; 236 PetscFunctionReturn(PETSC_SUCCESS); 237 } 238 239 /* 240 TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize) 241 is reached. 242 */ 243 static PetscErrorCode TSEventRecorderResize(TSEvent event) 244 { 245 PetscReal *time; 246 PetscInt *stepnum; 247 PetscInt *nevents; 248 PetscInt **eventidx; 249 PetscInt i, fact = 2; 250 251 PetscFunctionBegin; 252 253 /* Create large arrays */ 254 PetscCall(PetscMalloc1(fact * event->recsize, &time)); 255 PetscCall(PetscMalloc1(fact * event->recsize, &stepnum)); 256 PetscCall(PetscMalloc1(fact * event->recsize, &nevents)); 257 PetscCall(PetscMalloc1(fact * event->recsize, &eventidx)); 258 for (i = 0; i < fact * event->recsize; i++) PetscCall(PetscMalloc1(event->nevents, &eventidx[i])); 259 260 /* Copy over data */ 261 PetscCall(PetscArraycpy(time, event->recorder.time, event->recsize)); 262 PetscCall(PetscArraycpy(stepnum, event->recorder.stepnum, event->recsize)); 263 PetscCall(PetscArraycpy(nevents, event->recorder.nevents, event->recsize)); 264 for (i = 0; i < event->recsize; i++) PetscCall(PetscArraycpy(eventidx[i], event->recorder.eventidx[i], event->recorder.nevents[i])); 265 266 /* Destroy old arrays */ 267 for (i = 0; i < event->recsize; i++) PetscCall(PetscFree(event->recorder.eventidx[i])); 268 PetscCall(PetscFree(event->recorder.eventidx)); 269 PetscCall(PetscFree(event->recorder.nevents)); 270 PetscCall(PetscFree(event->recorder.stepnum)); 271 PetscCall(PetscFree(event->recorder.time)); 272 273 /* Set pointers */ 274 event->recorder.time = time; 275 event->recorder.stepnum = stepnum; 276 event->recorder.nevents = nevents; 277 event->recorder.eventidx = eventidx; 278 279 /* Double size */ 280 event->recsize *= fact; 281 282 PetscFunctionReturn(PETSC_SUCCESS); 283 } 284 285 /* 286 Helper routine to handle user post-events and recording 287 */ 288 static PetscErrorCode TSPostEvent(TS ts, PetscReal t, Vec U) 289 { 290 TSEvent event = ts->event; 291 PetscBool terminate = PETSC_FALSE; 292 PetscBool restart = PETSC_FALSE; 293 PetscInt i, ctr, stepnum; 294 PetscBool inflag[2], outflag[2]; 295 PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 296 297 PetscFunctionBegin; 298 if (event->postevent) { 299 PetscObjectState state_prev, state_post; 300 PetscCall(PetscObjectStateGet((PetscObject)U, &state_prev)); 301 PetscCallBack("Event post-event processing", (*event->postevent)(ts, event->nevents_zero, event->events_zero, t, U, forwardsolve, event->ctx)); 302 PetscCall(PetscObjectStateGet((PetscObject)U, &state_post)); 303 if (state_prev != state_post) restart = PETSC_TRUE; 304 } 305 306 /* Handle termination events and step restart */ 307 for (i = 0; i < event->nevents_zero; i++) 308 if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE; 309 inflag[0] = restart; 310 inflag[1] = terminate; 311 PetscCall(MPIU_Allreduce(inflag, outflag, 2, MPIU_BOOL, MPI_LOR, ((PetscObject)ts)->comm)); 312 restart = outflag[0]; 313 terminate = outflag[1]; 314 if (restart) PetscCall(TSRestartStep(ts)); 315 if (terminate) PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_EVENT)); 316 event->status = terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP; 317 318 /* Reset event indicator function values as states might get changed by the post-event callback */ 319 if (event->postevent) { 320 PetscCall(VecLockReadPush(U)); 321 PetscCallBack("Event indicator", (*event->indicator)(ts, t, U, event->fvalue, event->ctx)); 322 PetscCall(VecLockReadPop(U)); 323 } 324 325 /* Cache current time and event residual functions */ 326 event->ptime_prev = t; 327 for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 328 329 /* Record the event in the event recorder */ 330 PetscCall(TSGetStepNumber(ts, &stepnum)); 331 ctr = event->recorder.ctr; 332 if (ctr == event->recsize) PetscCall(TSEventRecorderResize(event)); 333 event->recorder.time[ctr] = t; 334 event->recorder.stepnum[ctr] = stepnum; 335 event->recorder.nevents[ctr] = event->nevents_zero; 336 for (i = 0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 337 event->recorder.ctr++; 338 PetscFunctionReturn(PETSC_SUCCESS); 339 } 340 341 /* Uses Anderson-Bjorck variant of regula falsi method */ 342 static inline PetscReal TSEventComputeStepSize(PetscReal tleft, PetscReal t, PetscReal tright, PetscScalar fleft, PetscScalar f, PetscScalar fright, PetscInt side, PetscReal dt) 343 { 344 PetscReal new_dt, scal = 1.0; 345 if (PetscRealPart(fleft) * PetscRealPart(f) < 0) { 346 if (side == 1) { 347 scal = (PetscRealPart(fright) - PetscRealPart(f)) / PetscRealPart(fright); 348 if (scal < PETSC_SMALL) scal = 0.5; 349 } 350 new_dt = (scal * PetscRealPart(fleft) * t - PetscRealPart(f) * tleft) / (scal * PetscRealPart(fleft) - PetscRealPart(f)) - tleft; 351 } else { 352 if (side == -1) { 353 scal = (PetscRealPart(fleft) - PetscRealPart(f)) / PetscRealPart(fleft); 354 if (scal < PETSC_SMALL) scal = 0.5; 355 } 356 new_dt = (PetscRealPart(f) * tright - scal * PetscRealPart(fright) * t) / (PetscRealPart(f) - scal * PetscRealPart(fright)) - t; 357 } 358 return PetscMin(dt, new_dt); 359 } 360 361 static PetscErrorCode TSEventDetection(TS ts) 362 { 363 TSEvent event = ts->event; 364 PetscReal t; 365 PetscInt i; 366 PetscInt fvalue_sign, fvalueprev_sign; 367 PetscInt in, out; 368 369 PetscFunctionBegin; 370 PetscCall(TSGetTime(ts, &t)); 371 for (i = 0; i < event->nevents; i++) { 372 if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 373 if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 374 event->status = TSEVENT_LOCATED_INTERVAL; 375 if (event->monitor) { 376 PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " interval detected due to zero value (tol=%g) [%g - %g]\n", event->iterctr, i, (double)event->vtol[i], (double)event->ptime_prev, (double)t)); 377 } 378 continue; 379 } 380 if (PetscAbsScalar(event->fvalue_prev[i]) < event->vtol[i]) continue; /* avoid duplicative detection if the previous endpoint is an event location */ 381 fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 382 fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 383 if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign)) { 384 if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 385 event->status = TSEVENT_LOCATED_INTERVAL; 386 if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " interval detected due to sign change [%g - %g]\n", event->iterctr, i, (double)event->ptime_prev, (double)t)); 387 } 388 } 389 in = (PetscInt)event->status; 390 PetscCall(MPIU_Allreduce(&in, &out, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts))); 391 event->status = (TSEventStatus)out; 392 PetscFunctionReturn(PETSC_SUCCESS); 393 } 394 395 static PetscErrorCode TSEventLocation(TS ts, PetscReal *dt) 396 { 397 TSEvent event = ts->event; 398 PetscReal diff = PetscAbsReal((event->ptime_right - event->ptime_prev) / 2); 399 PetscInt i; 400 PetscReal t; 401 PetscInt fvalue_sign, fvalueprev_sign; 402 PetscInt rollback = 0, in[2], out[2]; 403 404 PetscFunctionBegin; 405 PetscCall(TSGetTime(ts, &t)); 406 event->nevents_zero = 0; 407 for (i = 0; i < event->nevents; i++) { 408 if (event->zerocrossing[i]) { 409 if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i] || *dt < event->timestep_min || PetscAbsReal(*dt) < diff * event->vtol[i]) { /* stopping criteria */ 410 event->status = TSEVENT_ZERO; 411 event->fvalue_right[i] = event->fvalue[i]; 412 continue; 413 } 414 /* Compute new time step */ 415 *dt = TSEventComputeStepSize(event->ptime_prev, t, event->ptime_right, event->fvalue_prev[i], event->fvalue[i], event->fvalue_right[i], event->side[i], *dt); 416 fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 417 fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 418 switch (event->direction[i]) { 419 case -1: 420 if (fvalue_sign < 0) { 421 rollback = 1; 422 event->fvalue_right[i] = event->fvalue[i]; 423 event->side[i] = 1; 424 } 425 break; 426 case 1: 427 if (fvalue_sign > 0) { 428 rollback = 1; 429 event->fvalue_right[i] = event->fvalue[i]; 430 event->side[i] = 1; 431 } 432 break; 433 case 0: 434 if (fvalue_sign != fvalueprev_sign) { /* trigger rollback only when there is a sign change */ 435 rollback = 1; 436 event->fvalue_right[i] = event->fvalue[i]; 437 event->side[i] = 1; 438 } 439 break; 440 } 441 if (event->status == TSEVENT_PROCESSING) event->side[i] = -1; 442 } 443 } 444 in[0] = (PetscInt)event->status; 445 in[1] = rollback; 446 PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)ts))); 447 event->status = (TSEventStatus)out[0]; 448 rollback = out[1]; 449 /* If rollback is true, the status will be overwritten so that an event at the endtime of current time step will be postponed to guarantee correct order */ 450 if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 451 if (event->status == TSEVENT_ZERO) { 452 for (i = 0; i < event->nevents; i++) { 453 if (event->zerocrossing[i]) { 454 if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i] || *dt < event->timestep_min || PetscAbsReal(*dt) < diff * event->vtol[i]) { /* stopping criteria */ 455 event->events_zero[event->nevents_zero++] = i; 456 if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Event %" PetscInt_FMT " zero crossing located at time %g\n", event->iterctr, i, (double)t)); 457 event->zerocrossing[i] = PETSC_FALSE; 458 } 459 } 460 event->side[i] = 0; 461 } 462 } 463 PetscFunctionReturn(PETSC_SUCCESS); 464 } 465 466 PetscErrorCode TSEventHandler(TS ts) 467 { 468 TSEvent event; 469 PetscReal t; 470 Vec U; 471 PetscInt i; 472 PetscReal dt, dt_min, dt_reset = 0.0; 473 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 476 if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS); 477 event = ts->event; 478 479 PetscCall(TSGetTime(ts, &t)); 480 PetscCall(TSGetTimeStep(ts, &dt)); 481 PetscCall(TSGetSolution(ts, &U)); 482 483 if (event->status == TSEVENT_NONE) { 484 event->timestep_prev = dt; 485 event->ptime_end = t; 486 } 487 if (event->status == TSEVENT_RESET_NEXTSTEP) { 488 /* user has specified a PostEventInterval dt */ 489 dt = event->timestep_posteventinterval; 490 if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 491 PetscReal maxdt = ts->max_time - t; 492 dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt); 493 } 494 PetscCall(TSSetTimeStep(ts, dt)); 495 event->status = TSEVENT_NONE; 496 } 497 498 PetscCall(VecLockReadPush(U)); 499 PetscCallBack("Event indicator", (*event->indicator)(ts, t, U, event->fvalue, event->ctx)); 500 PetscCall(VecLockReadPop(U)); 501 502 /* Detect the events */ 503 PetscCall(TSEventDetection(ts)); 504 505 /* Locate the events */ 506 if (event->status == TSEVENT_LOCATED_INTERVAL || event->status == TSEVENT_PROCESSING) { 507 /* Approach the zero crosing by setting a new step size */ 508 PetscCall(TSEventLocation(ts, &dt)); 509 /* Roll back when new events are detected */ 510 if (event->status == TSEVENT_LOCATED_INTERVAL) { 511 PetscCall(TSRollBack(ts)); 512 PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); 513 event->iterctr++; 514 } 515 PetscCall(MPIU_Allreduce(&dt, &dt_min, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)ts))); 516 if (dt_reset > 0.0 && dt_reset < dt_min) dt_min = dt_reset; 517 PetscCall(TSSetTimeStep(ts, dt_min)); 518 /* Found the zero crossing */ 519 if (event->status == TSEVENT_ZERO) { 520 PetscCall(TSPostEvent(ts, t, U)); 521 dt = event->ptime_end - t; 522 if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */ 523 dt = event->timestep_prev; 524 event->status = TSEVENT_NONE; 525 } 526 if (event->timestep_postevent) { /* user has specified a PostEvent dt*/ 527 dt = event->timestep_postevent; 528 } 529 if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 530 PetscReal maxdt = ts->max_time - t; 531 dt = dt > maxdt ? maxdt : (PetscIsCloseAtTol(dt, maxdt, 10 * PETSC_MACHINE_EPSILON, 0) ? maxdt : dt); 532 } 533 PetscCall(TSSetTimeStep(ts, dt)); 534 event->iterctr = 0; 535 } 536 /* Have not found the zero crosing yet */ 537 if (event->status == TSEVENT_PROCESSING) { 538 if (event->monitor) PetscCall(PetscViewerASCIIPrintf(event->monitor, "TSEvent: iter %" PetscInt_FMT " - Stepping forward as no event detected in interval [%g - %g]\n", event->iterctr, (double)event->ptime_prev, (double)t)); 539 event->iterctr++; 540 } 541 } 542 if (event->status == TSEVENT_LOCATED_INTERVAL) { /* The step has been rolled back */ 543 event->status = TSEVENT_PROCESSING; 544 event->ptime_right = t; 545 } else { 546 for (i = 0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 547 event->ptime_prev = t; 548 } 549 PetscFunctionReturn(PETSC_SUCCESS); 550 } 551 552 PetscErrorCode TSAdjointEventHandler(TS ts) 553 { 554 TSEvent event; 555 PetscReal t; 556 Vec U; 557 PetscInt ctr; 558 PetscBool forwardsolve = PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 559 560 PetscFunctionBegin; 561 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 562 if (!ts->event) PetscFunctionReturn(PETSC_SUCCESS); 563 event = ts->event; 564 565 PetscCall(TSGetTime(ts, &t)); 566 PetscCall(TSGetSolution(ts, &U)); 567 568 ctr = event->recorder.ctr - 1; 569 if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 570 /* Call the user post-event function */ 571 if (event->postevent) { 572 PetscCallBack("Event post-event processing", (*event->postevent)(ts, event->recorder.nevents[ctr], event->recorder.eventidx[ctr], t, U, forwardsolve, event->ctx)); 573 event->recorder.ctr--; 574 } 575 } 576 PetscCall(PetscBarrier((PetscObject)ts)); 577 PetscFunctionReturn(PETSC_SUCCESS); 578 } 579 580 /*@ 581 TSGetNumEvents - Get the numbers of events currently set to be detected 582 583 Logically Collective 584 585 Input Parameter: 586 . ts - the `TS` context 587 588 Output Parameter: 589 . nevents - the number of events 590 591 Level: intermediate 592 593 .seealso: [](ch_ts), `TSEvent`, `TSSetEventHandler()` 594 @*/ 595 PetscErrorCode TSGetNumEvents(TS ts, PetscInt *nevents) 596 { 597 PetscFunctionBegin; 598 *nevents = ts->event->nevents; 599 PetscFunctionReturn(PETSC_SUCCESS); 600 } 601