1 /* 2 Code for Timestepping with implicit backwards Euler. 3 */ 4 #include "src/ts/tsimpl.h" /*I "petscts.h" I*/ 5 6 typedef struct { 7 Vec update; /* work vector where new solution is formed */ 8 Vec func; /* work vector where F(t[i],u[i]) is stored */ 9 Vec rhs; /* work vector for RHS; vec_sol/dt */ 10 11 /* information used for Pseudo-timestepping */ 12 13 PetscErrorCode (*dt)(TS,PetscReal*,void*); /* compute next timestep, and related context */ 14 void *dtctx; 15 PetscErrorCode (*verify)(TS,Vec,void*,PetscReal*,int*); /* verify previous timestep and related context */ 16 void *verifyctx; 17 18 PetscReal initial_fnorm,fnorm; /* original and current norm of F(u) */ 19 PetscReal fnorm_previous; 20 21 PetscReal dt_increment; /* scaling that dt is incremented each time-step */ 22 PetscTruth increment_dt_from_initial_dt; 23 } TS_Pseudo; 24 25 /* ------------------------------------------------------------------------------*/ 26 27 #undef __FUNCT__ 28 #define __FUNCT__ "TSPseudoComputeTimeStep" 29 /*@ 30 TSPseudoComputeTimeStep - Computes the next timestep for a currently running 31 pseudo-timestepping process. 32 33 Collective on TS 34 35 Input Parameter: 36 . ts - timestep context 37 38 Output Parameter: 39 . dt - newly computed timestep 40 41 Level: advanced 42 43 Notes: 44 The routine to be called here to compute the timestep should be 45 set by calling TSPseudoSetTimeStep(). 46 47 .keywords: timestep, pseudo, compute 48 49 .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 50 @*/ 51 PetscErrorCode TSPseudoComputeTimeStep(TS ts,PetscReal *dt) 52 { 53 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 54 PetscErrorCode ierr; 55 56 PetscFunctionBegin; 57 ierr = PetscLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 58 ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx);CHKERRQ(ierr); 59 ierr = PetscLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0);CHKERRQ(ierr); 60 PetscFunctionReturn(0); 61 } 62 63 64 /* ------------------------------------------------------------------------------*/ 65 #undef __FUNCT__ 66 #define __FUNCT__ "TSPseudoDefaultVerifyTimeStep" 67 /*@C 68 TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 69 70 Collective on TS 71 72 Input Parameters: 73 + ts - the timestep context 74 . dtctx - unused timestep context 75 - update - latest solution vector 76 77 Output Parameters: 78 + newdt - the timestep to use for the next step 79 - flag - flag indicating whether the last time step was acceptable 80 81 Level: advanced 82 83 Note: 84 This routine always returns a flag of 1, indicating an acceptable 85 timestep. 86 87 .keywords: timestep, pseudo, default, verify 88 89 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 90 @*/ 91 PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,PetscReal *newdt,int *flag) 92 { 93 PetscFunctionBegin; 94 *flag = 1; 95 PetscFunctionReturn(0); 96 } 97 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "TSPseudoVerifyTimeStep" 101 /*@ 102 TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 103 104 Collective on TS 105 106 Input Parameters: 107 + ts - timestep context 108 - update - latest solution vector 109 110 Output Parameters: 111 + dt - newly computed timestep (if it had to shrink) 112 - flag - indicates if current timestep was ok 113 114 Level: advanced 115 116 Notes: 117 The routine to be called here to compute the timestep should be 118 set by calling TSPseudoSetVerifyTimeStep(). 119 120 .keywords: timestep, pseudo, verify 121 122 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 123 @*/ 124 PetscErrorCode TSPseudoVerifyTimeStep(TS ts,Vec update,PetscReal *dt,int *flag) 125 { 126 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 127 PetscErrorCode ierr; 128 129 PetscFunctionBegin; 130 if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);} 131 132 ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag);CHKERRQ(ierr); 133 134 PetscFunctionReturn(0); 135 } 136 137 /* --------------------------------------------------------------------------------*/ 138 139 #undef __FUNCT__ 140 #define __FUNCT__ "TSStep_Pseudo" 141 static PetscErrorCode TSStep_Pseudo(TS ts,int *steps,PetscReal *ptime) 142 { 143 Vec sol = ts->vec_sol; 144 PetscErrorCode ierr; 145 int i,max_steps = ts->max_steps,its,ok,lits; 146 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 147 PetscReal current_time_step; 148 149 PetscFunctionBegin; 150 *steps = -ts->steps; 151 152 ierr = VecCopy(sol,pseudo->update);CHKERRQ(ierr); 153 for (i=0; i<max_steps && ts->ptime < ts->max_time; i++) { 154 ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 155 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 156 current_time_step = ts->time_step; 157 while (PETSC_TRUE) { 158 ts->ptime += current_time_step; 159 ierr = SNESSolve(ts->snes,pseudo->update);CHKERRQ(ierr); 160 ierr = SNESGetNumberLinearIterations(ts->snes,&lits);CHKERRQ(ierr); 161 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 162 ts->nonlinear_its += its; ts->linear_its += lits; 163 ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok);CHKERRQ(ierr); 164 if (ok) break; 165 ts->ptime -= current_time_step; 166 current_time_step = ts->time_step; 167 } 168 ierr = VecCopy(pseudo->update,sol);CHKERRQ(ierr); 169 ts->steps++; 170 } 171 ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 172 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 173 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 174 175 *steps += ts->steps; 176 *ptime = ts->ptime; 177 PetscFunctionReturn(0); 178 } 179 180 /*------------------------------------------------------------*/ 181 #undef __FUNCT__ 182 #define __FUNCT__ "TSDestroy_Pseudo" 183 static PetscErrorCode TSDestroy_Pseudo(TS ts) 184 { 185 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 186 PetscErrorCode ierr; 187 188 PetscFunctionBegin; 189 if (pseudo->update) {ierr = VecDestroy(pseudo->update);CHKERRQ(ierr);} 190 if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 191 if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 192 ierr = PetscFree(pseudo);CHKERRQ(ierr); 193 PetscFunctionReturn(0); 194 } 195 196 197 /*------------------------------------------------------------*/ 198 199 /* 200 This defines the nonlinear equation that is to be solved with SNES 201 202 (U^{n+1} - U^{n})/dt - F(U^{n+1}) 203 */ 204 #undef __FUNCT__ 205 #define __FUNCT__ "TSPseudoFunction" 206 PetscErrorCode TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 207 { 208 TS ts = (TS) ctx; 209 PetscScalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 210 PetscErrorCode ierr; 211 int i,n; 212 213 PetscFunctionBegin; 214 /* apply user provided function */ 215 ierr = TSComputeRHSFunction(ts,ts->ptime,x,y);CHKERRQ(ierr); 216 /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 217 ierr = VecGetArray(ts->vec_sol,&un);CHKERRQ(ierr); 218 ierr = VecGetArray(x,&unp1);CHKERRQ(ierr); 219 ierr = VecGetArray(y,&Funp1);CHKERRQ(ierr); 220 ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr); 221 for (i=0; i<n; i++) { 222 Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 223 } 224 ierr = VecRestoreArray(ts->vec_sol,&un);CHKERRQ(ierr); 225 ierr = VecRestoreArray(x,&unp1);CHKERRQ(ierr); 226 ierr = VecRestoreArray(y,&Funp1);CHKERRQ(ierr); 227 228 PetscFunctionReturn(0); 229 } 230 231 /* 232 This constructs the Jacobian needed for SNES 233 234 J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 235 */ 236 #undef __FUNCT__ 237 #define __FUNCT__ "TSPseudoJacobian" 238 PetscErrorCode TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 239 { 240 TS ts = (TS) ctx; 241 PetscErrorCode ierr; 242 PetscScalar mone = -1.0,mdt = 1.0/ts->time_step; 243 244 PetscFunctionBegin; 245 /* construct users Jacobian */ 246 ierr = TSComputeRHSJacobian(ts,ts->ptime,x,AA,BB,str);CHKERRQ(ierr); 247 248 /* shift and scale Jacobian */ 249 ierr = MatScale(&mone,*AA);CHKERRQ(ierr); 250 ierr = MatShift(&mdt,*AA);CHKERRQ(ierr); 251 if (*BB != *AA && *str != SAME_PRECONDITIONER) { 252 ierr = MatScale(&mone,*BB);CHKERRQ(ierr); 253 ierr = MatShift(&mdt,*BB);CHKERRQ(ierr); 254 } 255 256 PetscFunctionReturn(0); 257 } 258 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "TSSetUp_Pseudo" 262 static PetscErrorCode TSSetUp_Pseudo(TS ts) 263 { 264 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 265 PetscErrorCode ierr; 266 267 PetscFunctionBegin; 268 /* ierr = SNESSetFromOptions(ts->snes);CHKERRQ(ierr); */ 269 ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 270 ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 271 ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 272 ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 /*------------------------------------------------------------*/ 276 277 #undef __FUNCT__ 278 #define __FUNCT__ "TSPseudoDefaultMonitor" 279 PetscErrorCode TSPseudoDefaultMonitor(TS ts,int step,PetscReal ptime,Vec v,void *ctx) 280 { 281 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 282 PetscErrorCode ierr; 283 284 PetscFunctionBegin; 285 ierr = (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "TSSetFromOptions_Pseudo" 291 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 292 { 293 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 294 PetscErrorCode ierr; 295 PetscTruth flg; 296 297 PetscFunctionBegin; 298 299 ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 300 ierr = PetscOptionsName("-ts_monitor","Monitor convergence","TSPseudoDefaultMonitor",&flg);CHKERRQ(ierr); 301 if (flg) { 302 ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 303 } 304 ierr = PetscOptionsName("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",&flg);CHKERRQ(ierr); 305 if (flg) { 306 ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 307 } 308 ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 309 ierr = PetscOptionsTail();CHKERRQ(ierr); 310 311 PetscFunctionReturn(0); 312 } 313 314 #undef __FUNCT__ 315 #define __FUNCT__ "TSView_Pseudo" 316 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 317 { 318 PetscFunctionBegin; 319 PetscFunctionReturn(0); 320 } 321 322 /* ----------------------------------------------------------------------------- */ 323 #undef __FUNCT__ 324 #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 325 /*@C 326 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 327 last timestep. 328 329 Collective on TS 330 331 Input Parameters: 332 + ts - timestep context 333 . dt - user-defined function to verify timestep 334 - ctx - [optional] user-defined context for private data 335 for the timestep verification routine (may be PETSC_NULL) 336 337 Level: advanced 338 339 Calling sequence of func: 340 . func (TS ts,Vec update,void *ctx,PetscReal *newdt,int *flag); 341 342 . update - latest solution vector 343 . ctx - [optional] timestep context 344 . newdt - the timestep to use for the next step 345 . flag - flag indicating whether the last time step was acceptable 346 347 Notes: 348 The routine set here will be called by TSPseudoVerifyTimeStep() 349 during the timestepping process. 350 351 .keywords: timestep, pseudo, set, verify 352 353 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 354 @*/ 355 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,int*),void* ctx) 356 { 357 PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,int *),void *); 358 359 PetscFunctionBegin; 360 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 361 362 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 363 if (f) { 364 ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 365 } 366 PetscFunctionReturn(0); 367 } 368 369 #undef __FUNCT__ 370 #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 371 /*@ 372 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 373 dt when using the TSPseudoDefaultTimeStep() routine. 374 375 Collective on TS 376 377 Input Parameters: 378 + ts - the timestep context 379 - inc - the scaling factor >= 1.0 380 381 Options Database Key: 382 $ -ts_pseudo_increment <increment> 383 384 Level: advanced 385 386 .keywords: timestep, pseudo, set, increment 387 388 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 389 @*/ 390 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 391 { 392 PetscErrorCode ierr,(*f)(TS,PetscReal); 393 394 PetscFunctionBegin; 395 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 396 397 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C",(void (**)(void))&f);CHKERRQ(ierr); 398 if (f) { 399 ierr = (*f)(ts,inc);CHKERRQ(ierr); 400 } 401 PetscFunctionReturn(0); 402 } 403 404 #undef __FUNCT__ 405 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 406 /*@ 407 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 408 is computed via the formula 409 $ dt = initial_dt*initial_fnorm/current_fnorm 410 rather than the default update, 411 $ dt = current_dt*previous_fnorm/current_fnorm. 412 413 Collective on TS 414 415 Input Parameter: 416 . ts - the timestep context 417 418 Options Database Key: 419 $ -ts_pseudo_increment_dt_from_initial_dt 420 421 Level: advanced 422 423 .keywords: timestep, pseudo, set, increment 424 425 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 426 @*/ 427 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 428 { 429 PetscErrorCode ierr,(*f)(TS); 430 431 PetscFunctionBegin; 432 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 433 434 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C",(void (**)(void))&f);CHKERRQ(ierr); 435 if (f) { 436 ierr = (*f)(ts);CHKERRQ(ierr); 437 } 438 PetscFunctionReturn(0); 439 } 440 441 442 #undef __FUNCT__ 443 #define __FUNCT__ "TSPseudoSetTimeStep" 444 /*@C 445 TSPseudoSetTimeStep - Sets the user-defined routine to be 446 called at each pseudo-timestep to update the timestep. 447 448 Collective on TS 449 450 Input Parameters: 451 + ts - timestep context 452 . dt - function to compute timestep 453 - ctx - [optional] user-defined context for private data 454 required by the function (may be PETSC_NULL) 455 456 Level: intermediate 457 458 Calling sequence of func: 459 . func (TS ts,PetscReal *newdt,void *ctx); 460 461 . newdt - the newly computed timestep 462 . ctx - [optional] timestep context 463 464 Notes: 465 The routine set here will be called by TSPseudoComputeTimeStep() 466 during the timestepping process. 467 468 .keywords: timestep, pseudo, set 469 470 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 471 @*/ 472 PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 473 { 474 PetscErrorCode ierr,(*f)(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *); 475 476 PetscFunctionBegin; 477 PetscValidHeaderSpecific(ts,TS_COOKIE,1); 478 479 ierr = PetscObjectQueryFunction((PetscObject)ts,"TSPseudoSetTimeStep_C",(void (**)(void))&f);CHKERRQ(ierr); 480 if (f) { 481 ierr = (*f)(ts,dt,ctx);CHKERRQ(ierr); 482 } 483 PetscFunctionReturn(0); 484 } 485 486 /* ----------------------------------------------------------------------------- */ 487 488 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,int*); /* force argument to next function to not be extern C*/ 489 EXTERN_C_BEGIN 490 #undef __FUNCT__ 491 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 492 PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 493 { 494 TS_Pseudo *pseudo; 495 496 PetscFunctionBegin; 497 pseudo = (TS_Pseudo*)ts->data; 498 pseudo->verify = dt; 499 pseudo->verifyctx = ctx; 500 PetscFunctionReturn(0); 501 } 502 EXTERN_C_END 503 504 EXTERN_C_BEGIN 505 #undef __FUNCT__ 506 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 507 PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 508 { 509 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 510 511 PetscFunctionBegin; 512 pseudo->dt_increment = inc; 513 PetscFunctionReturn(0); 514 } 515 EXTERN_C_END 516 517 EXTERN_C_BEGIN 518 #undef __FUNCT__ 519 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 520 PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 521 { 522 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 523 524 PetscFunctionBegin; 525 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 526 PetscFunctionReturn(0); 527 } 528 EXTERN_C_END 529 530 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 531 EXTERN_C_BEGIN 532 #undef __FUNCT__ 533 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 534 PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 535 { 536 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 537 538 PetscFunctionBegin; 539 pseudo->dt = dt; 540 pseudo->dtctx = ctx; 541 PetscFunctionReturn(0); 542 } 543 EXTERN_C_END 544 545 /* ----------------------------------------------------------------------------- */ 546 547 EXTERN_C_BEGIN 548 #undef __FUNCT__ 549 #define __FUNCT__ "TSCreate_Pseudo" 550 PetscErrorCode TSCreate_Pseudo(TS ts) 551 { 552 TS_Pseudo *pseudo; 553 PetscErrorCode ierr; 554 555 PetscFunctionBegin; 556 ts->ops->destroy = TSDestroy_Pseudo; 557 ts->ops->view = TSView_Pseudo; 558 559 if (ts->problem_type == TS_LINEAR) { 560 SETERRQ(PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 561 } 562 if (!ts->A) { 563 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set Jacobian"); 564 } 565 566 ts->ops->setup = TSSetUp_Pseudo; 567 ts->ops->step = TSStep_Pseudo; 568 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 569 570 /* create the required nonlinear solver context */ 571 ierr = SNESCreate(ts->comm,&ts->snes);CHKERRQ(ierr); 572 573 ierr = PetscNew(TS_Pseudo,&pseudo);CHKERRQ(ierr); 574 PetscLogObjectMemory(ts,sizeof(TS_Pseudo)); 575 576 ierr = PetscMemzero(pseudo,sizeof(TS_Pseudo));CHKERRQ(ierr); 577 ts->data = (void*)pseudo; 578 579 pseudo->dt_increment = 1.1; 580 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 581 pseudo->dt = TSPseudoDefaultTimeStep; 582 583 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 584 "TSPseudoSetVerifyTimeStep_Pseudo", 585 TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 586 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 587 "TSPseudoSetTimeStepIncrement_Pseudo", 588 TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 589 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 590 "TSPseudoIncrementDtFromInitialDt_Pseudo", 591 TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 592 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 593 "TSPseudoSetTimeStep_Pseudo", 594 TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 595 PetscFunctionReturn(0); 596 } 597 EXTERN_C_END 598 599 #undef __FUNCT__ 600 #define __FUNCT__ "TSPseudoDefaultTimeStep" 601 /*@C 602 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 603 Use with TSPseudoSetTimeStep(). 604 605 Collective on TS 606 607 Input Parameters: 608 . ts - the timestep context 609 . dtctx - unused timestep context 610 611 Output Parameter: 612 . newdt - the timestep to use for the next step 613 614 Level: advanced 615 616 .keywords: timestep, pseudo, default 617 618 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 619 @*/ 620 PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 621 { 622 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 623 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 624 PetscErrorCode ierr; 625 626 PetscFunctionBegin; 627 ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 628 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 629 if (pseudo->initial_fnorm == 0.0) { 630 /* first time through so compute initial function norm */ 631 pseudo->initial_fnorm = pseudo->fnorm; 632 fnorm_previous = pseudo->fnorm; 633 } 634 if (pseudo->fnorm == 0.0) { 635 *newdt = 1.e12*inc*ts->time_step; 636 } else if (pseudo->increment_dt_from_initial_dt) { 637 *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 638 } else { 639 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 640 } 641 pseudo->fnorm_previous = pseudo->fnorm; 642 PetscFunctionReturn(0); 643 } 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661