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