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