1 2 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3 4 PetscClassId TSADAPT_CLASSID; 5 6 static PetscFunctionList TSAdaptList; 7 static PetscBool TSAdaptPackageInitialized; 8 static PetscBool TSAdaptRegisterAllCalled; 9 10 PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt); 11 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 12 PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt); 13 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 14 PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt); 15 PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt); 16 17 /*@C 18 TSAdaptRegister - adds a TSAdapt implementation 19 20 Not Collective 21 22 Input Parameters: 23 + name_scheme - name of user-defined adaptivity scheme 24 - routine_create - routine to create method context 25 26 Notes: 27 TSAdaptRegister() may be called multiple times to add several user-defined families. 28 29 Sample usage: 30 .vb 31 TSAdaptRegister("my_scheme",MySchemeCreate); 32 .ve 33 34 Then, your scheme can be chosen with the procedural interface via 35 $ TSAdaptSetType(ts,"my_scheme") 36 or at runtime via the option 37 $ -ts_adapt_type my_scheme 38 39 Level: advanced 40 41 .seealso: `TSAdaptRegisterAll()` 42 @*/ 43 PetscErrorCode TSAdaptRegister(const char sname[], PetscErrorCode (*function)(TSAdapt)) 44 { 45 PetscFunctionBegin; 46 PetscCall(TSAdaptInitializePackage()); 47 PetscCall(PetscFunctionListAdd(&TSAdaptList, sname, function)); 48 PetscFunctionReturn(0); 49 } 50 51 /*@C 52 TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 53 54 Not Collective 55 56 Level: advanced 57 58 .seealso: `TSAdaptRegisterDestroy()` 59 @*/ 60 PetscErrorCode TSAdaptRegisterAll(void) 61 { 62 PetscFunctionBegin; 63 if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 64 TSAdaptRegisterAllCalled = PETSC_TRUE; 65 PetscCall(TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None)); 66 PetscCall(TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic)); 67 PetscCall(TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP)); 68 PetscCall(TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL)); 69 PetscCall(TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE)); 70 PetscCall(TSAdaptRegister(TSADAPTHISTORY, TSAdaptCreate_History)); 71 PetscFunctionReturn(0); 72 } 73 74 /*@C 75 TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 76 called from PetscFinalize(). 77 78 Level: developer 79 80 .seealso: `PetscFinalize()` 81 @*/ 82 PetscErrorCode TSAdaptFinalizePackage(void) 83 { 84 PetscFunctionBegin; 85 PetscCall(PetscFunctionListDestroy(&TSAdaptList)); 86 TSAdaptPackageInitialized = PETSC_FALSE; 87 TSAdaptRegisterAllCalled = PETSC_FALSE; 88 PetscFunctionReturn(0); 89 } 90 91 /*@C 92 TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 93 called from TSInitializePackage(). 94 95 Level: developer 96 97 .seealso: `PetscInitialize()` 98 @*/ 99 PetscErrorCode TSAdaptInitializePackage(void) 100 { 101 PetscFunctionBegin; 102 if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 103 TSAdaptPackageInitialized = PETSC_TRUE; 104 PetscCall(PetscClassIdRegister("TSAdapt", &TSADAPT_CLASSID)); 105 PetscCall(TSAdaptRegisterAll()); 106 PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage)); 107 PetscFunctionReturn(0); 108 } 109 110 /*@C 111 TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 112 113 Logicially Collective on TSAdapt 114 115 Input Parameters: 116 + adapt - the TS adapter, most likely obtained with TSGetAdapt() 117 - type - either TSADAPTBASIC or TSADAPTNONE 118 119 Options Database: 120 . -ts_adapt_type <basic or dsp or none> - to set the adapter type 121 122 Level: intermediate 123 124 .seealso: `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()` 125 @*/ 126 PetscErrorCode TSAdaptSetType(TSAdapt adapt, TSAdaptType type) 127 { 128 PetscBool match; 129 PetscErrorCode (*r)(TSAdapt); 130 131 PetscFunctionBegin; 132 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 133 PetscValidCharPointer(type, 2); 134 PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match)); 135 if (match) PetscFunctionReturn(0); 136 PetscCall(PetscFunctionListFind(TSAdaptList, type, &r)); 137 PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSAdapt type \"%s\" given", type); 138 PetscTryTypeMethod(adapt, destroy); 139 PetscCall(PetscMemzero(adapt->ops, sizeof(struct _TSAdaptOps))); 140 PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type)); 141 PetscCall((*r)(adapt)); 142 PetscFunctionReturn(0); 143 } 144 145 /*@C 146 TSAdaptGetType - gets the TS adapter method type (as a string). 147 148 Not Collective 149 150 Input Parameter: 151 . adapt - The TS adapter, most likely obtained with TSGetAdapt() 152 153 Output Parameter: 154 . type - The name of TS adapter method 155 156 Level: intermediate 157 158 .seealso `TSAdaptSetType()` 159 @*/ 160 PetscErrorCode TSAdaptGetType(TSAdapt adapt, TSAdaptType *type) 161 { 162 PetscFunctionBegin; 163 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 164 PetscValidPointer(type, 2); 165 *type = ((PetscObject)adapt)->type_name; 166 PetscFunctionReturn(0); 167 } 168 169 PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt, const char prefix[]) 170 { 171 PetscFunctionBegin; 172 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 173 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix)); 174 PetscFunctionReturn(0); 175 } 176 177 /*@C 178 TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 179 180 Collective on PetscViewer 181 182 Input Parameters: 183 + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 184 some related function before a call to TSAdaptLoad(). 185 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 186 HDF5 file viewer, obtained from PetscViewerHDF5Open() 187 188 Level: intermediate 189 190 Notes: 191 The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 192 193 Notes for advanced users: 194 Most users should not need to know the details of the binary storage 195 format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 196 But for anyone who's interested, the standard binary matrix storage 197 format is 198 .vb 199 has not yet been determined 200 .ve 201 202 .seealso: `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()` 203 @*/ 204 PetscErrorCode TSAdaptLoad(TSAdapt adapt, PetscViewer viewer) 205 { 206 PetscBool isbinary; 207 char type[256]; 208 209 PetscFunctionBegin; 210 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 211 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 212 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 213 PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 214 215 PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 216 PetscCall(TSAdaptSetType(adapt, type)); 217 PetscTryTypeMethod(adapt, load, viewer); 218 PetscFunctionReturn(0); 219 } 220 221 PetscErrorCode TSAdaptView(TSAdapt adapt, PetscViewer viewer) 222 { 223 PetscBool iascii, isbinary, isnone, isglee; 224 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 227 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt), &viewer)); 228 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 229 PetscCheckSameComm(adapt, 1, viewer, 2); 230 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 231 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 232 if (iascii) { 233 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer)); 234 PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTNONE, &isnone)); 235 PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTGLEE, &isglee)); 236 if (!isnone) { 237 if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer, " always accepting steps\n")); 238 PetscCall(PetscViewerASCIIPrintf(viewer, " safety factor %g\n", (double)adapt->safety)); 239 PetscCall(PetscViewerASCIIPrintf(viewer, " extra safety factor after step rejection %g\n", (double)adapt->reject_safety)); 240 PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest increase %g\n", (double)adapt->clip[1])); 241 PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest decrease %g\n", (double)adapt->clip[0])); 242 PetscCall(PetscViewerASCIIPrintf(viewer, " maximum allowed timestep %g\n", (double)adapt->dt_max)); 243 PetscCall(PetscViewerASCIIPrintf(viewer, " minimum allowed timestep %g\n", (double)adapt->dt_min)); 244 PetscCall(PetscViewerASCIIPrintf(viewer, " maximum solution absolute value to be ignored %g\n", (double)adapt->ignore_max)); 245 } 246 if (isglee) { 247 if (adapt->glee_use_local) { 248 PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses local error control\n")); 249 } else { 250 PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses global error control\n")); 251 } 252 } 253 PetscCall(PetscViewerASCIIPushTab(viewer)); 254 PetscTryTypeMethod(adapt, view, viewer); 255 PetscCall(PetscViewerASCIIPopTab(viewer)); 256 } else if (isbinary) { 257 char type[256]; 258 259 /* need to save FILE_CLASS_ID for adapt class */ 260 PetscCall(PetscStrncpy(type, ((PetscObject)adapt)->type_name, 256)); 261 PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 262 } else PetscTryTypeMethod(adapt, view, viewer); 263 PetscFunctionReturn(0); 264 } 265 266 /*@ 267 TSAdaptReset - Resets a TSAdapt context. 268 269 Collective on TS 270 271 Input Parameter: 272 . adapt - the TSAdapt context obtained from TSAdaptCreate() 273 274 Level: developer 275 276 .seealso: `TSAdaptCreate()`, `TSAdaptDestroy()` 277 @*/ 278 PetscErrorCode TSAdaptReset(TSAdapt adapt) 279 { 280 PetscFunctionBegin; 281 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 282 PetscTryTypeMethod(adapt, reset); 283 PetscFunctionReturn(0); 284 } 285 286 PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 287 { 288 PetscFunctionBegin; 289 if (!*adapt) PetscFunctionReturn(0); 290 PetscValidHeaderSpecific(*adapt, TSADAPT_CLASSID, 1); 291 if (--((PetscObject)(*adapt))->refct > 0) { 292 *adapt = NULL; 293 PetscFunctionReturn(0); 294 } 295 296 PetscCall(TSAdaptReset(*adapt)); 297 298 PetscTryTypeMethod((*adapt), destroy); 299 PetscCall(PetscViewerDestroy(&(*adapt)->monitor)); 300 PetscCall(PetscHeaderDestroy(adapt)); 301 PetscFunctionReturn(0); 302 } 303 304 /*@ 305 TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 306 307 Collective on TSAdapt 308 309 Input Parameters: 310 + adapt - adaptive controller context 311 - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 312 313 Options Database Keys: 314 . -ts_adapt_monitor - to turn on monitoring 315 316 Level: intermediate 317 318 .seealso: `TSAdaptChoose()` 319 @*/ 320 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt, PetscBool flg) 321 { 322 PetscFunctionBegin; 323 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 324 PetscValidLogicalCollectiveBool(adapt, flg, 2); 325 if (flg) { 326 if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt), "stdout", &adapt->monitor)); 327 } else { 328 PetscCall(PetscViewerDestroy(&adapt->monitor)); 329 } 330 PetscFunctionReturn(0); 331 } 332 333 /*@C 334 TSAdaptSetCheckStage - Set a callback to check convergence for a stage 335 336 Logically collective on TSAdapt 337 338 Input Parameters: 339 + adapt - adaptive controller context 340 - func - stage check function 341 342 Arguments of func: 343 $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 344 345 + adapt - adaptive controller context 346 . ts - time stepping context 347 - accept - pending choice of whether to accept, can be modified by this routine 348 349 Level: advanced 350 351 .seealso: `TSAdaptChoose()` 352 @*/ 353 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt, TS, PetscReal, Vec, PetscBool *)) 354 { 355 PetscFunctionBegin; 356 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 357 adapt->checkstage = func; 358 PetscFunctionReturn(0); 359 } 360 361 /*@ 362 TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 363 any error or stability condition not meeting the prescribed goal. 364 365 Logically collective on TSAdapt 366 367 Input Parameters: 368 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 369 - flag - whether to always accept steps 370 371 Options Database Keys: 372 . -ts_adapt_always_accept - to always accept steps 373 374 Level: intermediate 375 376 .seealso: `TSAdapt`, `TSAdaptChoose()` 377 @*/ 378 PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag) 379 { 380 PetscFunctionBegin; 381 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 382 PetscValidLogicalCollectiveBool(adapt, flag, 2); 383 adapt->always_accept = flag; 384 PetscFunctionReturn(0); 385 } 386 387 /*@ 388 TSAdaptSetSafety - Set safety factors 389 390 Logically collective on TSAdapt 391 392 Input Parameters: 393 + adapt - adaptive controller context 394 . safety - safety factor relative to target error/stability goal 395 - reject_safety - extra safety factor to apply if the last step was rejected 396 397 Options Database Keys: 398 + -ts_adapt_safety <safety> - to set safety factor 399 - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 400 401 Level: intermediate 402 403 .seealso: `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()` 404 @*/ 405 PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety) 406 { 407 PetscFunctionBegin; 408 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 409 PetscValidLogicalCollectiveReal(adapt, safety, 2); 410 PetscValidLogicalCollectiveReal(adapt, reject_safety, 3); 411 PetscCheck(safety == PETSC_DEFAULT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety); 412 PetscCheck(safety == PETSC_DEFAULT || safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be less than one", (double)safety); 413 PetscCheck(reject_safety == PETSC_DEFAULT || reject_safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be non negative", (double)reject_safety); 414 PetscCheck(reject_safety == PETSC_DEFAULT || reject_safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be less than one", (double)reject_safety); 415 if (safety != PETSC_DEFAULT) adapt->safety = safety; 416 if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 417 PetscFunctionReturn(0); 418 } 419 420 /*@ 421 TSAdaptGetSafety - Get safety factors 422 423 Not Collective 424 425 Input Parameter: 426 . adapt - adaptive controller context 427 428 Output Parameters: 429 . safety - safety factor relative to target error/stability goal 430 + reject_safety - extra safety factor to apply if the last step was rejected 431 432 Level: intermediate 433 434 .seealso: `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()` 435 @*/ 436 PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety) 437 { 438 PetscFunctionBegin; 439 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 440 if (safety) PetscValidRealPointer(safety, 2); 441 if (reject_safety) PetscValidRealPointer(reject_safety, 3); 442 if (safety) *safety = adapt->safety; 443 if (reject_safety) *reject_safety = adapt->reject_safety; 444 PetscFunctionReturn(0); 445 } 446 447 /*@ 448 TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components. 449 450 Logically collective on TSAdapt 451 452 Input Parameters: 453 + adapt - adaptive controller context 454 - max_ignore - threshold for solution components that are ignored during error estimation 455 456 Options Database Keys: 457 . -ts_adapt_max_ignore <max_ignore> - to set the threshold 458 459 Level: intermediate 460 461 .seealso: `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()` 462 @*/ 463 PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore) 464 { 465 PetscFunctionBegin; 466 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 467 PetscValidLogicalCollectiveReal(adapt, max_ignore, 2); 468 adapt->ignore_max = max_ignore; 469 PetscFunctionReturn(0); 470 } 471 472 /*@ 473 TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value). 474 475 Not Collective 476 477 Input Parameter: 478 . adapt - adaptive controller context 479 480 Output Parameter: 481 . max_ignore - threshold for solution components that are ignored during error estimation 482 483 Level: intermediate 484 485 .seealso: `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()` 486 @*/ 487 PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore) 488 { 489 PetscFunctionBegin; 490 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 491 PetscValidRealPointer(max_ignore, 2); 492 *max_ignore = adapt->ignore_max; 493 PetscFunctionReturn(0); 494 } 495 496 /*@ 497 TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 498 499 Logically collective on TSAdapt 500 501 Input Parameters: 502 + adapt - adaptive controller context 503 . low - admissible decrease factor 504 - high - admissible increase factor 505 506 Options Database Keys: 507 . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 508 509 Level: intermediate 510 511 .seealso: `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()` 512 @*/ 513 PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high) 514 { 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 517 PetscValidLogicalCollectiveReal(adapt, low, 2); 518 PetscValidLogicalCollectiveReal(adapt, high, 3); 519 PetscCheck(low == PETSC_DEFAULT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low); 520 PetscCheck(low == PETSC_DEFAULT || low <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be less than one", (double)low); 521 PetscCheck(high == PETSC_DEFAULT || high >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Increase factor %g must be greater than one", (double)high); 522 if (low != PETSC_DEFAULT) adapt->clip[0] = low; 523 if (high != PETSC_DEFAULT) adapt->clip[1] = high; 524 PetscFunctionReturn(0); 525 } 526 527 /*@ 528 TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 529 530 Not Collective 531 532 Input Parameter: 533 . adapt - adaptive controller context 534 535 Output Parameters: 536 + low - optional, admissible decrease factor 537 - high - optional, admissible increase factor 538 539 Level: intermediate 540 541 .seealso: `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()` 542 @*/ 543 PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high) 544 { 545 PetscFunctionBegin; 546 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 547 if (low) PetscValidRealPointer(low, 2); 548 if (high) PetscValidRealPointer(high, 3); 549 if (low) *low = adapt->clip[0]; 550 if (high) *high = adapt->clip[1]; 551 PetscFunctionReturn(0); 552 } 553 554 /*@ 555 TSAdaptSetScaleSolveFailed - Scale step by this factor if solve fails 556 557 Logically collective on TSAdapt 558 559 Input Parameters: 560 + adapt - adaptive controller context 561 - scale - scale 562 563 Options Database Keys: 564 . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 565 566 Level: intermediate 567 568 .seealso: `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()` 569 @*/ 570 PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale) 571 { 572 PetscFunctionBegin; 573 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 574 PetscValidLogicalCollectiveReal(adapt, scale, 2); 575 PetscCheck(scale == PETSC_DEFAULT || scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale); 576 PetscCheck(scale == PETSC_DEFAULT || scale <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be less than one", (double)scale); 577 if (scale != PETSC_DEFAULT) adapt->scale_solve_failed = scale; 578 PetscFunctionReturn(0); 579 } 580 581 /*@ 582 TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 583 584 Not Collective 585 586 Input Parameter: 587 . adapt - adaptive controller context 588 589 Output Parameter: 590 . scale - scale factor 591 592 Level: intermediate 593 594 .seealso: `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()` 595 @*/ 596 PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale) 597 { 598 PetscFunctionBegin; 599 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 600 if (scale) PetscValidRealPointer(scale, 2); 601 if (scale) *scale = adapt->scale_solve_failed; 602 PetscFunctionReturn(0); 603 } 604 605 /*@ 606 TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 607 608 Logically collective on TSAdapt 609 610 Input Parameters: 611 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 612 . hmin - minimum time step 613 - hmax - maximum time step 614 615 Options Database Keys: 616 + -ts_adapt_dt_min <min> - to set minimum time step 617 - -ts_adapt_dt_max <max> - to set maximum time step 618 619 Level: intermediate 620 621 .seealso: `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()` 622 @*/ 623 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax) 624 { 625 PetscFunctionBegin; 626 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 627 PetscValidLogicalCollectiveReal(adapt, hmin, 2); 628 PetscValidLogicalCollectiveReal(adapt, hmax, 3); 629 PetscCheck(hmin == PETSC_DEFAULT || hmin >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmin); 630 PetscCheck(hmax == PETSC_DEFAULT || hmax >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmax); 631 if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 632 if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 633 hmin = adapt->dt_min; 634 hmax = adapt->dt_max; 635 PetscCheck(hmax > hmin, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Maximum time step %g must greater than minimum time step %g", (double)hmax, (double)hmin); 636 PetscFunctionReturn(0); 637 } 638 639 /*@ 640 TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 641 642 Not Collective 643 644 Input Parameter: 645 . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 646 647 Output Parameters: 648 + hmin - minimum time step 649 - hmax - maximum time step 650 651 Level: intermediate 652 653 .seealso: `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()` 654 @*/ 655 PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax) 656 { 657 PetscFunctionBegin; 658 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 659 if (hmin) PetscValidRealPointer(hmin, 2); 660 if (hmax) PetscValidRealPointer(hmax, 3); 661 if (hmin) *hmin = adapt->dt_min; 662 if (hmax) *hmax = adapt->dt_max; 663 PetscFunctionReturn(0); 664 } 665 666 /* 667 TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 668 669 Collective on TSAdapt 670 671 Input Parameter: 672 . adapt - the TSAdapt context 673 674 Options Database Keys: 675 + -ts_adapt_type <type> - algorithm to use for adaptivity 676 . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 677 . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 678 . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 679 . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 680 . -ts_adapt_dt_min <min> - minimum timestep to use 681 . -ts_adapt_dt_max <max> - maximum timestep to use 682 . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 683 . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 684 - -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver 685 686 Level: advanced 687 688 Notes: 689 This function is automatically called by TSSetFromOptions() 690 691 .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`, 692 `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()` 693 */ 694 PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems *PetscOptionsObject) 695 { 696 char type[256] = TSADAPTBASIC; 697 PetscReal safety, reject_safety, clip[2], scale, hmin, hmax; 698 PetscBool set, flg; 699 PetscInt two; 700 701 PetscFunctionBegin; 702 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 703 /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 704 * function can only be called from inside TSSetFromOptions() */ 705 PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options"); 706 PetscCall(PetscOptionsFList("-ts_adapt_type", "Algorithm to use for adaptivity", "TSAdaptSetType", TSAdaptList, ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type, type, sizeof(type), &flg)); 707 if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type)); 708 709 PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set)); 710 if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg)); 711 712 safety = adapt->safety; 713 reject_safety = adapt->reject_safety; 714 PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set)); 715 PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg)); 716 if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety)); 717 718 two = 2; 719 clip[0] = adapt->clip[0]; 720 clip[1] = adapt->clip[1]; 721 PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set)); 722 PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip"); 723 if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1])); 724 725 hmin = adapt->dt_min; 726 hmax = adapt->dt_max; 727 PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set)); 728 PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg)); 729 if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax)); 730 731 PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set)); 732 PetscCall(PetscOptionsBool("-ts_adapt_glee_use_local", "GLEE adaptor uses local error estimation for step control", "", adapt->glee_use_local, &adapt->glee_use_local, &set)); 733 734 PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set)); 735 if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale)); 736 737 PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL)); 738 PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported"); 739 740 PetscCall(PetscOptionsInt("-ts_adapt_time_step_increase_delay", "Number of timesteps to delay increasing the time step after it has been decreased due to failed solver", "TSAdaptSetTimeStepIncreaseDelay", adapt->timestepjustdecreased_delay, &adapt->timestepjustdecreased_delay, NULL)); 741 742 PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set)); 743 if (set) PetscCall(TSAdaptSetMonitor(adapt, flg)); 744 745 PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject); 746 PetscOptionsHeadEnd(); 747 PetscFunctionReturn(0); 748 } 749 750 /*@ 751 TSAdaptCandidatesClear - clear any previously set candidate schemes 752 753 Logically collective on TSAdapt 754 755 Input Parameter: 756 . adapt - adaptive controller 757 758 Level: developer 759 760 .seealso: `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 761 @*/ 762 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 763 { 764 PetscFunctionBegin; 765 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 766 PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates))); 767 PetscFunctionReturn(0); 768 } 769 770 /*@C 771 TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 772 773 Logically collective on TSAdapt 774 775 Input Parameters: 776 + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 777 . name - name of the candidate scheme to add 778 . order - order of the candidate scheme 779 . stageorder - stage order of the candidate scheme 780 . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 781 . cost - relative measure of the amount of work required for the candidate scheme 782 - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 783 784 Note: 785 This routine is not available in Fortran. 786 787 Level: developer 788 789 .seealso: `TSAdaptCandidatesClear()`, `TSAdaptChoose()` 790 @*/ 791 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse) 792 { 793 PetscInt c; 794 795 PetscFunctionBegin; 796 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 797 PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order); 798 if (inuse) { 799 PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 800 adapt->candidates.inuse_set = PETSC_TRUE; 801 } 802 /* first slot if this is the current scheme, otherwise the next available slot */ 803 c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 804 805 adapt->candidates.name[c] = name; 806 adapt->candidates.order[c] = order; 807 adapt->candidates.stageorder[c] = stageorder; 808 adapt->candidates.ccfl[c] = ccfl; 809 adapt->candidates.cost[c] = cost; 810 adapt->candidates.n++; 811 PetscFunctionReturn(0); 812 } 813 814 /*@C 815 TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 816 817 Not Collective 818 819 Input Parameter: 820 . adapt - time step adaptivity context 821 822 Output Parameters: 823 + n - number of candidate schemes, always at least 1 824 . order - the order of each candidate scheme 825 . stageorder - the stage order of each candidate scheme 826 . ccfl - the CFL coefficient of each scheme 827 - cost - the relative cost of each scheme 828 829 Level: developer 830 831 Note: 832 The current scheme is always returned in the first slot 833 834 .seealso: `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 835 @*/ 836 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost) 837 { 838 PetscFunctionBegin; 839 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 840 if (n) *n = adapt->candidates.n; 841 if (order) *order = adapt->candidates.order; 842 if (stageorder) *stageorder = adapt->candidates.stageorder; 843 if (ccfl) *ccfl = adapt->candidates.ccfl; 844 if (cost) *cost = adapt->candidates.cost; 845 PetscFunctionReturn(0); 846 } 847 848 /*@C 849 TSAdaptChoose - choose which method and step size to use for the next step 850 851 Collective on TSAdapt 852 853 Input Parameters: 854 + adapt - adaptive contoller 855 . ts - time stepper 856 - h - current step size 857 858 Output Parameters: 859 + next_sc - optional, scheme to use for the next step 860 . next_h - step size to use for the next step 861 - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 862 863 Note: 864 The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 865 being retried after an initial rejection. 866 867 Level: developer 868 869 .seealso: `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()` 870 @*/ 871 PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept) 872 { 873 PetscInt ncandidates = adapt->candidates.n; 874 PetscInt scheme = 0; 875 PetscReal wlte = -1.0; 876 PetscReal wltea = -1.0; 877 PetscReal wlter = -1.0; 878 879 PetscFunctionBegin; 880 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 881 PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 882 if (next_sc) PetscValidIntPointer(next_sc, 4); 883 PetscValidRealPointer(next_h, 5); 884 PetscValidBoolPointer(accept, 6); 885 if (next_sc) *next_sc = 0; 886 887 /* Do not mess with adaptivity while handling events*/ 888 if (ts->event && ts->event->status != TSEVENT_NONE) { 889 *next_h = h; 890 *accept = PETSC_TRUE; 891 PetscFunctionReturn(0); 892 } 893 894 PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter); 895 PetscCheck(scheme >= 0 && (ncandidates <= 0 || scheme < ncandidates), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Chosen scheme %" PetscInt_FMT " not in valid range 0..%" PetscInt_FMT, scheme, ncandidates - 1); 896 PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h); 897 if (next_sc) *next_sc = scheme; 898 899 if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 900 /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 901 PetscReal t = ts->ptime + ts->time_step, h = *next_h; 902 PetscReal tend = t + h, tmax, hmax; 903 PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 904 PetscReal b = adapt->matchstepfac[1]; 905 906 if (ts->tspan) { 907 if (PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) /* hit a span time point */ 908 if (ts->tspan->spanctr + 1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr + 1]; 909 else tmax = ts->max_time; /* hit the last span time point */ 910 else tmax = ts->tspan->span_times[ts->tspan->spanctr]; 911 } else tmax = ts->max_time; 912 hmax = tmax - t; 913 if (t < tmax && tend > tmax) *next_h = hmax; 914 if (t < tmax && tend < tmax && h * b > hmax) *next_h = hmax / 2; 915 if (t < tmax && tend < tmax && h * a > hmax) *next_h = hmax; 916 /* if step size is changed to match a span time point */ 917 if (ts->tspan && h != *next_h && !adapt->dt_span_cached) adapt->dt_span_cached = h; 918 /* reset time step after a span time point */ 919 if (ts->tspan && h == *next_h && adapt->dt_span_cached && PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) { 920 *next_h = adapt->dt_span_cached; 921 adapt->dt_span_cached = 0; 922 } 923 } 924 if (adapt->monitor) { 925 const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 926 PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 927 if (wlte < 0) { 928 PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s %s %" PetscInt_FMT ":%s step %3" PetscInt_FMT " %s t=%-11g+%10.3e dt=%-10.3e\n", ((PetscObject)adapt)->type_name, ((PetscObject)ts)->type_name, scheme, sc_name, ts->steps, *accept ? "accepted" : "rejected", 929 (double)ts->ptime, (double)h, (double)*next_h)); 930 } else { 931 PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s %s %" PetscInt_FMT ":%s step %3" PetscInt_FMT " %s t=%-11g+%10.3e dt=%-10.3e wlte=%5.3g wltea=%5.3g wlter=%5.3g\n", ((PetscObject)adapt)->type_name, ((PetscObject)ts)->type_name, scheme, sc_name, ts->steps, *accept ? "accepted" : "rejected", 932 (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter)); 933 } 934 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 935 } 936 PetscFunctionReturn(0); 937 } 938 939 /*@ 940 TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 941 before increasing the time step. 942 943 Logicially Collective on TSAdapt 944 945 Input Parameters: 946 + adapt - adaptive controller context 947 - cnt - the number of timesteps 948 949 Options Database Key: 950 . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 951 952 Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 953 The successful use of this option is problem dependent 954 955 Developer Note: there is no theory to support this option 956 957 Level: advanced 958 959 .seealso: 960 @*/ 961 PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt) 962 { 963 PetscFunctionBegin; 964 adapt->timestepjustdecreased_delay = cnt; 965 PetscFunctionReturn(0); 966 } 967 968 /*@ 969 TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails or solution vector is infeasible) 970 971 Collective on TSAdapt 972 973 Input Parameters: 974 + adapt - adaptive controller context 975 . ts - time stepper 976 . t - Current simulation time 977 - Y - Current solution vector 978 979 Output Parameter: 980 . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 981 982 Level: developer 983 984 .seealso: 985 @*/ 986 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept) 987 { 988 SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 989 990 PetscFunctionBegin; 991 PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 992 PetscValidHeaderSpecific(ts, TS_CLASSID, 2); 993 PetscValidBoolPointer(accept, 5); 994 995 if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason)); 996 if (snesreason < 0) { 997 *accept = PETSC_FALSE; 998 if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 999 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 1000 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures)); 1001 if (adapt->monitor) { 1002 PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 1003 PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected t=%-11g+%10.3e, nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed\n", ((PetscObject)adapt)->type_name, ts->steps, 1004 (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures)); 1005 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 1006 } 1007 } 1008 } else { 1009 *accept = PETSC_TRUE; 1010 PetscCall(TSFunctionDomainError(ts, t, Y, accept)); 1011 if (*accept && adapt->checkstage) { 1012 PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept)); 1013 if (!*accept) { 1014 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by user function provided by TSSetFunctionDomainError()\n", ts->steps)); 1015 if (adapt->monitor) { 1016 PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 1017 PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected by user function provided by TSSetFunctionDomainError()\n", ((PetscObject)adapt)->type_name, ts->steps)); 1018 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 1019 } 1020 } 1021 } 1022 } 1023 1024 if (!(*accept) && !ts->reason) { 1025 PetscReal dt, new_dt; 1026 PetscCall(TSGetTimeStep(ts, &dt)); 1027 new_dt = dt * adapt->scale_solve_failed; 1028 PetscCall(TSSetTimeStep(ts, new_dt)); 1029 adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 1030 if (adapt->monitor) { 1031 PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 1032 PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected (%s) t=%-11g+%10.3e retrying with dt=%-10.3e\n", ((PetscObject)adapt)->type_name, ts->steps, SNESConvergedReasons[snesreason], (double)ts->ptime, (double)dt, (double)new_dt)); 1033 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel)); 1034 } 1035 } 1036 PetscFunctionReturn(0); 1037 } 1038 1039 /*@ 1040 TSAdaptCreate - create an adaptive controller context for time stepping 1041 1042 Collective 1043 1044 Input Parameter: 1045 . comm - The communicator 1046 1047 Output Parameter: 1048 . adapt - new TSAdapt object 1049 1050 Level: developer 1051 1052 Notes: 1053 TSAdapt creation is handled by TS, so users should not need to call this function. 1054 1055 .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()` 1056 @*/ 1057 PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt) 1058 { 1059 TSAdapt adapt; 1060 1061 PetscFunctionBegin; 1062 PetscValidPointer(inadapt, 2); 1063 *inadapt = NULL; 1064 PetscCall(TSAdaptInitializePackage()); 1065 1066 PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView)); 1067 1068 adapt->always_accept = PETSC_FALSE; 1069 adapt->safety = 0.9; 1070 adapt->reject_safety = 0.5; 1071 adapt->clip[0] = 0.1; 1072 adapt->clip[1] = 10.; 1073 adapt->dt_min = 1e-20; 1074 adapt->dt_max = 1e+20; 1075 adapt->ignore_max = -1.0; 1076 adapt->glee_use_local = PETSC_TRUE; 1077 adapt->scale_solve_failed = 0.25; 1078 /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1079 to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1080 adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1081 adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 1082 adapt->wnormtype = NORM_2; 1083 adapt->timestepjustdecreased_delay = 0; 1084 1085 *inadapt = adapt; 1086 PetscFunctionReturn(0); 1087 } 1088