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