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