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_GLEE(TSAdapt); 11 PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt); 12 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 13 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(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(TSADAPTGLEE ,TSAdaptCreate_GLEE);CHKERRQ(ierr); 71 ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 72 ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr); 73 ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);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 TSCreate_GLLE() 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 ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 146 if (match) PetscFunctionReturn(0); 147 ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 148 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 149 if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 150 ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 151 ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 152 ierr = (*r)(adapt);CHKERRQ(ierr); 153 PetscFunctionReturn(0); 154 } 155 156 /*@C 157 TSAdaptGetType - gets the TS adapter method type (as a string). 158 159 Not Collective 160 161 Input Parameter: 162 . adapt - The TS adapter, most likely obtained with TSGetAdapt() 163 164 Output Parameter: 165 . type - The name of TS adapter method 166 167 Level: intermediate 168 169 .keywords: TSAdapt, get, type 170 .seealso TSAdaptSetType() 171 @*/ 172 PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 173 { 174 PetscFunctionBegin; 175 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 176 PetscValidPointer(type,2); 177 *type = ((PetscObject)adapt)->type_name; 178 PetscFunctionReturn(0); 179 } 180 181 PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 182 { 183 PetscErrorCode ierr; 184 185 PetscFunctionBegin; 186 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 187 ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 188 PetscFunctionReturn(0); 189 } 190 191 /*@C 192 TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 193 194 Collective on PetscViewer 195 196 Input Parameters: 197 + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 198 some related function before a call to TSAdaptLoad(). 199 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 200 HDF5 file viewer, obtained from PetscViewerHDF5Open() 201 202 Level: intermediate 203 204 Notes: 205 The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 206 207 Notes for advanced users: 208 Most users should not need to know the details of the binary storage 209 format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 210 But for anyone who's interested, the standard binary matrix storage 211 format is 212 .vb 213 has not yet been determined 214 .ve 215 216 .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 217 @*/ 218 PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 219 { 220 PetscErrorCode ierr; 221 PetscBool isbinary; 222 char type[256]; 223 224 PetscFunctionBegin; 225 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 226 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 227 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 228 if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 229 230 ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 231 ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 232 if (adapt->ops->load) { 233 ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 234 } 235 PetscFunctionReturn(0); 236 } 237 238 PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 239 { 240 PetscErrorCode ierr; 241 PetscBool iascii,isbinary; 242 243 PetscFunctionBegin; 244 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 245 if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 246 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 247 PetscCheckSameComm(adapt,1,viewer,2); 248 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 249 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 250 if (iascii) { 251 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 252 ierr = PetscViewerASCIIPrintf(viewer," number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 253 if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 254 ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 255 ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 256 if (adapt->ops->view) { 257 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 258 ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 259 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 260 } 261 } else if (isbinary) { 262 char type[256]; 263 264 /* need to save FILE_CLASS_ID for adapt class */ 265 ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 266 ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 267 } else if (adapt->ops->view) { 268 ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 269 } 270 PetscFunctionReturn(0); 271 } 272 273 /*@ 274 TSAdaptReset - Resets a TSAdapt context. 275 276 Collective on TS 277 278 Input Parameter: 279 . adapt - the TSAdapt context obtained from TSAdaptCreate() 280 281 Level: developer 282 283 .seealso: TSAdaptCreate(), TSAdaptDestroy() 284 @*/ 285 PetscErrorCode TSAdaptReset(TSAdapt adapt) 286 { 287 PetscErrorCode ierr; 288 289 PetscFunctionBegin; 290 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 291 if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 292 PetscFunctionReturn(0); 293 } 294 295 PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 296 { 297 PetscErrorCode ierr; 298 299 PetscFunctionBegin; 300 if (!*adapt) PetscFunctionReturn(0); 301 PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 302 if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 303 304 ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 305 306 if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 307 ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 308 ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 309 PetscFunctionReturn(0); 310 } 311 312 /*@ 313 TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 314 315 Collective on TSAdapt 316 317 Input Arguments: 318 + adapt - adaptive controller context 319 - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 320 321 Options Database Keys: 322 . -ts_adapt_monitor 323 324 Level: intermediate 325 326 .seealso: TSAdaptChoose() 327 @*/ 328 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 329 { 330 PetscErrorCode ierr; 331 332 PetscFunctionBegin; 333 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 334 PetscValidLogicalCollectiveBool(adapt,flg,2); 335 if (flg) { 336 if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 337 } else { 338 ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 339 } 340 PetscFunctionReturn(0); 341 } 342 343 /*@C 344 TSAdaptSetCheckStage - Set a callback to check convergence for a stage 345 346 Logically collective on TSAdapt 347 348 Input Arguments: 349 + adapt - adaptive controller context 350 - func - stage check function 351 352 Arguments of func: 353 $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 354 355 + adapt - adaptive controller context 356 . ts - time stepping context 357 - accept - pending choice of whether to accept, can be modified by this routine 358 359 Level: advanced 360 361 .seealso: TSAdaptChoose() 362 @*/ 363 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 364 { 365 366 PetscFunctionBegin; 367 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 368 adapt->checkstage = func; 369 PetscFunctionReturn(0); 370 } 371 372 /*@ 373 TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 374 any error or stability condition not meeting the prescribed goal. 375 376 377 Input Arguments: 378 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 379 - flag - whether to always accept steps 380 381 Options Database Keys: 382 . -ts_adapt_always_accept 383 384 Level: intermediate 385 386 .seealso: TSAdapt, TSAdaptChoose() 387 @*/ 388 PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 389 { 390 PetscFunctionBegin; 391 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 392 PetscValidLogicalCollectiveBool(adapt,flag,2); 393 adapt->always_accept = flag; 394 PetscFunctionReturn(0); 395 } 396 397 /*@ 398 TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 399 400 Logically Collective 401 402 Input Arguments: 403 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 404 . hmin - minimum time step 405 - hmax - maximum time step 406 407 Options Database Keys: 408 + -ts_adapt_dt_min - minimum time step 409 - -ts_adapt_dt_max - maximum time step 410 411 Level: intermediate 412 413 .seealso: TSAdapt, TSAdaptChoose() 414 @*/ 415 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 416 { 417 418 PetscFunctionBegin; 419 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 420 PetscValidLogicalCollectiveReal(adapt,hmin,2); 421 PetscValidLogicalCollectiveReal(adapt,hmax,3); 422 if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin); 423 if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax); 424 if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 425 if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 426 #if defined(PETSC_USE_DEBUG) 427 hmin = adapt->dt_min; 428 hmax = adapt->dt_max; 429 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); 430 #endif 431 PetscFunctionReturn(0); 432 } 433 434 /* 435 TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 436 437 Collective on TSAdapt 438 439 Input Parameter: 440 . adapt - the TSAdapt context 441 442 Options Database Keys: 443 + -ts_adapt_type <type> - basic 444 . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 445 . -ts_adapt_dt_min <min> - minimum timestep to use 446 . -ts_adapt_dt_max <max> - maximum timestep to use 447 . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 448 - -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 449 450 Level: advanced 451 452 Notes: 453 This function is automatically called by TSSetFromOptions() 454 455 .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 456 457 .seealso: TSGetType() 458 */ 459 PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 460 { 461 PetscErrorCode ierr; 462 char type[256] = TSADAPTBASIC; 463 PetscReal hmin,hmax; 464 PetscBool set,flg; 465 466 PetscFunctionBegin; 467 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 468 /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 469 * function can only be called from inside TSSetFromOptions() */ 470 ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 471 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); 472 if (flg || !((PetscObject)adapt)->type_name) { 473 ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 474 } 475 ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 476 if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 477 ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin=adapt->dt_min,&hmin,&set);CHKERRQ(ierr); 478 ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax=adapt->dt_max,&hmax,&flg);CHKERRQ(ierr); 479 if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 480 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); 481 ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 482 if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 483 ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 484 if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 485 if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 486 ierr = PetscOptionsTail();CHKERRQ(ierr); 487 PetscFunctionReturn(0); 488 } 489 490 /*@ 491 TSAdaptCandidatesClear - clear any previously set candidate schemes 492 493 Logically Collective 494 495 Input Argument: 496 . adapt - adaptive controller 497 498 Level: developer 499 500 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 501 @*/ 502 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 503 { 504 PetscErrorCode ierr; 505 506 PetscFunctionBegin; 507 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 508 ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 509 PetscFunctionReturn(0); 510 } 511 512 /*@C 513 TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 514 515 Logically Collective 516 517 Input Arguments: 518 + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 519 . name - name of the candidate scheme to add 520 . order - order of the candidate scheme 521 . stageorder - stage order of the candidate scheme 522 . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 523 . cost - relative measure of the amount of work required for the candidate scheme 524 - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 525 526 Note: 527 This routine is not available in Fortran. 528 529 Level: developer 530 531 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 532 @*/ 533 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 534 { 535 PetscInt c; 536 537 PetscFunctionBegin; 538 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 539 if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 540 if (inuse) { 541 if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 542 adapt->candidates.inuse_set = PETSC_TRUE; 543 } 544 /* first slot if this is the current scheme, otherwise the next available slot */ 545 c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 546 547 adapt->candidates.name[c] = name; 548 adapt->candidates.order[c] = order; 549 adapt->candidates.stageorder[c] = stageorder; 550 adapt->candidates.ccfl[c] = ccfl; 551 adapt->candidates.cost[c] = cost; 552 adapt->candidates.n++; 553 PetscFunctionReturn(0); 554 } 555 556 /*@C 557 TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 558 559 Not Collective 560 561 Input Arguments: 562 . adapt - time step adaptivity context 563 564 Output Arguments: 565 + n - number of candidate schemes, always at least 1 566 . order - the order of each candidate scheme 567 . stageorder - the stage order of each candidate scheme 568 . ccfl - the CFL coefficient of each scheme 569 - cost - the relative cost of each scheme 570 571 Level: developer 572 573 Note: 574 The current scheme is always returned in the first slot 575 576 .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 577 @*/ 578 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 579 { 580 PetscFunctionBegin; 581 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 582 if (n) *n = adapt->candidates.n; 583 if (order) *order = adapt->candidates.order; 584 if (stageorder) *stageorder = adapt->candidates.stageorder; 585 if (ccfl) *ccfl = adapt->candidates.ccfl; 586 if (cost) *cost = adapt->candidates.cost; 587 PetscFunctionReturn(0); 588 } 589 590 /*@C 591 TSAdaptChoose - choose which method and step size to use for the next step 592 593 Logically Collective 594 595 Input Arguments: 596 + adapt - adaptive contoller 597 - h - current step size 598 599 Output Arguments: 600 + next_sc - optional, scheme to use for the next step 601 . next_h - step size to use for the next step 602 - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 603 604 Note: 605 The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 606 being retried after an initial rejection. 607 608 Level: developer 609 610 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 611 @*/ 612 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 613 { 614 PetscErrorCode ierr; 615 PetscInt ncandidates = adapt->candidates.n; 616 PetscInt scheme = 0; 617 PetscReal wlte = -1.0; 618 PetscReal wltea = -1.0; 619 PetscReal wlter = -1.0; 620 621 PetscFunctionBegin; 622 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 623 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 624 if (next_sc) PetscValidIntPointer(next_sc,4); 625 PetscValidPointer(next_h,5); 626 PetscValidIntPointer(accept,6); 627 if (next_sc) *next_sc = 0; 628 629 /* Do not mess with adaptivity while handling events*/ 630 if (ts->event && ts->event->status != TSEVENT_NONE) { 631 *next_h = h; 632 *accept = PETSC_TRUE; 633 PetscFunctionReturn(0); 634 } 635 636 ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 637 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); 638 if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 639 if (next_sc) *next_sc = scheme; 640 641 if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 642 /* Reduce time step if it overshoots max time */ 643 if (ts->ptime + ts->time_step + *next_h >= ts->max_time) { 644 PetscReal next_dt = ts->max_time - (ts->ptime + ts->time_step); 645 if (next_dt > PETSC_SMALL) *next_h = next_dt; 646 else ts->reason = TS_CONVERGED_TIME; 647 } 648 } 649 650 if (adapt->monitor) { 651 const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 652 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 653 if (wlte < 0) { 654 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); 655 } else { 656 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); 657 } 658 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 659 } 660 PetscFunctionReturn(0); 661 } 662 663 /*@ 664 TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 665 666 Collective 667 668 Input Arguments: 669 + adapt - adaptive controller context 670 . ts - time stepper 671 . t - Current simulation time 672 - Y - Current solution vector 673 674 Output Arguments: 675 . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 676 677 Level: developer 678 679 .seealso: 680 @*/ 681 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 682 { 683 PetscErrorCode ierr; 684 SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 685 686 PetscFunctionBegin; 687 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 688 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 689 PetscValidIntPointer(accept,3); 690 691 if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 692 if (snesreason < 0) { 693 *accept = PETSC_FALSE; 694 if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 695 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 696 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); 697 if (adapt->monitor) { 698 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 699 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); 700 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 701 } 702 } 703 } else { 704 *accept = PETSC_TRUE; 705 ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 706 if(*accept && adapt->checkstage) { 707 ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 708 } 709 } 710 711 if(!(*accept) && !ts->reason) { 712 PetscReal dt,new_dt; 713 ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 714 new_dt = dt * adapt->scale_solve_failed; 715 ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 716 if (adapt->monitor) { 717 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 718 ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected t=%-11g+%10.3e retrying with dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)dt,(double)new_dt);CHKERRQ(ierr); 719 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 720 } 721 } 722 PetscFunctionReturn(0); 723 } 724 725 /*@ 726 TSAdaptCreate - create an adaptive controller context for time stepping 727 728 Collective on MPI_Comm 729 730 Input Parameter: 731 . comm - The communicator 732 733 Output Parameter: 734 . adapt - new TSAdapt object 735 736 Level: developer 737 738 Notes: 739 TSAdapt creation is handled by TS, so users should not need to call this function. 740 741 .keywords: TSAdapt, create 742 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 743 @*/ 744 PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 745 { 746 PetscErrorCode ierr; 747 TSAdapt adapt; 748 749 PetscFunctionBegin; 750 PetscValidPointer(inadapt,1); 751 *inadapt = NULL; 752 ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 753 754 ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 755 756 adapt->always_accept = PETSC_FALSE; 757 adapt->dt_min = 1e-20; 758 #if defined(PETSC_USE_REAL_SINGLE) 759 adapt->dt_max = 1e35; 760 #else 761 adapt->dt_max = 1e50; 762 #endif 763 adapt->scale_solve_failed = 0.25; 764 adapt->wnormtype = NORM_2; 765 ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr); 766 767 *inadapt = adapt; 768 PetscFunctionReturn(0); 769 } 770