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