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