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