1 2 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 3 4 static PetscFList TSAdaptList; 5 static PetscBool TSAdaptPackageInitialized; 6 static PetscBool TSAdaptRegisterAllCalled; 7 static PetscClassId TSADAPT_CLASSID; 8 9 EXTERN_C_BEGIN 10 PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 11 PetscErrorCode TSAdaptCreate_None(TSAdapt); 12 PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 13 EXTERN_C_END 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "TSAdaptRegister" 17 /*@C 18 TSAdaptRegister - see TSAdaptRegisterDynamic() 19 20 Level: advanced 21 @*/ 22 PetscErrorCode TSAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSAdapt)) 23 { 24 PetscErrorCode ierr; 25 char fullname[PETSC_MAX_PATH_LEN]; 26 27 PetscFunctionBegin; 28 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 29 ierr = PetscFListAdd(&TSAdaptList,sname,fullname,(void(*)(void))function);CHKERRQ(ierr); 30 PetscFunctionReturn(0); 31 } 32 33 #undef __FUNCT__ 34 #define __FUNCT__ "TSAdaptRegisterAll" 35 /*@C 36 TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 37 38 Not Collective 39 40 Level: advanced 41 42 .keywords: TSAdapt, register, all 43 44 .seealso: TSAdaptRegisterDestroy() 45 @*/ 46 PetscErrorCode TSAdaptRegisterAll(const char path[]) 47 { 48 PetscErrorCode ierr; 49 50 PetscFunctionBegin; 51 ierr = TSAdaptRegisterDynamic(TSADAPTBASIC,path,"TSAdaptCreate_Basic",TSAdaptCreate_Basic);CHKERRQ(ierr); 52 ierr = TSAdaptRegisterDynamic(TSADAPTNONE, path,"TSAdaptCreate_None", TSAdaptCreate_None);CHKERRQ(ierr); 53 ierr = TSAdaptRegisterDynamic(TSADAPTCFL, path,"TSAdaptCreate_CFL", TSAdaptCreate_CFL);CHKERRQ(ierr); 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "TSAdaptFinalizePackage" 59 /*@C 60 TSFinalizePackage - This function destroys everything in the TS package. It is 61 called from PetscFinalize(). 62 63 Level: developer 64 65 .keywords: Petsc, destroy, package 66 .seealso: PetscFinalize() 67 @*/ 68 PetscErrorCode TSAdaptFinalizePackage(void) 69 { 70 PetscFunctionBegin; 71 TSAdaptPackageInitialized = PETSC_FALSE; 72 TSAdaptRegisterAllCalled = PETSC_FALSE; 73 TSAdaptList = PETSC_NULL; 74 PetscFunctionReturn(0); 75 } 76 77 #undef __FUNCT__ 78 #define __FUNCT__ "TSAdaptInitializePackage" 79 /*@C 80 TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 81 called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 82 TSCreate_GL() when using static libraries. 83 84 Input Parameter: 85 path - The dynamic library path, or PETSC_NULL 86 87 Level: developer 88 89 .keywords: TSAdapt, initialize, package 90 .seealso: PetscInitialize() 91 @*/ 92 PetscErrorCode TSAdaptInitializePackage(const char path[]) 93 { 94 PetscErrorCode ierr; 95 96 PetscFunctionBegin; 97 if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 98 TSAdaptPackageInitialized = PETSC_TRUE; 99 ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 100 ierr = TSAdaptRegisterAll(path);CHKERRQ(ierr); 101 ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 102 PetscFunctionReturn(0); 103 } 104 105 #undef __FUNCT__ 106 #define __FUNCT__ "TSAdaptRegisterDestroy" 107 /*@C 108 TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister()/TSAdaptRegisterDynamic(). 109 110 Not Collective 111 112 Level: advanced 113 114 .keywords: TSAdapt, register, destroy 115 .seealso: TSAdaptRegister(), TSAdaptRegisterAll(), TSAdaptRegisterDynamic() 116 @*/ 117 PetscErrorCode TSAdaptRegisterDestroy(void) 118 { 119 PetscErrorCode ierr; 120 121 PetscFunctionBegin; 122 ierr = PetscFListDestroy(&TSAdaptList);CHKERRQ(ierr); 123 TSAdaptRegisterAllCalled = PETSC_FALSE; 124 PetscFunctionReturn(0); 125 } 126 127 128 #undef __FUNCT__ 129 #define __FUNCT__ "TSAdaptSetType" 130 PetscErrorCode TSAdaptSetType(TSAdapt adapt,const TSAdaptType type) 131 { 132 PetscErrorCode ierr,(*r)(TSAdapt); 133 134 PetscFunctionBegin; 135 ierr = PetscFListFind(TSAdaptList,((PetscObject)adapt)->comm,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr); 136 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 137 if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);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__ "TSAdaptView" 156 PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 157 { 158 PetscErrorCode ierr; 159 PetscBool iascii; 160 161 PetscFunctionBegin; 162 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 163 if (iascii) { 164 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");CHKERRQ(ierr); 165 ierr = PetscViewerASCIIPrintf(viewer,"number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 166 if (adapt->ops->view) { 167 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 168 ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 169 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 170 } 171 } else { 172 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name); 173 } 174 PetscFunctionReturn(0); 175 } 176 177 #undef __FUNCT__ 178 #define __FUNCT__ "TSAdaptDestroy" 179 PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 180 { 181 PetscErrorCode ierr; 182 183 PetscFunctionBegin; 184 if (!*adapt) PetscFunctionReturn(0); 185 PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 186 if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);} 187 if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 188 ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 189 ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 190 PetscFunctionReturn(0); 191 } 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "TSAdaptSetMonitor" 195 /*@ 196 TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 197 198 Collective on TSAdapt 199 200 Input Arguments: 201 + adapt - adaptive controller context 202 - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 203 204 Level: intermediate 205 206 .seealso: TSAdaptChoose() 207 @*/ 208 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 209 { 210 PetscErrorCode ierr; 211 212 PetscFunctionBegin; 213 if (flg) { 214 if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(((PetscObject)adapt)->comm,"stdout",&adapt->monitor);CHKERRQ(ierr);} 215 } else { 216 ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 217 } 218 PetscFunctionReturn(0); 219 } 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "TSAdaptSetCheckStage" 223 /*@C 224 TSAdaptSetCheckStage - set a callback to check convergence for a stage 225 226 Logically collective on TSAdapt 227 228 Input Arguments: 229 + adapt - adaptive controller context 230 - func - stage check function 231 232 Arguments of func: 233 $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 234 235 + adapt - adaptive controller context 236 . ts - time stepping context 237 - accept - pending choice of whether to accept, can be modified by this routine 238 239 Level: advanced 240 241 .seealso: TSAdaptChoose() 242 @*/ 243 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*)) 244 { 245 246 PetscFunctionBegin; 247 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 248 adapt->ops->checkstage = func; 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNCT__ 253 #define __FUNCT__ "TSAdaptSetStepLimits" 254 /*@ 255 TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 256 257 Logically Collective 258 259 Input Arguments: 260 + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 261 . hmin - minimum time step 262 - hmax - maximum time step 263 264 Options Database Keys: 265 + -ts_adapt_dt_min - minimum time step 266 - -ts_adapt_dt_max - maximum time step 267 268 Level: intermediate 269 270 .seealso: TSAdapt 271 @*/ 272 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 273 { 274 275 PetscFunctionBegin; 276 if (hmin != PETSC_DECIDE) adapt->dt_min = hmin; 277 if (hmax != PETSC_DECIDE) adapt->dt_max = hmax; 278 PetscFunctionReturn(0); 279 } 280 281 #undef __FUNCT__ 282 #define __FUNCT__ "TSAdaptSetFromOptions" 283 /*@ 284 TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 285 286 Collective on TSAdapt 287 288 Input Parameter: 289 . adapt - the TSAdapt context 290 291 Options Database Keys: 292 . -ts_adapt_type <type> - basic 293 294 Level: advanced 295 296 Notes: 297 This function is automatically called by TSSetFromOptions() 298 299 .keywords: TS, TSGetAdapt(), TSAdaptSetType() 300 301 .seealso: TSGetType() 302 @*/ 303 PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt) 304 { 305 PetscErrorCode ierr; 306 char type[256] = TSADAPTBASIC; 307 PetscBool set,flg; 308 309 PetscFunctionBegin; 310 /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 311 * function can only be called from inside TSSetFromOptions_GL() */ 312 ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr); 313 ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList, 314 ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof type,&flg);CHKERRQ(ierr); 315 if (flg || !((PetscObject)adapt)->type_name) { 316 ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 317 } 318 if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);} 319 ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,PETSC_NULL);CHKERRQ(ierr); 320 ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,PETSC_NULL);CHKERRQ(ierr); 321 ierr = PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,PETSC_NULL);CHKERRQ(ierr); 322 ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 323 if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 324 ierr = PetscOptionsTail();CHKERRQ(ierr); 325 PetscFunctionReturn(0); 326 } 327 328 #undef __FUNCT__ 329 #define __FUNCT__ "TSAdaptCandidatesClear" 330 /*@ 331 TSAdaptCandidatesClear - clear any previously set candidate schemes 332 333 Logically Collective 334 335 Input Argument: 336 . adapt - adaptive controller 337 338 Level: developer 339 340 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 341 @*/ 342 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 343 { 344 PetscErrorCode ierr; 345 346 PetscFunctionBegin; 347 ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 348 PetscFunctionReturn(0); 349 } 350 351 #undef __FUNCT__ 352 #define __FUNCT__ "TSAdaptCandidateAdd" 353 /*@C 354 TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 355 356 Logically Collective 357 358 Input Arguments: 359 + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 360 . name - name of the candidate scheme to add 361 . order - order of the candidate scheme 362 . stageorder - stage order of the candidate scheme 363 . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 364 . cost - relative measure of the amount of work required for the candidate scheme 365 - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 366 367 Note: 368 This routine is not available in Fortran. 369 370 Level: developer 371 372 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 373 @*/ 374 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 375 { 376 PetscInt c; 377 378 PetscFunctionBegin; 379 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 380 if (order < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 381 if (inuse) { 382 if (adapt->candidates.inuse_set) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()"); 383 adapt->candidates.inuse_set = PETSC_TRUE; 384 } 385 /* first slot if this is the current scheme, otherwise the next available slot */ 386 c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 387 adapt->candidates.name[c] = name; 388 adapt->candidates.order[c] = order; 389 adapt->candidates.stageorder[c] = stageorder; 390 adapt->candidates.ccfl[c] = ccfl; 391 adapt->candidates.cost[c] = cost; 392 adapt->candidates.n++; 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "TSAdaptCandidatesGet" 398 /*@C 399 TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 400 401 Not Collective 402 403 Input Arguments: 404 . adapt - time step adaptivity context 405 406 Output Arguments: 407 + n - number of candidate schemes, always at least 1 408 . order - the order of each candidate scheme 409 . stageorder - the stage order of each candidate scheme 410 . ccfl - the CFL coefficient of each scheme 411 - cost - the relative cost of each scheme 412 413 Level: developer 414 415 Note: 416 The current scheme is always returned in the first slot 417 418 .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 419 @*/ 420 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 421 { 422 PetscFunctionBegin; 423 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 424 if (n) *n = adapt->candidates.n; 425 if (order) *order = adapt->candidates.order; 426 if (stageorder) *stageorder = adapt->candidates.stageorder; 427 if (ccfl) *ccfl = adapt->candidates.ccfl; 428 if (cost) *cost = adapt->candidates.cost; 429 PetscFunctionReturn(0); 430 } 431 432 #undef __FUNCT__ 433 #define __FUNCT__ "TSAdaptChoose" 434 /*@C 435 TSAdaptChoose - choose which method and step size to use for the next step 436 437 Logically Collective 438 439 Input Arguments: 440 + adapt - adaptive contoller 441 - h - current step size 442 443 Output Arguments: 444 + next_sc - scheme to use for the next step 445 . next_h - step size to use for the next step 446 - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 447 448 Note: 449 The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 450 being retried after an initial rejection. 451 452 Level: developer 453 454 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 455 @*/ 456 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 457 { 458 PetscErrorCode ierr; 459 PetscReal wlte = -1.0; 460 461 PetscFunctionBegin; 462 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 463 PetscValidHeaderSpecific(ts,TS_CLASSID,2); 464 PetscValidIntPointer(next_sc,4); 465 PetscValidPointer(next_h,5); 466 PetscValidIntPointer(accept,6); 467 if (adapt->candidates.n < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n); 468 if (!adapt->candidates.inuse_set) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n); 469 ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr); 470 if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1); 471 if (!(*next_h > 0.)) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h); 472 473 if (adapt->monitor) { 474 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 475 if (wlte < 0) { 476 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); 477 } else { 478 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); 479 } 480 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 481 } 482 PetscFunctionReturn(0); 483 } 484 485 #undef __FUNCT__ 486 #define __FUNCT__ "TSAdaptCheckStage" 487 /*@ 488 TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 489 490 Collective 491 492 Input Arguments: 493 + adapt - adaptive controller context 494 - ts - time stepper 495 496 Output Arguments: 497 . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 498 499 Level: developer 500 501 .seealso: 502 @*/ 503 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept) 504 { 505 PetscErrorCode ierr; 506 SNES snes; 507 SNESConvergedReason snesreason; 508 509 PetscFunctionBegin; 510 *accept = PETSC_TRUE; 511 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 512 ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr); 513 if (snesreason < 0) { 514 PetscReal dt,new_dt; 515 *accept = PETSC_FALSE; 516 ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 517 if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 518 ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 519 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); 520 if (adapt->monitor) { 521 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 522 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); 523 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 524 } 525 } else { 526 new_dt = dt*adapt->scale_solve_failed; 527 ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 528 if (adapt->monitor) { 529 ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 530 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); 531 ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 532 } 533 } 534 } 535 if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);} 536 PetscFunctionReturn(0); 537 } 538 539 540 541 #undef __FUNCT__ 542 #define __FUNCT__ "TSAdaptCreate" 543 /*@ 544 TSAdaptCreate - create an adaptive controller context for time stepping 545 546 Collective on MPI_Comm 547 548 Input Parameter: 549 . comm - The communicator 550 551 Output Parameter: 552 . adapt - new TSAdapt object 553 554 Level: developer 555 556 Notes: 557 TSAdapt creation is handled by TS, so users should not need to call this function. 558 559 .keywords: TSAdapt, create 560 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 561 @*/ 562 PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 563 { 564 PetscErrorCode ierr; 565 TSAdapt adapt; 566 567 PetscFunctionBegin; 568 *inadapt = 0; 569 ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,0,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 570 571 adapt->dt_min = 1e-20; 572 adapt->dt_max = 1e50; 573 adapt->scale_solve_failed = 0.25; 574 575 *inadapt = adapt; 576 PetscFunctionReturn(0); 577 } 578