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_DSP(TSAdapt); 13 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 14 PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt); 15 PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt); 16 17 /*@C 18 TSAdaptRegister - adds a TSAdapt implementation 19 20 Not Collective 21 22 Input Parameters: 23 + name_scheme - name of user-defined adaptivity scheme 24 - routine_create - routine to create method context 25 26 Notes: 27 TSAdaptRegister() may be called multiple times to add several user-defined families. 28 29 Sample usage: 30 .vb 31 TSAdaptRegister("my_scheme",MySchemeCreate); 32 .ve 33 34 Then, your scheme can be chosen with the procedural interface via 35 $ TSAdaptSetType(ts,"my_scheme") 36 or at runtime via the option 37 $ -ts_adapt_type my_scheme 38 39 Level: advanced 40 41 .keywords: TSAdapt, register 42 43 .seealso: TSAdaptRegisterAll() 44 @*/ 45 PetscErrorCode TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt)) 46 { 47 PetscErrorCode ierr; 48 49 PetscFunctionBegin; 50 ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 51 ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr); 52 PetscFunctionReturn(0); 53 } 54 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(TSADAPTDSP, TSAdaptCreate_DSP);CHKERRQ(ierr); 76 ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 77 ierr = TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE);CHKERRQ(ierr); 78 ierr = TSAdaptRegister(TSADAPTHISTORY,TSAdaptCreate_History);CHKERRQ(ierr); 79 PetscFunctionReturn(0); 80 } 81 82 /*@C 83 TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 84 called from PetscFinalize(). 85 86 Level: developer 87 88 .keywords: Petsc, destroy, package 89 .seealso: PetscFinalize() 90 @*/ 91 PetscErrorCode TSAdaptFinalizePackage(void) 92 { 93 PetscErrorCode ierr; 94 95 PetscFunctionBegin; 96 ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 97 TSAdaptPackageInitialized = PETSC_FALSE; 98 TSAdaptRegisterAllCalled = PETSC_FALSE; 99 PetscFunctionReturn(0); 100 } 101 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 TSAdaptCreate() 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 /*@C 126 TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 127 128 Logicially Collective on TSAdapt 129 130 Input Parameter: 131 + adapt - the TS adapter, most likely obtained with TSGetAdapt() 132 - type - either TSADAPTBASIC or TSADAPTNONE 133 134 Options Database: 135 . -ts_adapt_type <basic or dsp or none> - to set the adapter type 136 137 Level: intermediate 138 139 .keywords: TSAdapt, create 140 141 .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType() 142 @*/ 143 PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 144 { 145 PetscBool match; 146 PetscErrorCode ierr,(*r)(TSAdapt); 147 148 PetscFunctionBegin; 149 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 150 PetscValidCharPointer(type,2); 151 ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 152 if (match) PetscFunctionReturn(0); 153 ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 154 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 155 if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 156 ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 157 ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 158 ierr = (*r)(adapt);CHKERRQ(ierr); 159 PetscFunctionReturn(0); 160 } 161 162 /*@C 163 TSAdaptGetType - gets the TS adapter method type (as a string). 164 165 Not Collective 166 167 Input Parameter: 168 . adapt - The TS adapter, most likely obtained with TSGetAdapt() 169 170 Output Parameter: 171 . type - The name of TS adapter method 172 173 Level: intermediate 174 175 .keywords: TSAdapt, get, type 176 .seealso TSAdaptSetType() 177 @*/ 178 PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type) 179 { 180 PetscFunctionBegin; 181 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 182 PetscValidPointer(type,2); 183 *type = ((PetscObject)adapt)->type_name; 184 PetscFunctionReturn(0); 185 } 186 187 PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 188 { 189 PetscErrorCode ierr; 190 191 PetscFunctionBegin; 192 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 193 ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 /*@C 198 TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 199 200 Collective on PetscViewer 201 202 Input Parameters: 203 + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 204 some related function before a call to TSAdaptLoad(). 205 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 206 HDF5 file viewer, obtained from PetscViewerHDF5Open() 207 208 Level: intermediate 209 210 Notes: 211 The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 212 213 Notes for advanced users: 214 Most users should not need to know the details of the binary storage 215 format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 216 But for anyone who's interested, the standard binary matrix storage 217 format is 218 .vb 219 has not yet been determined 220 .ve 221 222 .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 223 @*/ 224 PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 225 { 226 PetscErrorCode ierr; 227 PetscBool isbinary; 228 char type[256]; 229 230 PetscFunctionBegin; 231 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 232 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 233 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 234 if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 235 236 ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 237 ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 238 if (adapt->ops->load) { 239 ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 240 } 241 PetscFunctionReturn(0); 242 } 243 244 PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 245 { 246 PetscErrorCode ierr; 247 PetscBool iascii,isbinary,isnone; 248 249 PetscFunctionBegin; 250 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 251 if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 252 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 253 PetscCheckSameComm(adapt,1,viewer,2); 254 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 255 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 256 if (iascii) { 257 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 258 ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);CHKERRQ(ierr); 259 if (!isnone) { 260 if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 261 ierr = PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety);CHKERRQ(ierr); 262 ierr = PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety);CHKERRQ(ierr); 263 ierr = PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1]);CHKERRQ(ierr); 264 ierr = PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0]);CHKERRQ(ierr); 265 ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 266 ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 267 } 268 if (adapt->ops->view) { 269 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 270 ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 271 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 272 } 273 } else if (isbinary) { 274 char type[256]; 275 276 /* need to save FILE_CLASS_ID for adapt class */ 277 ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 278 ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 279 } else if (adapt->ops->view) { 280 ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 281 } 282 PetscFunctionReturn(0); 283 } 284 285 /*@ 286 TSAdaptReset - Resets a TSAdapt context. 287 288 Collective on TS 289 290 Input Parameter: 291 . adapt - the TSAdapt context obtained from TSAdaptCreate() 292 293 Level: developer 294 295 .seealso: TSAdaptCreate(), TSAdaptDestroy() 296 @*/ 297 PetscErrorCode TSAdaptReset(TSAdapt adapt) 298 { 299 PetscErrorCode ierr; 300 301 PetscFunctionBegin; 302 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 303 if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 304 PetscFunctionReturn(0); 305 } 306 307 PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 308 { 309 PetscErrorCode ierr; 310 311 PetscFunctionBegin; 312 if (!*adapt) PetscFunctionReturn(0); 313 PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 314 if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 315 316 ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 317 318 if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 319 ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 320 ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 321 PetscFunctionReturn(0); 322 } 323 324 /*@ 325 TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 326 327 Collective on TSAdapt 328 329 Input Arguments: 330 + adapt - adaptive controller context 331 - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 332 333 Options Database Keys: 334 . -ts_adapt_monitor - to turn on monitoring 335 336 Level: intermediate 337 338 .seealso: TSAdaptChoose() 339 @*/ 340 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 341 { 342 PetscErrorCode ierr; 343 344 PetscFunctionBegin; 345 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 346 PetscValidLogicalCollectiveBool(adapt,flg,2); 347 if (flg) { 348 if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 349 } else { 350 ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 351 } 352 PetscFunctionReturn(0); 353 } 354 355 /*@C 356 TSAdaptSetCheckStage - Set a callback to check convergence for a stage 357 358 Logically collective on TSAdapt 359 360 Input Arguments: 361 + adapt - adaptive controller context 362 - func - stage check function 363 364 Arguments of func: 365 $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 366 367 + adapt - adaptive controller context 368 . ts - time stepping context 369 - accept - pending choice of whether to accept, can be modified by this routine 370 371 Level: advanced 372 373 .seealso: TSAdaptChoose() 374 @*/ 375 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 376 { 377 378 PetscFunctionBegin; 379 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 380 adapt->checkstage = func; 381 PetscFunctionReturn(0); 382 } 383 384 /*@ 385 TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 386 any error or stability condition not meeting the prescribed goal. 387 388 Logically collective on TSAdapt 389 390 Input Arguments: 391 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 392 - flag - whether to always accept steps 393 394 Options Database Keys: 395 . -ts_adapt_always_accept - to always accept steps 396 397 Level: intermediate 398 399 .seealso: TSAdapt, TSAdaptChoose() 400 @*/ 401 PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 402 { 403 PetscFunctionBegin; 404 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 405 PetscValidLogicalCollectiveBool(adapt,flag,2); 406 adapt->always_accept = flag; 407 PetscFunctionReturn(0); 408 } 409 410 /*@ 411 TSAdaptSetSafety - Set safety factors 412 413 Logically collective on TSAdapt 414 415 Input Arguments: 416 + adapt - adaptive controller context 417 . safety - safety factor relative to target error/stability goal 418 - reject_safety - extra safety factor to apply if the last step was rejected 419 420 Options Database Keys: 421 + -ts_adapt_safety <safety> - to set safety factor 422 - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor 423 424 Level: intermediate 425 426 .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose() 427 @*/ 428 PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety) 429 { 430 PetscFunctionBegin; 431 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 432 PetscValidLogicalCollectiveReal(adapt,safety,2); 433 PetscValidLogicalCollectiveReal(adapt,reject_safety,3); 434 if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety); 435 if (safety != PETSC_DEFAULT && safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be less than one",(double)safety); 436 if (reject_safety != PETSC_DEFAULT && reject_safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be non negative",(double)reject_safety); 437 if (reject_safety != PETSC_DEFAULT && reject_safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be less than one",(double)reject_safety); 438 if (safety != PETSC_DEFAULT) adapt->safety = safety; 439 if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety; 440 PetscFunctionReturn(0); 441 } 442 443 /*@ 444 TSAdaptGetSafety - Get safety factors 445 446 Not Collective 447 448 Input Arguments: 449 . adapt - adaptive controller context 450 451 Ouput Arguments: 452 . safety - safety factor relative to target error/stability goal 453 + reject_safety - extra safety factor to apply if the last step was rejected 454 455 Level: intermediate 456 457 .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose() 458 @*/ 459 PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety) 460 { 461 PetscFunctionBegin; 462 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 463 if (safety) PetscValidRealPointer(safety,2); 464 if (reject_safety) PetscValidRealPointer(reject_safety,3); 465 if (safety) *safety = adapt->safety; 466 if (reject_safety) *reject_safety = adapt->reject_safety; 467 PetscFunctionReturn(0); 468 } 469 470 /*@ 471 TSAdaptSetClip - Sets the admissible decrease/increase factor in step size 472 473 Logically collective on TSAdapt 474 475 Input Arguments: 476 + adapt - adaptive controller context 477 . low - admissible decrease factor 478 - high - admissible increase factor 479 480 Options Database Keys: 481 . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors 482 483 Level: intermediate 484 485 .seealso: TSAdaptChoose(), TSAdaptGetClip() 486 @*/ 487 PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 488 { 489 PetscFunctionBegin; 490 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 491 PetscValidLogicalCollectiveReal(adapt,low,2); 492 PetscValidLogicalCollectiveReal(adapt,high,3); 493 if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 494 if (low != PETSC_DEFAULT && low > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 495 if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high); 496 if (low != PETSC_DEFAULT) adapt->clip[0] = low; 497 if (high != PETSC_DEFAULT) adapt->clip[1] = high; 498 PetscFunctionReturn(0); 499 } 500 501 /*@ 502 TSAdaptGetClip - Gets the admissible decrease/increase factor in step size 503 504 Not Collective 505 506 Input Arguments: 507 . adapt - adaptive controller context 508 509 Ouput Arguments: 510 + low - optional, admissible decrease factor 511 - high - optional, admissible increase factor 512 513 Level: intermediate 514 515 .seealso: TSAdaptChoose(), TSAdaptSetClip() 516 @*/ 517 PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 518 { 519 PetscFunctionBegin; 520 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 521 if (low) PetscValidRealPointer(low,2); 522 if (high) PetscValidRealPointer(high,3); 523 if (low) *low = adapt->clip[0]; 524 if (high) *high = adapt->clip[1]; 525 PetscFunctionReturn(0); 526 } 527 528 /*@ 529 TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller 530 531 Logically collective on TSAdapt 532 533 Input Arguments: 534 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 535 . hmin - minimum time step 536 - hmax - maximum time step 537 538 Options Database Keys: 539 + -ts_adapt_dt_min <min> - to set minimum time step 540 - -ts_adapt_dt_max <max> - to set maximum time step 541 542 Level: intermediate 543 544 .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose() 545 @*/ 546 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 547 { 548 549 PetscFunctionBegin; 550 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 551 PetscValidLogicalCollectiveReal(adapt,hmin,2); 552 PetscValidLogicalCollectiveReal(adapt,hmax,3); 553 if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin); 554 if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax); 555 if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 556 if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 557 hmin = adapt->dt_min; 558 hmax = adapt->dt_max; 559 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); 560 PetscFunctionReturn(0); 561 } 562 563 /*@ 564 TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller 565 566 Not Collective 567 568 Input Arguments: 569 . adapt - time step adaptivity context, usually gotten with TSGetAdapt() 570 571 Output Arguments: 572 + hmin - minimum time step 573 - hmax - maximum time step 574 575 Level: intermediate 576 577 .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose() 578 @*/ 579 PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax) 580 { 581 582 PetscFunctionBegin; 583 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 584 if (hmin) PetscValidRealPointer(hmin,2); 585 if (hmax) PetscValidRealPointer(hmax,3); 586 if (hmin) *hmin = adapt->dt_min; 587 if (hmax) *hmax = adapt->dt_max; 588 PetscFunctionReturn(0); 589 } 590 591 /* 592 TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 593 594 Collective on TSAdapt 595 596 Input Parameter: 597 . adapt - the TSAdapt context 598 599 Options Database Keys: 600 + -ts_adapt_type <type> - algorithm to use for adaptivity 601 . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 602 . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal 603 . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected 604 . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors 605 . -ts_adapt_dt_min <min> - minimum timestep to use 606 . -ts_adapt_dt_max <max> - maximum timestep to use 607 . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 608 . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 609 - -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver 610 611 Level: advanced 612 613 Notes: 614 This function is automatically called by TSSetFromOptions() 615 616 .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 617 618 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(), 619 TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor() 620 */ 621 PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 622 { 623 PetscErrorCode ierr; 624 char type[256] = TSADAPTBASIC; 625 PetscReal safety,reject_safety,clip[2],hmin,hmax; 626 PetscBool set,flg; 627 PetscInt two; 628 629 PetscFunctionBegin; 630 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 631 /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 632 * function can only be called from inside TSSetFromOptions() */ 633 ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 634 ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 635 if (flg || !((PetscObject)adapt)->type_name) { 636 ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 637 } 638 639 ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 640 if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 641 642 safety = adapt->safety; reject_safety = adapt->reject_safety; 643 ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr); 644 ierr = PetscOptionsReal("-ts_adapt_reject_safety","Extra safety factor to apply if the last step was rejected","TSAdaptSetSafety",reject_safety,&reject_safety,&flg);CHKERRQ(ierr); 645 if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);} 646 647 two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1]; 648 ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr); 649 if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip"); 650 if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);} 651 652 hmin = adapt->dt_min; hmax = adapt->dt_max; 653 ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr); 654 ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr); 655 if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 656 657 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); 658 659 ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 660 if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 661 662 ierr = PetscOptionsInt("-ts_adapt_time_step_increase_delay","Number of timesteps to delay increasing the time step after it has been decreased due to failed solver","TSAdaptSetTimeStepIncreaseDelay",adapt->timestepjustdecreased_delay,&adapt->timestepjustdecreased_delay,NULL);CHKERRQ(ierr); 663 664 ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 665 if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 666 667 if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 668 ierr = PetscOptionsTail();CHKERRQ(ierr); 669 PetscFunctionReturn(0); 670 } 671 672 /*@ 673 TSAdaptCandidatesClear - clear any previously set candidate schemes 674 675 Logically collective on TSAdapt 676 677 Input Argument: 678 . adapt - adaptive controller 679 680 Level: developer 681 682 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 683 @*/ 684 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 685 { 686 PetscErrorCode ierr; 687 688 PetscFunctionBegin; 689 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 690 ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 691 PetscFunctionReturn(0); 692 } 693 694 /*@C 695 TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 696 697 Logically collective on TSAdapt 698 699 Input Arguments: 700 + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 701 . name - name of the candidate scheme to add 702 . order - order of the candidate scheme 703 . stageorder - stage order of the candidate scheme 704 . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 705 . cost - relative measure of the amount of work required for the candidate scheme 706 - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 707 708 Note: 709 This routine is not available in Fortran. 710 711 Level: developer 712 713 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 714 @*/ 715 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 716 { 717 PetscInt c; 718 719 PetscFunctionBegin; 720 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 721 if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 722 if (inuse) { 723 if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 724 adapt->candidates.inuse_set = PETSC_TRUE; 725 } 726 /* first slot if this is the current scheme, otherwise the next available slot */ 727 c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 728 729 adapt->candidates.name[c] = name; 730 adapt->candidates.order[c] = order; 731 adapt->candidates.stageorder[c] = stageorder; 732 adapt->candidates.ccfl[c] = ccfl; 733 adapt->candidates.cost[c] = cost; 734 adapt->candidates.n++; 735 PetscFunctionReturn(0); 736 } 737 738 /*@C 739 TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 740 741 Not Collective 742 743 Input Arguments: 744 . adapt - time step adaptivity context 745 746 Output Arguments: 747 + n - number of candidate schemes, always at least 1 748 . order - the order of each candidate scheme 749 . stageorder - the stage order of each candidate scheme 750 . ccfl - the CFL coefficient of each scheme 751 - cost - the relative cost of each scheme 752 753 Level: developer 754 755 Note: 756 The current scheme is always returned in the first slot 757 758 .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 759 @*/ 760 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 761 { 762 PetscFunctionBegin; 763 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 764 if (n) *n = adapt->candidates.n; 765 if (order) *order = adapt->candidates.order; 766 if (stageorder) *stageorder = adapt->candidates.stageorder; 767 if (ccfl) *ccfl = adapt->candidates.ccfl; 768 if (cost) *cost = adapt->candidates.cost; 769 PetscFunctionReturn(0); 770 } 771 772 /*@C 773 TSAdaptChoose - choose which method and step size to use for the next step 774 775 Collective on TSAdapt 776 777 Input Arguments: 778 + adapt - adaptive contoller 779 - h - current step size 780 781 Output Arguments: 782 + next_sc - optional, scheme to use for the next step 783 . next_h - step size to use for the next step 784 - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 785 786 Note: 787 The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 788 being retried after an initial rejection. 789 790 Level: developer 791 792 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 793 @*/ 794 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 795 { 796 PetscErrorCode ierr; 797 PetscInt ncandidates = adapt->candidates.n; 798 PetscInt scheme = 0; 799 PetscReal wlte = -1.0; 800 PetscReal wltea = -1.0; 801 PetscReal wlter = -1.0; 802 803 PetscFunctionBegin; 804 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 805 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 806 if (next_sc) PetscValidIntPointer(next_sc,4); 807 PetscValidPointer(next_h,5); 808 PetscValidIntPointer(accept,6); 809 if (next_sc) *next_sc = 0; 810 811 /* Do not mess with adaptivity while handling events*/ 812 if (ts->event && ts->event->status != TSEVENT_NONE) { 813 *next_h = h; 814 *accept = PETSC_TRUE; 815 PetscFunctionReturn(0); 816 } 817 818 ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 819 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); 820 if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 821 if (next_sc) *next_sc = scheme; 822 823 if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 824 /* Increase/reduce step size if end time of next step is close to or overshoots max time */ 825 PetscReal t = ts->ptime + ts->time_step, h = *next_h; 826 PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t; 827 PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]); 828 PetscReal b = adapt->matchstepfac[1]; 829 if (t < tmax && tend > tmax) *next_h = hmax; 830 if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2; 831 if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax; 832 } 833 834 if (adapt->monitor) { 835 const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 836 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 837 if (wlte < 0) { 838 ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h);CHKERRQ(ierr); 839 } else { 840 ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e wlte=%5.3g wltea=%5.3g wlter=%5.3g\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h,(double)wlte,(double)wltea,(double)wlter);CHKERRQ(ierr); 841 } 842 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 843 } 844 PetscFunctionReturn(0); 845 } 846 847 /*@ 848 TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver 849 before increasing the time step. 850 851 Logicially Collective on TSAdapt 852 853 Input Arguments: 854 + adapt - adaptive controller context 855 - cnt - the number of timesteps 856 857 Options Database Key: 858 . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase 859 860 Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0. 861 The successful use of this option is problem dependent 862 863 Developer Note: there is no theory to support this option 864 865 Level: advanced 866 867 .seealso: 868 @*/ 869 PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt) 870 { 871 PetscFunctionBegin; 872 adapt->timestepjustdecreased_delay = cnt; 873 PetscFunctionReturn(0); 874 } 875 876 877 /*@ 878 TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 879 880 Collective on TSAdapt 881 882 Input Arguments: 883 + adapt - adaptive controller context 884 . ts - time stepper 885 . t - Current simulation time 886 - Y - Current solution vector 887 888 Output Arguments: 889 . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 890 891 Level: developer 892 893 .seealso: 894 @*/ 895 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 896 { 897 PetscErrorCode ierr; 898 SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 899 900 PetscFunctionBegin; 901 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 902 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 903 PetscValidIntPointer(accept,3); 904 905 if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 906 if (snesreason < 0) { 907 *accept = PETSC_FALSE; 908 if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 909 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 910 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); 911 if (adapt->monitor) { 912 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 913 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); 914 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 915 } 916 } 917 } else { 918 *accept = PETSC_TRUE; 919 ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 920 if(*accept && adapt->checkstage) { 921 ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 922 } 923 } 924 925 if(!(*accept) && !ts->reason) { 926 PetscReal dt,new_dt; 927 ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 928 new_dt = dt * adapt->scale_solve_failed; 929 ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 930 adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay; 931 if (adapt->monitor) { 932 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 933 ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected (%s) t=%-11g+%10.3e retrying with dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,SNESConvergedReasons[snesreason],(double)ts->ptime,(double)dt,(double)new_dt);CHKERRQ(ierr); 934 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 935 } 936 } 937 PetscFunctionReturn(0); 938 } 939 940 /*@ 941 TSAdaptCreate - create an adaptive controller context for time stepping 942 943 Collective on MPI_Comm 944 945 Input Parameter: 946 . comm - The communicator 947 948 Output Parameter: 949 . adapt - new TSAdapt object 950 951 Level: developer 952 953 Notes: 954 TSAdapt creation is handled by TS, so users should not need to call this function. 955 956 .keywords: TSAdapt, create 957 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 958 @*/ 959 PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 960 { 961 PetscErrorCode ierr; 962 TSAdapt adapt; 963 964 PetscFunctionBegin; 965 PetscValidPointer(inadapt,1); 966 *inadapt = NULL; 967 ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 968 969 ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 970 971 adapt->always_accept = PETSC_FALSE; 972 adapt->safety = 0.9; 973 adapt->reject_safety = 0.5; 974 adapt->clip[0] = 0.1; 975 adapt->clip[1] = 10.; 976 adapt->dt_min = 1e-20; 977 adapt->dt_max = 1e+20; 978 adapt->scale_solve_failed = 0.25; 979 /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case 980 to prevent from situations were unreasonably small time steps are taken in order to match the final time */ 981 adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */ 982 adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */ 983 adapt->wnormtype = NORM_2; 984 adapt->timestepjustdecreased_delay = 0; 985 986 *inadapt = adapt; 987 PetscFunctionReturn(0); 988 } 989