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