1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "TSEventInitialize" 5 /* 6 TSEventInitialize - Initializes TSEvent for TSSolve 7 */ 8 PetscErrorCode TSEventInitialize(TSEvent event,TS ts,PetscReal t,Vec U) 9 { 10 PetscErrorCode ierr; 11 12 PetscFunctionBegin; 13 if (!event) PetscFunctionReturn(0); 14 PetscValidPointer(event,1); 15 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 16 PetscValidHeaderSpecific(U,VEC_CLASSID,4); 17 event->ptime_prev = t; 18 event->iterctr=0; 19 ierr = (*event->eventhandler)(ts,t,U,event->fvalue_prev,event->ctx);CHKERRQ(ierr); 20 /* Initialize the event recorder */ 21 event->recorder.ctr = 0; 22 PetscFunctionReturn(0); 23 } 24 25 #undef __FUNCT__ 26 #define __FUNCT__ "TSEventDestroy" 27 PetscErrorCode TSEventDestroy(TSEvent *event) 28 { 29 PetscErrorCode ierr; 30 PetscInt i; 31 32 PetscFunctionBegin; 33 PetscValidPointer(event,1); 34 if (!*event) PetscFunctionReturn(0); 35 ierr = PetscFree((*event)->fvalue);CHKERRQ(ierr); 36 ierr = PetscFree((*event)->fvalue_prev);CHKERRQ(ierr); 37 ierr = PetscFree((*event)->fvalue_right);CHKERRQ(ierr); 38 ierr = PetscFree((*event)->zerocrossing);CHKERRQ(ierr); 39 ierr = PetscFree((*event)->side);CHKERRQ(ierr); 40 ierr = PetscFree((*event)->direction);CHKERRQ(ierr); 41 ierr = PetscFree((*event)->terminate);CHKERRQ(ierr); 42 ierr = PetscFree((*event)->events_zero);CHKERRQ(ierr); 43 ierr = PetscFree((*event)->vtol);CHKERRQ(ierr); 44 45 46 for (i=0; i < (*event)->recsize; i++) { 47 ierr = PetscFree((*event)->recorder.eventidx[i]);CHKERRQ(ierr); 48 } 49 ierr = PetscFree((*event)->recorder.eventidx);CHKERRQ(ierr); 50 ierr = PetscFree((*event)->recorder.nevents);CHKERRQ(ierr); 51 ierr = PetscFree((*event)->recorder.stepnum);CHKERRQ(ierr); 52 ierr = PetscFree((*event)->recorder.time);CHKERRQ(ierr); 53 54 ierr = PetscViewerDestroy(&(*event)->monitor);CHKERRQ(ierr); 55 ierr = PetscFree(*event);CHKERRQ(ierr); 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "TSSetEventTolerances" 61 /*@ 62 TSSetEventTolerances - Set tolerances for event zero crossings when using event handler 63 64 Logically Collective 65 66 Input Arguments: 67 + ts - time integration context 68 . tol - scalar tolerance, PETSC_DECIDE to leave current value 69 - vtol - array of tolerances or NULL, used in preference to tol if present 70 71 Options Database Keys: 72 . -ts_event_tol <tol> tolerance for event zero crossing 73 74 Notes: 75 Must call TSSetEventHandler() before setting the tolerances. 76 77 The size of vtol is equal to the number of events. 78 79 Level: beginner 80 81 .seealso: TS, TSEvent, TSSetEventHandler() 82 @*/ 83 PetscErrorCode TSSetEventTolerances(TS ts,PetscReal tol,PetscReal vtol[]) 84 { 85 TSEvent event; 86 PetscInt i; 87 88 PetscFunctionBegin; 89 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 90 if (vtol) PetscValidRealPointer(vtol,3); 91 if (!ts->event) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set the events first by calling TSSetEventHandler()"); 92 93 event = ts->event; 94 if (vtol) { 95 for (i=0; i < event->nevents; i++) event->vtol[i] = vtol[i]; 96 } else { 97 if (tol != PETSC_DECIDE || tol != PETSC_DEFAULT) { 98 for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 99 } 100 } 101 PetscFunctionReturn(0); 102 } 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "TSSetEventHandler" 106 /*@C 107 TSSetEventHandler - Sets a monitoring function used for detecting events 108 109 Logically Collective on TS 110 111 Input Parameters: 112 + ts - the TS context obtained from TSCreate() 113 . nevents - number of local events 114 . direction - direction of zero crossing to be detected. -1 => Zero crossing in negative direction, 115 +1 => Zero crossing in positive direction, 0 => both ways (one for each event) 116 . terminate - flag to indicate whether time stepping should be terminated after 117 event is detected (one for each event) 118 . eventhandler - event monitoring routine 119 . postevent - [optional] post-event function 120 - ctx - [optional] user-defined context for private data for the 121 event monitor and post event routine (use NULL if no 122 context is desired) 123 124 Calling sequence of eventhandler: 125 PetscErrorCode EventHandler(TS ts,PetscReal t,Vec U,PetscScalar fvalue[],void* ctx) 126 127 Input Parameters: 128 + ts - the TS context 129 . t - current time 130 . U - current iterate 131 - ctx - [optional] context passed with eventhandler 132 133 Output parameters: 134 . fvalue - function value of events at time t 135 136 Calling sequence of postevent: 137 PetscErrorCode PostEvent(TS ts,PetscInt nevents_zero, PetscInt events_zero[], PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 138 139 Input Parameters: 140 + ts - the TS context 141 . nevents_zero - number of local events whose event function is zero 142 . events_zero - indices of local events which have reached zero 143 . t - current time 144 . U - current solution 145 . forwardsolve - Flag to indicate whether TS is doing a forward solve (1) or adjoint solve (0) 146 - ctx - the context passed with eventhandler 147 148 Level: intermediate 149 150 .keywords: TS, event, set, monitor 151 152 .seealso: TSCreate(), TSSetTimeStep(), TSSetConvergedReason() 153 @*/ 154 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) 155 { 156 PetscErrorCode ierr; 157 TSEvent event; 158 PetscInt i; 159 PetscBool flg; 160 #if defined PETSC_USE_REAL_SINGLE 161 PetscReal tol=1e-4; 162 #else 163 PetscReal to=1e-6; 164 #endif 165 166 PetscFunctionBegin; 167 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 168 PetscValidIntPointer(direction,2); 169 PetscValidIntPointer(terminate,3); 170 171 ierr = PetscNewLog(ts,&event);CHKERRQ(ierr); 172 ierr = PetscMalloc1(nevents,&event->fvalue);CHKERRQ(ierr); 173 ierr = PetscMalloc1(nevents,&event->fvalue_prev);CHKERRQ(ierr); 174 ierr = PetscMalloc1(nevents,&event->fvalue_right);CHKERRQ(ierr); 175 ierr = PetscMalloc1(nevents,&event->zerocrossing);CHKERRQ(ierr); 176 ierr = PetscMalloc1(nevents,&event->side);CHKERRQ(ierr); 177 ierr = PetscMalloc1(nevents,&event->direction);CHKERRQ(ierr); 178 ierr = PetscMalloc1(nevents,&event->terminate);CHKERRQ(ierr); 179 ierr = PetscMalloc1(nevents,&event->vtol);CHKERRQ(ierr); 180 for (i=0; i < nevents; i++) { 181 event->direction[i] = direction[i]; 182 event->terminate[i] = terminate[i]; 183 event->zerocrossing[i] = PETSC_FALSE; 184 event->side[i] = 0; 185 } 186 ierr = PetscMalloc1(nevents,&event->events_zero);CHKERRQ(ierr); 187 event->nevents = nevents; 188 event->eventhandler = eventhandler; 189 event->postevent = postevent; 190 event->ctx = ctx; 191 192 event->recsize = 12; /* Initial size of the recorder */ 193 ierr = PetscOptionsBegin(((PetscObject)ts)->comm,"","TS Event options","");CHKERRQ(ierr); 194 { 195 ierr = PetscOptionsReal("-ts_event_tol","Scalar event tolerance for zero crossing check","",tol,&tol,NULL);CHKERRQ(ierr); 196 ierr = PetscOptionsName("-ts_event_monitor","Print choices made by event handler","",&flg);CHKERRQ(ierr); 197 ierr = PetscOptionsInt("-ts_event_recorder_initial_size","Initial size of event recorder","",event->recsize,&event->recsize,NULL);CHKERRQ(ierr); 198 } 199 PetscOptionsEnd(); 200 201 ierr = PetscMalloc1(event->recsize,&event->recorder.time);CHKERRQ(ierr); 202 ierr = PetscMalloc1(event->recsize,&event->recorder.stepnum);CHKERRQ(ierr); 203 ierr = PetscMalloc1(event->recsize,&event->recorder.nevents);CHKERRQ(ierr); 204 ierr = PetscMalloc1(event->recsize,&event->recorder.eventidx);CHKERRQ(ierr); 205 for (i=0; i < event->recsize; i++) { 206 ierr = PetscMalloc1(event->nevents,&event->recorder.eventidx[i]);CHKERRQ(ierr); 207 } 208 209 for (i=0; i < event->nevents; i++) event->vtol[i] = tol; 210 if (flg) {ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,"stdout",&event->monitor);CHKERRQ(ierr);} 211 212 ierr = TSEventDestroy(&ts->event);CHKERRQ(ierr); 213 ts->event = event; 214 PetscFunctionReturn(0); 215 } 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "TSEventRecorderResize" 219 /* 220 TSEventRecorderResize - Resizes (2X) the event recorder arrays whenever the recording limit (event->recsize) 221 is reached. 222 */ 223 PetscErrorCode TSEventRecorderResize(TSEvent event) 224 { 225 PetscErrorCode ierr; 226 PetscReal *time; 227 PetscInt *stepnum; 228 PetscInt *nevents; 229 PetscInt **eventidx; 230 PetscInt i,fact=2; 231 232 PetscFunctionBegin; 233 234 /* Create large arrays */ 235 ierr = PetscMalloc1(fact*event->recsize,&time);CHKERRQ(ierr); 236 ierr = PetscMalloc1(fact*event->recsize,&stepnum);CHKERRQ(ierr); 237 ierr = PetscMalloc1(fact*event->recsize,&nevents);CHKERRQ(ierr); 238 ierr = PetscMalloc1(fact*event->recsize,&eventidx);CHKERRQ(ierr); 239 for (i=0; i < fact*event->recsize; i++) { 240 ierr = PetscMalloc1(event->nevents,&eventidx[i]);CHKERRQ(ierr); 241 } 242 243 /* Copy over data */ 244 ierr = PetscMemcpy(time,event->recorder.time,event->recsize*sizeof(PetscReal));CHKERRQ(ierr); 245 ierr = PetscMemcpy(stepnum,event->recorder.stepnum,event->recsize*sizeof(PetscInt));CHKERRQ(ierr); 246 ierr = PetscMemcpy(nevents,event->recorder.nevents,event->recsize*sizeof(PetscInt));CHKERRQ(ierr); 247 for(i=0; i < event->recsize; i++) { 248 ierr = PetscMemcpy(eventidx[i],event->recorder.eventidx[i],event->recorder.nevents[i]*sizeof(PetscInt));CHKERRQ(ierr); 249 } 250 251 /* Destroy old arrays */ 252 for (i=0; i < event->recsize; i++) { 253 ierr = PetscFree(event->recorder.eventidx[i]);CHKERRQ(ierr); 254 } 255 ierr = PetscFree(event->recorder.eventidx);CHKERRQ(ierr); 256 ierr = PetscFree(event->recorder.nevents);CHKERRQ(ierr); 257 ierr = PetscFree(event->recorder.stepnum);CHKERRQ(ierr); 258 ierr = PetscFree(event->recorder.time);CHKERRQ(ierr); 259 260 /* Set pointers */ 261 event->recorder.time = time; 262 event->recorder.stepnum = stepnum; 263 event->recorder.nevents = nevents; 264 event->recorder.eventidx = eventidx; 265 266 /* Double size */ 267 event->recsize *= fact; 268 269 PetscFunctionReturn(0); 270 } 271 272 /* 273 Helper rutine to handle user postenvents and recording 274 */ 275 #undef __FUNCT__ 276 #define __FUNCT__ "TSPostEvent" 277 static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U) 278 { 279 PetscErrorCode ierr; 280 TSEvent event = ts->event; 281 PetscBool terminate = PETSC_FALSE; 282 PetscInt i,ctr,stepnum; 283 PetscBool ts_terminate; 284 PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ 285 286 PetscFunctionBegin; 287 if (event->postevent) { 288 ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 289 } 290 291 /* Handle termination events */ 292 for (i=0; i < event->nevents_zero; i++) { 293 terminate = (PetscBool)(terminate || event->terminate[event->events_zero[i]]); 294 } 295 ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); 296 if (ts_terminate) { 297 ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr); 298 event->status = TSEVENT_NONE; 299 } else { 300 event->status = TSEVENT_RESET_NEXTSTEP; 301 } 302 303 event->ptime_prev = t; 304 /* Reset event residual functions as states might get changed by the postevent callback */ 305 if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);} 306 /* Cache current event residual functions */ 307 for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; 308 309 /* Record the event in the event recorder */ 310 ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr); 311 ctr = event->recorder.ctr; 312 if (ctr == event->recsize) { 313 ierr = TSEventRecorderResize(event);CHKERRQ(ierr); 314 } 315 event->recorder.time[ctr] = t; 316 event->recorder.stepnum[ctr] = stepnum; 317 event->recorder.nevents[ctr] = event->nevents_zero; 318 for (i=0; i < event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; 319 event->recorder.ctr++; 320 PetscFunctionReturn(0); 321 } 322 323 /* Uses Anderson-Bjorck variant of regular falsi method */ 324 PETSC_STATIC_INLINE PetscReal TSEventComputeStepSize(PetscReal tleft,PetscReal t, PetscReal tright, PetscScalar fleft, PetscScalar f, PetscScalar fright, PetscInt side,PetscReal dt) 325 { 326 PetscReal scal=1.0; 327 if(PetscRealPart(fleft)*PetscRealPart(f) < 0) { 328 if(side == 1) { 329 scal = (PetscRealPart(fright) - PetscRealPart(f))/PetscRealPart(fright); 330 if(scal < PETSC_SMALL) scal = 0.5; 331 } 332 dt = PetscMin(dt,(scal*PetscRealPart(fleft)*t - PetscRealPart(f)*tleft)/(scal*PetscRealPart(fleft) - PetscRealPart(f)) - tleft); 333 } else { 334 if(side == -1) { 335 scal = (PetscRealPart(fleft) - PetscRealPart(f))/PetscRealPart(fleft); 336 if(scal < PETSC_SMALL) scal = 0.5; 337 } 338 dt = PetscMin(dt,(PetscRealPart(f)*tright - scal*PetscRealPart(fright)*t)/(PetscRealPart(f) - scal*PetscRealPart(fright)) - t); 339 } 340 return dt; 341 } 342 343 344 #undef __FUNCT__ 345 #define __FUNCT__ "TSEventHandler" 346 PetscErrorCode TSEventHandler(TS ts) 347 { 348 PetscErrorCode ierr; 349 TSEvent event; 350 PetscReal t; 351 Vec U; 352 PetscInt i; 353 PetscReal dt,dt_min; 354 PetscInt rollback=0,in[2],out[2]; 355 PetscInt fvalue_sign,fvalueprev_sign; 356 357 PetscFunctionBegin; 358 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 359 if (!ts->event) PetscFunctionReturn(0); 360 event = ts->event; 361 362 ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 363 ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 364 ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 365 366 if (event->status == TSEVENT_NONE) { 367 event->timestep_prev = dt; 368 } 369 370 if (event->status == TSEVENT_RESET_NEXTSTEP) { 371 /* Restore previous time step */ 372 dt = event->timestep_prev; 373 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 374 event->status = TSEVENT_NONE; 375 } 376 377 if (event->status == TSEVENT_NONE) { 378 event->ptime_end = t; 379 } 380 381 ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); 382 383 for (i=0; i < event->nevents; i++) { 384 if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 385 event->status = TSEVENT_ZERO; 386 event->fvalue_right[i] = event->fvalue[i]; 387 continue; 388 } 389 fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); 390 fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); 391 if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) { 392 switch (event->direction[i]) { 393 case -1: 394 if (fvalue_sign < 0) { 395 rollback = 1; 396 397 /* Compute new time step */ 398 dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); 399 400 if (event->monitor) { 401 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); 402 } 403 event->fvalue_right[i] = event->fvalue[i]; 404 event->side[i] = 1; 405 406 if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 407 event->status = TSEVENT_LOCATED_INTERVAL; 408 } 409 break; 410 case 1: 411 if (fvalue_sign > 0) { 412 rollback = 1; 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 417 if (event->monitor) { 418 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); 419 } 420 event->fvalue_right[i] = event->fvalue[i]; 421 event->side[i] = 1; 422 423 if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 424 event->status = TSEVENT_LOCATED_INTERVAL; 425 } 426 break; 427 case 0: 428 rollback = 1; 429 430 /* Compute new time step */ 431 dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); 432 433 if (event->monitor) { 434 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); 435 } 436 event->fvalue_right[i] = event->fvalue[i]; 437 event->side[i] = 1; 438 439 if(!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; 440 event->status = TSEVENT_LOCATED_INTERVAL; 441 442 break; 443 } 444 } 445 } 446 447 in[0] = event->status; in[1] = rollback; 448 ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 449 event->status = (TSEventStatus)out[0]; rollback = out[1]; 450 if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; 451 452 event->nevents_zero = 0; 453 if (event->status == TSEVENT_ZERO) { 454 for (i=0; i < event->nevents; i++) { 455 if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { 456 event->events_zero[event->nevents_zero++] = i; 457 if (event->monitor) { 458 ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);CHKERRQ(ierr); 459 } 460 event->zerocrossing[i] = PETSC_FALSE; 461 } 462 event->side[i] = 0; 463 } 464 ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr); 465 466 dt = event->ptime_end - t; 467 if (PetscAbsReal(dt) < PETSC_SMALL) dt += event->timestep_prev; /* XXX Should be done better */ 468 ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 469 event->iterctr = 0; 470 PetscFunctionReturn(0); 471 } 472 473 if (event->status == TSEVENT_LOCATED_INTERVAL) { 474 ierr = TSRollBack(ts);CHKERRQ(ierr); 475 ts->steps--; 476 ts->total_steps--; 477 ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr); 478 event->status = TSEVENT_PROCESSING; 479 480 event->ptime_right = t; 481 } else { 482 for(i=0; i < event->nevents; i++) { 483 if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) { 484 /* Compute new time step */ 485 dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); 486 487 event->side[i] = -1; 488 } 489 event->fvalue_prev[i] = event->fvalue[i]; 490 } 491 492 if (event->monitor && event->status == TSEVENT_PROCESSING) { 493 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); 494 } 495 event->ptime_prev = t; 496 } 497 498 if(event->status == TSEVENT_PROCESSING) event->iterctr++; 499 500 ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 501 ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr); 502 PetscFunctionReturn(0); 503 } 504 505 #undef __FUNCT__ 506 #define __FUNCT__ "TSAdjointEventHandler" 507 PetscErrorCode TSAdjointEventHandler(TS ts) 508 { 509 PetscErrorCode ierr; 510 TSEvent event; 511 PetscReal t; 512 Vec U; 513 PetscInt ctr; 514 PetscBool forwardsolve=PETSC_FALSE; /* Flag indicating that TS is doing an adjoint solve */ 515 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 518 if (!ts->event) PetscFunctionReturn(0); 519 event = ts->event; 520 521 ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 522 ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 523 524 ctr = event->recorder.ctr-1; 525 if (ctr >= 0 && PetscAbsReal(t - event->recorder.time[ctr]) < PETSC_SMALL) { 526 /* Call the user postevent function */ 527 if (event->postevent) { 528 ierr = (*event->postevent)(ts,event->recorder.nevents[ctr],event->recorder.eventidx[ctr],t,U,forwardsolve,event->ctx);CHKERRQ(ierr); 529 event->recorder.ctr--; 530 } 531 } 532 533 PetscBarrier((PetscObject)ts); 534 PetscFunctionReturn(0); 535 } 536