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