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