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 if (adapt->ops->destroy) PetscCall((*adapt->ops->destroy)(adapt)); 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 if (adapt->ops->load) PetscCall((*adapt->ops->load)(adapt,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 if (adapt->ops->view) { 254 PetscCall(PetscViewerASCIIPushTab(viewer)); 255 PetscCall((*adapt->ops->view)(adapt,viewer)); 256 PetscCall(PetscViewerASCIIPopTab(viewer)); 257 } 258 } else if (isbinary) { 259 char type[256]; 260 261 /* need to save FILE_CLASS_ID for adapt class */ 262 PetscCall(PetscStrncpy(type,((PetscObject)adapt)->type_name,256)); 263 PetscCall(PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR)); 264 } else if (adapt->ops->view) PetscCall((*adapt->ops->view)(adapt,viewer)); 265 PetscFunctionReturn(0); 266 } 267 268 /*@ 269 TSAdaptReset - Resets a TSAdapt context. 270 271 Collective on TS 272 273 Input Parameter: 274 . adapt - the TSAdapt context obtained from TSAdaptCreate() 275 276 Level: developer 277 278 .seealso: `TSAdaptCreate()`, `TSAdaptDestroy()` 279 @*/ 280 PetscErrorCode TSAdaptReset(TSAdapt adapt) 281 { 282 PetscFunctionBegin; 283 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 284 if (adapt->ops->reset) PetscCall((*adapt->ops->reset)(adapt)); 285 PetscFunctionReturn(0); 286 } 287 288 PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 289 { 290 PetscFunctionBegin; 291 if (!*adapt) PetscFunctionReturn(0); 292 PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 293 if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 294 295 PetscCall(TSAdaptReset(*adapt)); 296 297 if ((*adapt)->ops->destroy) PetscCall((*(*adapt)->ops->destroy)(*adapt)); 298 PetscCall(PetscViewerDestroy(&(*adapt)->monitor)); 299 PetscCall(PetscHeaderDestroy(adapt)); 300 PetscFunctionReturn(0); 301 } 302 303 /*@ 304 TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 305 306 Collective on TSAdapt 307 308 Input Parameters: 309 + adapt - adaptive controller context 310 - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 311 312 Options Database Keys: 313 . -ts_adapt_monitor - to turn on monitoring 314 315 Level: intermediate 316 317 .seealso: `TSAdaptChoose()` 318 @*/ 319 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 320 { 321 PetscFunctionBegin; 322 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 323 PetscValidLogicalCollectiveBool(adapt,flg,2); 324 if (flg) { 325 if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor)); 326 } else { 327 PetscCall(PetscViewerDestroy(&adapt->monitor)); 328 } 329 PetscFunctionReturn(0); 330 } 331 332 /*@C 333 TSAdaptSetCheckStage - Set a callback to check convergence for a stage 334 335 Logically collective on TSAdapt 336 337 Input Parameters: 338 + adapt - adaptive controller context 339 - func - stage check function 340 341 Arguments of func: 342 $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 343 344 + adapt - adaptive controller context 345 . ts - time stepping context 346 - accept - pending choice of whether to accept, can be modified by this routine 347 348 Level: advanced 349 350 .seealso: `TSAdaptChoose()` 351 @*/ 352 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 353 { 354 PetscFunctionBegin; 355 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 356 adapt->checkstage = func; 357 PetscFunctionReturn(0); 358 } 359 360 /*@ 361 TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 362 any error or stability condition not meeting the prescribed goal. 363 364 Logically collective on TSAdapt 365 366 Input Parameters: 367 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 368 - flag - whether to always accept steps 369 370 Options Database Keys: 371 . -ts_adapt_always_accept - to always accept steps 372 373 Level: intermediate 374 375 .seealso: `TSAdapt`, `TSAdaptChoose()` 376 @*/ 377 PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 378 { 379 PetscFunctionBegin; 380 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 381 PetscValidLogicalCollectiveBool(adapt,flag,2); 382 adapt->always_accept = flag; 383 PetscFunctionReturn(0); 384 } 385 386 /*@ 387 TSAdaptSetSafety - Set safety factors 388 389 Logically collective on TSAdapt 390 391 Input Parameters: 392 + adapt - adaptive controller context 393 . safety - safety factor relative to target error/stability goal 394 - reject_safety - extra safety factor to apply if the last step was rejected 395 396 Options Database Keys: 397 + -ts_adapt_safety <safety> - to set safety factor 398 - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 399 400 Level: intermediate 401 402 .seealso: `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()` 403 @*/ 404 PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 405 { 406 PetscFunctionBegin; 407 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 408 PetscValidLogicalCollectiveReal(adapt,safety,2); 409 PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 410 PetscCheck(safety == PETSC_DEFAULT || safety >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 411 PetscCheck(safety == PETSC_DEFAULT || safety <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be less than one",(double)safety); 412 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); 413 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); 414 if (safety != PETSC_DEFAULT) adapt->safety = safety; 415 if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 416 PetscFunctionReturn(0); 417 } 418 419 /*@ 420 TSAdaptGetSafety - Get safety factors 421 422 Not Collective 423 424 Input Parameter: 425 . adapt - adaptive controller context 426 427 Output Parameters: 428 . safety - safety factor relative to target error/stability goal 429 + reject_safety - extra safety factor to apply if the last step was rejected 430 431 Level: intermediate 432 433 .seealso: `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()` 434 @*/ 435 PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 436 { 437 PetscFunctionBegin; 438 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 439 if (safety) PetscValidRealPointer(safety,2); 440 if (reject_safety) PetscValidRealPointer(reject_safety,3); 441 if (safety) *safety = adapt->safety; 442 if (reject_safety) *reject_safety = adapt->reject_safety; 443 PetscFunctionReturn(0); 444 } 445 446 /*@ 447 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. 448 449 Logically collective on TSAdapt 450 451 Input Parameters: 452 + adapt - adaptive controller context 453 - max_ignore - threshold for solution components that are ignored during error estimation 454 455 Options Database Keys: 456 . -ts_adapt_max_ignore <max_ignore> - to set the threshold 457 458 Level: intermediate 459 460 .seealso: `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()` 461 @*/ 462 PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt,PetscReal max_ignore) 463 { 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 466 PetscValidLogicalCollectiveReal(adapt,max_ignore,2); 467 adapt->ignore_max = max_ignore; 468 PetscFunctionReturn(0); 469 } 470 471 /*@ 472 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). 473 474 Not Collective 475 476 Input Parameter: 477 . adapt - adaptive controller context 478 479 Output Parameter: 480 . max_ignore - threshold for solution components that are ignored during error estimation 481 482 Level: intermediate 483 484 .seealso: `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()` 485 @*/ 486 PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt,PetscReal *max_ignore) 487 { 488 PetscFunctionBegin; 489 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 490 PetscValidRealPointer(max_ignore,2); 491 *max_ignore = adapt->ignore_max; 492 PetscFunctionReturn(0); 493 } 494 495 /*@ 496 TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 497 498 Logically collective on TSAdapt 499 500 Input Parameters: 501 + adapt - adaptive controller context 502 . low - admissible decrease factor 503 - high - admissible increase factor 504 505 Options Database Keys: 506 . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 507 508 Level: intermediate 509 510 .seealso: `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()` 511 @*/ 512 PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 513 { 514 PetscFunctionBegin; 515 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 516 PetscValidLogicalCollectiveReal(adapt,low,2); 517 PetscValidLogicalCollectiveReal(adapt,high,3); 518 PetscCheck(low == PETSC_DEFAULT || low >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 519 PetscCheck(low == PETSC_DEFAULT || low <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 520 PetscCheck(high == PETSC_DEFAULT || high >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be greater than one",(double)high); 521 if (low != PETSC_DEFAULT) adapt->clip[0] = low; 522 if (high != PETSC_DEFAULT) adapt->clip[1] = high; 523 PetscFunctionReturn(0); 524 } 525 526 /*@ 527 TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 528 529 Not Collective 530 531 Input Parameter: 532 . adapt - adaptive controller context 533 534 Output Parameters: 535 + low - optional, admissible decrease factor 536 - high - optional, admissible increase factor 537 538 Level: intermediate 539 540 .seealso: `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()` 541 @*/ 542 PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 543 { 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 546 if (low) PetscValidRealPointer(low,2); 547 if (high) PetscValidRealPointer(high,3); 548 if (low) *low = adapt->clip[0]; 549 if (high) *high = adapt->clip[1]; 550 PetscFunctionReturn(0); 551 } 552 553 /*@ 554 TSAdaptSetScaleSolveFailed - Scale step by this factor if solve fails 555 556 Logically collective on TSAdapt 557 558 Input Parameters: 559 + adapt - adaptive controller context 560 - scale - scale 561 562 Options Database Keys: 563 . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails 564 565 Level: intermediate 566 567 .seealso: `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()` 568 @*/ 569 PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt,PetscReal scale) 570 { 571 PetscFunctionBegin; 572 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 573 PetscValidLogicalCollectiveReal(adapt,scale,2); 574 PetscCheck(scale == PETSC_DEFAULT || scale > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be positive",(double)scale); 575 PetscCheck(scale == PETSC_DEFAULT || scale <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scale factor %g must be less than one",(double)scale); 576 if (scale != PETSC_DEFAULT) adapt->scale_solve_failed = scale; 577 PetscFunctionReturn(0); 578 } 579 580 /*@ 581 TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size 582 583 Not Collective 584 585 Input Parameter: 586 . adapt - adaptive controller context 587 588 Output Parameter: 589 . scale - scale factor 590 591 Level: intermediate 592 593 .seealso: `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()` 594 @*/ 595 PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt,PetscReal *scale) 596 { 597 PetscFunctionBegin; 598 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 599 if (scale) PetscValidRealPointer(scale,2); 600 if (scale) *scale = adapt->scale_solve_failed; 601 PetscFunctionReturn(0); 602 } 603 604 /*@ 605 TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 606 607 Logically collective on TSAdapt 608 609 Input Parameters: 610 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 611 . hmin - minimum time step 612 - hmax - maximum time step 613 614 Options Database Keys: 615 + -ts_adapt_dt_min <min> - to set minimum time step 616 - -ts_adapt_dt_max <max> - to set maximum time step 617 618 Level: intermediate 619 620 .seealso: `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()` 621 @*/ 622 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 623 { 624 PetscFunctionBegin; 625 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 626 PetscValidLogicalCollectiveReal(adapt,hmin,2); 627 PetscValidLogicalCollectiveReal(adapt,hmax,3); 628 PetscCheck(hmin == PETSC_DEFAULT || hmin >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin); 629 PetscCheck(hmax == PETSC_DEFAULT || hmax >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax); 630 if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 631 if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 632 hmin = adapt->dt_min; 633 hmax = adapt->dt_max; 634 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); 635 PetscFunctionReturn(0); 636 } 637 638 /*@ 639 TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 640 641 Not Collective 642 643 Input Parameter: 644 . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 645 646 Output Parameters: 647 + hmin - minimum time step 648 - hmax - maximum time step 649 650 Level: intermediate 651 652 .seealso: `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()` 653 @*/ 654 PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 655 { 656 PetscFunctionBegin; 657 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 658 if (hmin) PetscValidRealPointer(hmin,2); 659 if (hmax) PetscValidRealPointer(hmax,3); 660 if (hmin) *hmin = adapt->dt_min; 661 if (hmax) *hmax = adapt->dt_max; 662 PetscFunctionReturn(0); 663 } 664 665 /* 666 TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 667 668 Collective on TSAdapt 669 670 Input Parameter: 671 . adapt - the TSAdapt context 672 673 Options Database Keys: 674 + -ts_adapt_type <type> - algorithm to use for adaptivity 675 . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 676 . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 677 . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 678 . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 679 . -ts_adapt_dt_min <min> - minimum timestep to use 680 . -ts_adapt_dt_max <max> - maximum timestep to use 681 . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 682 . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 683 - -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver 684 685 Level: advanced 686 687 Notes: 688 This function is automatically called by TSSetFromOptions() 689 690 .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`, 691 `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()` 692 */ 693 PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 694 { 695 char type[256] = TSADAPTBASIC; 696 PetscReal safety,reject_safety,clip[2],scale,hmin,hmax; 697 PetscBool set,flg; 698 PetscInt two; 699 700 PetscFunctionBegin; 701 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,2); 702 /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 703 * function can only be called from inside TSSetFromOptions() */ 704 PetscOptionsHeadBegin(PetscOptionsObject,"TS Adaptivity options"); 705 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)); 706 if (flg || !((PetscObject)adapt)->type_name) { 707 PetscCall(TSAdaptSetType(adapt,type)); 708 } 709 710 PetscCall(PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set)); 711 if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt,flg)); 712 713 safety = adapt->safety; 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; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 719 PetscCall(PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set)); 720 PetscCheck(!set || (two == 2),PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 721 if (set) PetscCall(TSAdaptSetClip(adapt,clip[0],clip[1])); 722 723 hmin = adapt->dt_min; hmax = adapt->dt_max; 724 PetscCall(PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set)); 725 PetscCall(PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg)); 726 if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt,hmin,hmax)); 727 728 PetscCall(PetscOptionsReal("-ts_adapt_max_ignore","Adaptor ignores (absolute) solution values smaller than this value","",adapt->ignore_max,&adapt->ignore_max,&set)); 729 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)); 730 731 PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","TSAdaptSetScaleSolveFailed",adapt->scale_solve_failed,&scale,&set)); 732 if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt,scale)); 733 734 PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL)); 735 PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY,PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 736 737 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)); 738 739 PetscCall(PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set)); 740 if (set) PetscCall(TSAdaptSetMonitor(adapt,flg)); 741 742 if (adapt->ops->setfromoptions) PetscCall((*adapt->ops->setfromoptions)(PetscOptionsObject,adapt)); 743 PetscOptionsHeadEnd(); 744 PetscFunctionReturn(0); 745 } 746 747 /*@ 748 TSAdaptCandidatesClear - clear any previously set candidate schemes 749 750 Logically collective on TSAdapt 751 752 Input Parameter: 753 . adapt - adaptive controller 754 755 Level: developer 756 757 .seealso: `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 758 @*/ 759 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 760 { 761 PetscFunctionBegin; 762 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 763 PetscCall(PetscMemzero(&adapt->candidates,sizeof(adapt->candidates))); 764 PetscFunctionReturn(0); 765 } 766 767 /*@C 768 TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 769 770 Logically collective on TSAdapt 771 772 Input Parameters: 773 + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 774 . name - name of the candidate scheme to add 775 . order - order of the candidate scheme 776 . stageorder - stage order of the candidate scheme 777 . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 778 . cost - relative measure of the amount of work required for the candidate scheme 779 - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 780 781 Note: 782 This routine is not available in Fortran. 783 784 Level: developer 785 786 .seealso: `TSAdaptCandidatesClear()`, `TSAdaptChoose()` 787 @*/ 788 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 789 { 790 PetscInt c; 791 792 PetscFunctionBegin; 793 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 794 PetscCheck(order >= 1,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %" PetscInt_FMT " must be a positive integer",order); 795 if (inuse) { 796 PetscCheck(!adapt->candidates.inuse_set,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 797 adapt->candidates.inuse_set = PETSC_TRUE; 798 } 799 /* first slot if this is the current scheme, otherwise the next available slot */ 800 c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 801 802 adapt->candidates.name[c] = name; 803 adapt->candidates.order[c] = order; 804 adapt->candidates.stageorder[c] = stageorder; 805 adapt->candidates.ccfl[c] = ccfl; 806 adapt->candidates.cost[c] = cost; 807 adapt->candidates.n++; 808 PetscFunctionReturn(0); 809 } 810 811 /*@C 812 TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 813 814 Not Collective 815 816 Input Parameter: 817 . adapt - time step adaptivity context 818 819 Output Parameters: 820 + n - number of candidate schemes, always at least 1 821 . order - the order of each candidate scheme 822 . stageorder - the stage order of each candidate scheme 823 . ccfl - the CFL coefficient of each scheme 824 - cost - the relative cost of each scheme 825 826 Level: developer 827 828 Note: 829 The current scheme is always returned in the first slot 830 831 .seealso: `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()` 832 @*/ 833 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 834 { 835 PetscFunctionBegin; 836 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 837 if (n) *n = adapt->candidates.n; 838 if (order) *order = adapt->candidates.order; 839 if (stageorder) *stageorder = adapt->candidates.stageorder; 840 if (ccfl) *ccfl = adapt->candidates.ccfl; 841 if (cost) *cost = adapt->candidates.cost; 842 PetscFunctionReturn(0); 843 } 844 845 /*@C 846 TSAdaptChoose - choose which method and step size to use for the next step 847 848 Collective on TSAdapt 849 850 Input Parameters: 851 + adapt - adaptive contoller 852 . ts - time stepper 853 - h - current step size 854 855 Output Parameters: 856 + next_sc - optional, scheme to use for the next step 857 . next_h - step size to use for the next step 858 - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 859 860 Note: 861 The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 862 being retried after an initial rejection. 863 864 Level: developer 865 866 .seealso: `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()` 867 @*/ 868 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 869 { 870 PetscInt ncandidates = adapt->candidates.n; 871 PetscInt scheme = 0; 872 PetscReal wlte = -1.0; 873 PetscReal wltea = -1.0; 874 PetscReal wlter = -1.0; 875 876 PetscFunctionBegin; 877 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 878 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 879 if (next_sc) PetscValidIntPointer(next_sc,4); 880 PetscValidRealPointer(next_h,5); 881 PetscValidBoolPointer(accept,6); 882 if (next_sc) *next_sc = 0; 883 884 /* Do not mess with adaptivity while handling events*/ 885 if (ts->event && ts->event->status != TSEVENT_NONE) { 886 *next_h = h; 887 *accept = PETSC_TRUE; 888 PetscFunctionReturn(0); 889 } 890 891 PetscCall((*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter)); 892 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); 893 PetscCheck(*next_h >= 0,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 894 if (next_sc) *next_sc = scheme; 895 896 if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 897 /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 898 PetscReal t = ts->ptime + ts->time_step, h = *next_h; 899 PetscReal tend = t + h, tmax, hmax; 900 PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 901 PetscReal b = adapt->matchstepfac[1]; 902 903 if (ts->tspan) { 904 if (PetscIsCloseAtTol(t,ts->tspan->span_times[ts->tspan->spanctr],10*PETSC_MACHINE_EPSILON,0)) /* hit a span time point */ 905 if (ts->tspan->spanctr+1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr+1]; 906 else tmax = ts->max_time; /* hit the last span time point */ 907 else tmax = ts->tspan->span_times[ts->tspan->spanctr]; 908 } else tmax = ts->max_time; 909 hmax = tmax - t; 910 if (t < tmax && tend > tmax) *next_h = hmax; 911 if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 912 if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 913 } 914 915 if (adapt->monitor) { 916 const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 917 PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 918 if (wlte < 0) { 919 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",(double)ts->ptime,(double)h,(double)*next_h)); 920 } else { 921 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",(double)ts->ptime,(double)h,(double)*next_h,(double)wlte,(double)wltea,(double)wlter)); 922 } 923 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 924 } 925 PetscFunctionReturn(0); 926 } 927 928 /*@ 929 TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 930 before increasing the time step. 931 932 Logicially Collective on TSAdapt 933 934 Input Parameters: 935 + adapt - adaptive controller context 936 - cnt - the number of timesteps 937 938 Options Database Key: 939 . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 940 941 Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 942 The successful use of this option is problem dependent 943 944 Developer Note: there is no theory to support this option 945 946 Level: advanced 947 948 .seealso: 949 @*/ 950 PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 951 { 952 PetscFunctionBegin; 953 adapt->timestepjustdecreased_delay = cnt; 954 PetscFunctionReturn(0); 955 } 956 957 /*@ 958 TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails or solution vector is infeasible) 959 960 Collective on TSAdapt 961 962 Input Parameters: 963 + adapt - adaptive controller context 964 . ts - time stepper 965 . t - Current simulation time 966 - Y - Current solution vector 967 968 Output Parameter: 969 . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 970 971 Level: developer 972 973 .seealso: 974 @*/ 975 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 976 { 977 SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 978 979 PetscFunctionBegin; 980 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 981 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 982 PetscValidBoolPointer(accept,5); 983 984 if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes,&snesreason)); 985 if (snesreason < 0) { 986 *accept = PETSC_FALSE; 987 if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 988 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 989 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)); 990 if (adapt->monitor) { 991 PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 992 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,(double)ts->ptime,(double)ts->time_step,ts->num_snes_failures)); 993 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 994 } 995 } 996 } else { 997 *accept = PETSC_TRUE; 998 PetscCall(TSFunctionDomainError(ts,t,Y,accept)); 999 if (*accept && adapt->checkstage) { 1000 PetscCall((*adapt->checkstage)(adapt,ts,t,Y,accept)); 1001 if (!*accept) { 1002 PetscCall(PetscInfo(ts,"Step=%" PetscInt_FMT ", solution rejected by user function provided by TSSetFunctionDomainError()\n",ts->steps)); 1003 if (adapt->monitor) { 1004 PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 1005 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)); 1006 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 1007 } 1008 } 1009 } 1010 } 1011 1012 if (!(*accept) && !ts->reason) { 1013 PetscReal dt,new_dt; 1014 PetscCall(TSGetTimeStep(ts,&dt)); 1015 new_dt = dt * adapt->scale_solve_failed; 1016 PetscCall(TSSetTimeStep(ts,new_dt)); 1017 adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 1018 if (adapt->monitor) { 1019 PetscCall(PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 1020 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)); 1021 PetscCall(PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel)); 1022 } 1023 } 1024 PetscFunctionReturn(0); 1025 } 1026 1027 /*@ 1028 TSAdaptCreate - create an adaptive controller context for time stepping 1029 1030 Collective 1031 1032 Input Parameter: 1033 . comm - The communicator 1034 1035 Output Parameter: 1036 . adapt - new TSAdapt object 1037 1038 Level: developer 1039 1040 Notes: 1041 TSAdapt creation is handled by TS, so users should not need to call this function. 1042 1043 .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()` 1044 @*/ 1045 PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 1046 { 1047 TSAdapt adapt; 1048 1049 PetscFunctionBegin; 1050 PetscValidPointer(inadapt,2); 1051 *inadapt = NULL; 1052 PetscCall(TSAdaptInitializePackage()); 1053 1054 PetscCall(PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView)); 1055 1056 adapt->always_accept = PETSC_FALSE; 1057 adapt->safety = 0.9; 1058 adapt->reject_safety = 0.5; 1059 adapt->clip[0] = 0.1; 1060 adapt->clip[1] = 10.; 1061 adapt->dt_min = 1e-20; 1062 adapt->dt_max = 1e+20; 1063 adapt->ignore_max = -1.0; 1064 adapt->glee_use_local = PETSC_TRUE; 1065 adapt->scale_solve_failed = 0.25; 1066 /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 1067 to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 1068 adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 1069 adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 1070 adapt->wnormtype = NORM_2; 1071 adapt->timestepjustdecreased_delay = 0; 1072 1073 *inadapt = adapt; 1074 PetscFunctionReturn(0); 1075 } 1076