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