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