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