1 /* 2 Code for Timestepping with implicit backwards Euler. 3 */ 4 #include <private/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 xdot; /* work vector for time derivative of state */ 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*,PetscBool *); /* 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 PetscBool 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,PetscBool *flag) 92 { 93 PetscFunctionBegin; 94 *flag = PETSC_TRUE; 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,PetscBool *flag) 125 { 126 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 127 PetscErrorCode ierr; 128 129 PetscFunctionBegin; 130 if (!pseudo->verify) {*flag = PETSC_TRUE; 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,PetscInt *steps,PetscReal *ptime) 142 { 143 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 144 Vec sol = ts->vec_sol, update = pseudo->update; 145 PetscErrorCode ierr; 146 PetscInt i,its,lits; 147 PetscBool ok; 148 PetscReal current_time_step; 149 150 PetscFunctionBegin; 151 *steps = -ts->steps; 152 *ptime = ts->ptime; 153 154 ierr = VecCopy(sol,update);CHKERRQ(ierr); 155 for (i=0; i<ts->max_steps && ts->ptime < ts->max_time; i++) { 156 157 ierr = TSPseudoComputeTimeStep(ts,&ts->time_step);CHKERRQ(ierr); 158 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 159 current_time_step = ts->time_step; 160 ierr = TSPreStep(ts);CHKERRQ(ierr); 161 while (PETSC_TRUE) { 162 ts->ptime += current_time_step; 163 ierr = SNESSolve(ts->snes,PETSC_NULL,update);CHKERRQ(ierr); 164 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 165 ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 166 ts->nonlinear_its += its; ts->linear_its += lits; 167 ierr = TSPseudoVerifyTimeStep(ts,update,&ts->time_step,&ok);CHKERRQ(ierr); 168 if (ok) break; 169 ts->ptime -= current_time_step; 170 current_time_step = ts->time_step; 171 } 172 ierr = VecCopy(update,sol);CHKERRQ(ierr); 173 ts->steps++; 174 ierr = TSPostStep(ts);CHKERRQ(ierr); 175 } 176 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 177 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 178 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 179 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 180 181 *steps += ts->steps; 182 *ptime = ts->ptime; 183 PetscFunctionReturn(0); 184 } 185 186 /*------------------------------------------------------------*/ 187 #undef __FUNCT__ 188 #define __FUNCT__ "TSReset_Pseudo" 189 static PetscErrorCode TSReset_Pseudo(TS ts) 190 { 191 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 ierr = VecDestroy(&pseudo->update);CHKERRQ(ierr); 196 ierr = VecDestroy(&pseudo->func);CHKERRQ(ierr); 197 ierr = VecDestroy(&pseudo->xdot);CHKERRQ(ierr); 198 PetscFunctionReturn(0); 199 } 200 201 #undef __FUNCT__ 202 #define __FUNCT__ "TSDestroy_Pseudo" 203 static PetscErrorCode TSDestroy_Pseudo(TS ts) 204 { 205 PetscErrorCode ierr; 206 207 PetscFunctionBegin; 208 ierr = TSReset_Pseudo(ts);CHKERRQ(ierr); 209 ierr = PetscFree(ts->data);CHKERRQ(ierr); 210 PetscFunctionReturn(0); 211 } 212 213 /*------------------------------------------------------------*/ 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "TSPseudoGetXdot" 217 /* 218 Compute Xdot = (X^{n+1}-X^n)/dt) = 0 219 */ 220 static PetscErrorCode TSPseudoGetXdot(TS ts,Vec X,Vec *Xdot) 221 { 222 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 223 PetscScalar mdt = 1.0/ts->time_step,*xnp1,*xn,*xdot; 224 PetscErrorCode ierr; 225 PetscInt i,n; 226 227 PetscFunctionBegin; 228 ierr = VecGetArray(ts->vec_sol,&xn);CHKERRQ(ierr); 229 ierr = VecGetArray(X,&xnp1);CHKERRQ(ierr); 230 ierr = VecGetArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 231 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 232 for (i=0; i<n; i++) { 233 xdot[i] = mdt*(xnp1[i] - xn[i]); 234 } 235 ierr = VecRestoreArray(ts->vec_sol,&xn);CHKERRQ(ierr); 236 ierr = VecRestoreArray(X,&xnp1);CHKERRQ(ierr); 237 ierr = VecRestoreArray(pseudo->xdot,&xdot);CHKERRQ(ierr); 238 *Xdot = pseudo->xdot; 239 PetscFunctionReturn(0); 240 } 241 242 #undef __FUNCT__ 243 #define __FUNCT__ "SNESTSFormFunction_Pseudo" 244 /* 245 The transient residual is 246 247 F(U^{n+1},(U^{n+1}-U^n)/dt) = 0 248 249 or for ODE, 250 251 (U^{n+1} - U^{n})/dt - F(U^{n+1}) = 0 252 253 This is the function that must be evaluated for transient simulation and for 254 finite difference Jacobians. On the first Newton step, this algorithm uses 255 a guess of U^{n+1} = U^n in which case the transient term vanishes and the 256 residual is actually the steady state residual. Pseudotransient 257 continuation as described in the literature is a linearly implicit 258 algorithm, it only takes this one Newton step with the steady state 259 residual, and then advances to the next time step. 260 */ 261 static PetscErrorCode SNESTSFormFunction_Pseudo(SNES snes,Vec X,Vec Y,TS ts) 262 { 263 Vec Xdot; 264 PetscErrorCode ierr; 265 266 PetscFunctionBegin; 267 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 268 ierr = TSComputeIFunction(ts,ts->ptime,X,Xdot,Y);CHKERRQ(ierr); 269 PetscFunctionReturn(0); 270 } 271 272 #undef __FUNCT__ 273 #define __FUNCT__ "SNESTSFormJacobian_Pseudo" 274 /* 275 This constructs the Jacobian needed for SNES. For DAE, this is 276 277 dF(X,Xdot)/dX + shift*dF(X,Xdot)/dXdot 278 279 and for ODE: 280 281 J = I/dt - J_{Frhs} where J_{Frhs} is the given Jacobian of Frhs. 282 */ 283 static PetscErrorCode SNESTSFormJacobian_Pseudo(SNES snes,Vec X,Mat *AA,Mat *BB,MatStructure *str,TS ts) 284 { 285 Vec Xdot; 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 ierr = TSPseudoGetXdot(ts,X,&Xdot);CHKERRQ(ierr); 290 ierr = TSComputeIJacobian(ts,ts->ptime,X,Xdot,1./ts->time_step,AA,BB,str);CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 295 #undef __FUNCT__ 296 #define __FUNCT__ "TSSetUp_Pseudo" 297 static PetscErrorCode TSSetUp_Pseudo(TS ts) 298 { 299 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 if (ts->problem_type != TS_NONLINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Only for nonlinear problems"); 304 305 ierr = VecDuplicate(ts->vec_sol,&pseudo->update);CHKERRQ(ierr); 306 ierr = VecDuplicate(ts->vec_sol,&pseudo->func);CHKERRQ(ierr); 307 ierr = VecDuplicate(ts->vec_sol,&pseudo->xdot);CHKERRQ(ierr); 308 PetscFunctionReturn(0); 309 } 310 /*------------------------------------------------------------*/ 311 312 #undef __FUNCT__ 313 #define __FUNCT__ "TSPseudoMonitorDefault" 314 PetscErrorCode TSPseudoMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 315 { 316 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 317 PetscErrorCode ierr; 318 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 319 320 PetscFunctionBegin; 321 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 322 ierr = PetscViewerASCIIPrintf(viewer,"TS %D dt %G time %G fnorm %G\n",step,ts->time_step,ptime,pseudo->fnorm);CHKERRQ(ierr); 323 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 324 PetscFunctionReturn(0); 325 } 326 327 #undef __FUNCT__ 328 #define __FUNCT__ "TSSetFromOptions_Pseudo" 329 static PetscErrorCode TSSetFromOptions_Pseudo(TS ts) 330 { 331 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 332 PetscErrorCode ierr; 333 PetscBool flg = PETSC_FALSE; 334 PetscViewer viewer; 335 336 PetscFunctionBegin; 337 ierr = PetscOptionsHead("Pseudo-timestepping options");CHKERRQ(ierr); 338 ierr = PetscOptionsBool("-ts_monitor","Monitor convergence","TSPseudoMonitorDefault",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 339 if (flg) { 340 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,"stdout",&viewer);CHKERRQ(ierr); 341 ierr = TSMonitorSet(ts,TSPseudoMonitorDefault,viewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 342 } 343 flg = PETSC_FALSE; 344 ierr = PetscOptionsBool("-ts_pseudo_increment_dt_from_initial_dt","Increase dt as a ratio from original dt","TSPseudoIncrementDtFromInitialDt",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 345 if (flg) { 346 ierr = TSPseudoIncrementDtFromInitialDt(ts);CHKERRQ(ierr); 347 } 348 ierr = PetscOptionsReal("-ts_pseudo_increment","Ratio to increase dt","TSPseudoSetTimeStepIncrement",pseudo->dt_increment,&pseudo->dt_increment,0);CHKERRQ(ierr); 349 ierr = PetscOptionsTail();CHKERRQ(ierr); 350 PetscFunctionReturn(0); 351 } 352 353 #undef __FUNCT__ 354 #define __FUNCT__ "TSView_Pseudo" 355 static PetscErrorCode TSView_Pseudo(TS ts,PetscViewer viewer) 356 { 357 PetscFunctionBegin; 358 PetscFunctionReturn(0); 359 } 360 361 /* ----------------------------------------------------------------------------- */ 362 #undef __FUNCT__ 363 #define __FUNCT__ "TSPseudoSetVerifyTimeStep" 364 /*@C 365 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 366 last timestep. 367 368 Logically Collective on TS 369 370 Input Parameters: 371 + ts - timestep context 372 . dt - user-defined function to verify timestep 373 - ctx - [optional] user-defined context for private data 374 for the timestep verification routine (may be PETSC_NULL) 375 376 Level: advanced 377 378 Calling sequence of func: 379 . func (TS ts,Vec update,void *ctx,PetscReal *newdt,PetscBool *flag); 380 381 . update - latest solution vector 382 . ctx - [optional] timestep context 383 . newdt - the timestep to use for the next step 384 . flag - flag indicating whether the last time step was acceptable 385 386 Notes: 387 The routine set here will be called by TSPseudoVerifyTimeStep() 388 during the timestepping process. 389 390 .keywords: timestep, pseudo, set, verify 391 392 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 393 @*/ 394 PetscErrorCode TSPseudoSetVerifyTimeStep(TS ts,PetscErrorCode (*dt)(TS,Vec,void*,PetscReal*,PetscBool *),void* ctx) 395 { 396 PetscErrorCode ierr; 397 398 PetscFunctionBegin; 399 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 400 ierr = PetscTryMethod(ts,"TSPseudoSetVerifyTimeStep_C",(TS,PetscErrorCode (*)(TS,Vec,void*,PetscReal *,PetscBool *),void *),(ts,dt,ctx));CHKERRQ(ierr); 401 PetscFunctionReturn(0); 402 } 403 404 #undef __FUNCT__ 405 #define __FUNCT__ "TSPseudoSetTimeStepIncrement" 406 /*@ 407 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 408 dt when using the TSPseudoDefaultTimeStep() routine. 409 410 Logically Collective on TS 411 412 Input Parameters: 413 + ts - the timestep context 414 - inc - the scaling factor >= 1.0 415 416 Options Database Key: 417 $ -ts_pseudo_increment <increment> 418 419 Level: advanced 420 421 .keywords: timestep, pseudo, set, increment 422 423 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 424 @*/ 425 PetscErrorCode TSPseudoSetTimeStepIncrement(TS ts,PetscReal inc) 426 { 427 PetscErrorCode ierr; 428 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 431 PetscValidLogicalCollectiveReal(ts,inc,2); 432 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStepIncrement_C",(TS,PetscReal),(ts,inc));CHKERRQ(ierr); 433 PetscFunctionReturn(0); 434 } 435 436 #undef __FUNCT__ 437 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt" 438 /*@ 439 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 440 is computed via the formula 441 $ dt = initial_dt*initial_fnorm/current_fnorm 442 rather than the default update, 443 $ dt = current_dt*previous_fnorm/current_fnorm. 444 445 Logically Collective on TS 446 447 Input Parameter: 448 . ts - the timestep context 449 450 Options Database Key: 451 $ -ts_pseudo_increment_dt_from_initial_dt 452 453 Level: advanced 454 455 .keywords: timestep, pseudo, set, increment 456 457 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 458 @*/ 459 PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS ts) 460 { 461 PetscErrorCode ierr; 462 463 PetscFunctionBegin; 464 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 465 ierr = PetscTryMethod(ts,"TSPseudoIncrementDtFromInitialDt_C",(TS),(ts));CHKERRQ(ierr); 466 PetscFunctionReturn(0); 467 } 468 469 470 #undef __FUNCT__ 471 #define __FUNCT__ "TSPseudoSetTimeStep" 472 /*@C 473 TSPseudoSetTimeStep - Sets the user-defined routine to be 474 called at each pseudo-timestep to update the timestep. 475 476 Logically Collective on TS 477 478 Input Parameters: 479 + ts - timestep context 480 . dt - function to compute timestep 481 - ctx - [optional] user-defined context for private data 482 required by the function (may be PETSC_NULL) 483 484 Level: intermediate 485 486 Calling sequence of func: 487 . func (TS ts,PetscReal *newdt,void *ctx); 488 489 . newdt - the newly computed timestep 490 . ctx - [optional] timestep context 491 492 Notes: 493 The routine set here will be called by TSPseudoComputeTimeStep() 494 during the timestepping process. 495 496 .keywords: timestep, pseudo, set 497 498 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 499 @*/ 500 PetscErrorCode TSPseudoSetTimeStep(TS ts,PetscErrorCode (*dt)(TS,PetscReal*,void*),void* ctx) 501 { 502 PetscErrorCode ierr; 503 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 506 ierr = PetscTryMethod(ts,"TSPseudoSetTimeStep_C",(TS,PetscErrorCode (*)(TS,PetscReal *,void *),void *),(ts,dt,ctx));CHKERRQ(ierr); 507 PetscFunctionReturn(0); 508 } 509 510 /* ----------------------------------------------------------------------------- */ 511 512 typedef PetscErrorCode (*FCN1)(TS,Vec,void*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 513 EXTERN_C_BEGIN 514 #undef __FUNCT__ 515 #define __FUNCT__ "TSPseudoSetVerifyTimeStep_Pseudo" 516 PetscErrorCode TSPseudoSetVerifyTimeStep_Pseudo(TS ts,FCN1 dt,void* ctx) 517 { 518 TS_Pseudo *pseudo; 519 520 PetscFunctionBegin; 521 pseudo = (TS_Pseudo*)ts->data; 522 pseudo->verify = dt; 523 pseudo->verifyctx = ctx; 524 PetscFunctionReturn(0); 525 } 526 EXTERN_C_END 527 528 EXTERN_C_BEGIN 529 #undef __FUNCT__ 530 #define __FUNCT__ "TSPseudoSetTimeStepIncrement_Pseudo" 531 PetscErrorCode TSPseudoSetTimeStepIncrement_Pseudo(TS ts,PetscReal inc) 532 { 533 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 534 535 PetscFunctionBegin; 536 pseudo->dt_increment = inc; 537 PetscFunctionReturn(0); 538 } 539 EXTERN_C_END 540 541 EXTERN_C_BEGIN 542 #undef __FUNCT__ 543 #define __FUNCT__ "TSPseudoIncrementDtFromInitialDt_Pseudo" 544 PetscErrorCode TSPseudoIncrementDtFromInitialDt_Pseudo(TS ts) 545 { 546 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 547 548 PetscFunctionBegin; 549 pseudo->increment_dt_from_initial_dt = PETSC_TRUE; 550 PetscFunctionReturn(0); 551 } 552 EXTERN_C_END 553 554 typedef PetscErrorCode (*FCN2)(TS,PetscReal*,void*); /* force argument to next function to not be extern C*/ 555 EXTERN_C_BEGIN 556 #undef __FUNCT__ 557 #define __FUNCT__ "TSPseudoSetTimeStep_Pseudo" 558 PetscErrorCode TSPseudoSetTimeStep_Pseudo(TS ts,FCN2 dt,void* ctx) 559 { 560 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 561 562 PetscFunctionBegin; 563 pseudo->dt = dt; 564 pseudo->dtctx = ctx; 565 PetscFunctionReturn(0); 566 } 567 EXTERN_C_END 568 569 /* ----------------------------------------------------------------------------- */ 570 /*MC 571 TSPSEUDO - Solve steady state ODE and DAE problems with pseudo time stepping 572 573 This method solves equations of the form 574 575 $ F(X,Xdot) = 0 576 577 for steady state using the iteration 578 579 $ [G'] S = -F(X,0) 580 $ X += S 581 582 where 583 584 $ G(Y) = F(Y,(Y-X)/dt) 585 586 This is linearly-implicit Euler with the residual always evaluated "at steady 587 state". See note below. 588 589 Options database keys: 590 + -ts_pseudo_increment <real> - ratio of increase dt 591 - -ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt 592 593 Level: beginner 594 595 References: 596 Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential-Algebraic Equations, 2003. 597 C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998. 598 599 Notes: 600 The residual computed by this method includes the transient term (Xdot is computed instead of 601 always being zero), but since the prediction from the last step is always the solution from the 602 last step, on the first Newton iteration we have 603 604 $ Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0 605 606 Therefore, the linear system solved by the first Newton iteration is equivalent to the one 607 described above and in the papers. If the user chooses to perform multiple Newton iterations, the 608 algorithm is no longer the one described in the referenced papers. 609 610 .seealso: TSCreate(), TS, TSSetType() 611 612 M*/ 613 EXTERN_C_BEGIN 614 #undef __FUNCT__ 615 #define __FUNCT__ "TSCreate_Pseudo" 616 PetscErrorCode TSCreate_Pseudo(TS ts) 617 { 618 TS_Pseudo *pseudo; 619 PetscErrorCode ierr; 620 621 PetscFunctionBegin; 622 ts->ops->reset = TSReset_Pseudo; 623 ts->ops->destroy = TSDestroy_Pseudo; 624 ts->ops->view = TSView_Pseudo; 625 626 ts->ops->setup = TSSetUp_Pseudo; 627 ts->ops->step = TSStep_Pseudo; 628 ts->ops->setfromoptions = TSSetFromOptions_Pseudo; 629 ts->ops->snesfunction = SNESTSFormFunction_Pseudo; 630 ts->ops->snesjacobian = SNESTSFormJacobian_Pseudo; 631 632 /* create the required nonlinear solver context */ 633 ts->problem_type = TS_NONLINEAR; 634 ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 635 636 ierr = PetscNewLog(ts,TS_Pseudo,&pseudo);CHKERRQ(ierr); 637 ts->data = (void*)pseudo; 638 639 pseudo->dt_increment = 1.1; 640 pseudo->increment_dt_from_initial_dt = PETSC_FALSE; 641 pseudo->dt = TSPseudoDefaultTimeStep; 642 643 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetVerifyTimeStep_C", 644 "TSPseudoSetVerifyTimeStep_Pseudo", 645 TSPseudoSetVerifyTimeStep_Pseudo);CHKERRQ(ierr); 646 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStepIncrement_C", 647 "TSPseudoSetTimeStepIncrement_Pseudo", 648 TSPseudoSetTimeStepIncrement_Pseudo);CHKERRQ(ierr); 649 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoIncrementDtFromInitialDt_C", 650 "TSPseudoIncrementDtFromInitialDt_Pseudo", 651 TSPseudoIncrementDtFromInitialDt_Pseudo);CHKERRQ(ierr); 652 ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPseudoSetTimeStep_C", 653 "TSPseudoSetTimeStep_Pseudo", 654 TSPseudoSetTimeStep_Pseudo);CHKERRQ(ierr); 655 PetscFunctionReturn(0); 656 } 657 EXTERN_C_END 658 659 #undef __FUNCT__ 660 #define __FUNCT__ "TSPseudoDefaultTimeStep" 661 /*@C 662 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 663 Use with TSPseudoSetTimeStep(). 664 665 Collective on TS 666 667 Input Parameters: 668 . ts - the timestep context 669 . dtctx - unused timestep context 670 671 Output Parameter: 672 . newdt - the timestep to use for the next step 673 674 Level: advanced 675 676 .keywords: timestep, pseudo, default 677 678 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 679 @*/ 680 PetscErrorCode TSPseudoDefaultTimeStep(TS ts,PetscReal* newdt,void* dtctx) 681 { 682 TS_Pseudo *pseudo = (TS_Pseudo*)ts->data; 683 PetscReal inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 684 PetscErrorCode ierr; 685 686 PetscFunctionBegin; 687 ierr = VecZeroEntries(pseudo->xdot);CHKERRQ(ierr); 688 ierr = TSComputeIFunction(ts,ts->ptime,ts->vec_sol,pseudo->xdot,pseudo->func);CHKERRQ(ierr); 689 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm);CHKERRQ(ierr); 690 if (pseudo->initial_fnorm == 0.0) { 691 /* first time through so compute initial function norm */ 692 pseudo->initial_fnorm = pseudo->fnorm; 693 fnorm_previous = pseudo->fnorm; 694 } 695 if (pseudo->fnorm == 0.0) { 696 *newdt = 1.e12*inc*ts->time_step; 697 } else if (pseudo->increment_dt_from_initial_dt) { 698 *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 699 } else { 700 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 701 } 702 pseudo->fnorm_previous = pseudo->fnorm; 703 PetscFunctionReturn(0); 704 } 705