petsc-3.10.5 2019-03-28
Report Typos and Errors

TSPSEUDO

Solve steady state ODE and DAE problems with pseudo time stepping This method solves equations of the form

   F(X,Xdot) = 0

for steady state using the iteration

   [G'] S = -F(X,0)
   X += S

where

   G(Y) = F(Y,(Y-X)/dt)

This is linearly-implicit Euler with the residual always evaluated "at steady state". See note below.

Options database keys

-ts_pseudo_increment <real> - ratio of increase dt
-ts_pseudo_increment_dt_from_initial_dt <truth> - Increase dt as a ratio from original dt
-ts_pseudo_fatol <atol> - stop iterating when the function norm is less than atol
-ts_pseudo_frtol <rtol> - stop iterating when the function norm divided by the initial function norm is less than rtol

References

1. - Todd S. Coffey and C. T. Kelley and David E. Keyes, Pseudotransient Continuation and Differential Algebraic Equations, 2003.
2. - C. T. Kelley and David E. Keyes, Convergence analysis of Pseudotransient Continuation, 1998.

Notes

The residual computed by this method includes the transient term (Xdot is computed instead of always being zero), but since the prediction from the last step is always the solution from the last step, on the first Newton iteration we have

 Xdot = (Xpredicted - Xold)/dt = (Xold-Xold)/dt = 0

Therefore, the linear system solved by the first Newton iteration is equivalent to the one described above and in the papers. If the user chooses to perform multiple Newton iterations, the algorithm is no longer the one described in the referenced papers.

See Also

TSCreate(), TS, TSSetType()

Level

beginner

Location

src/ts/impls/pseudo/posindep.c

Examples

src/ts/examples/tutorials/ex1.c.html
src/ts/examples/tutorials/ex24.c.html
src/ts/examples/tutorials/ex26.c.html
src/ts/examples/tutorials/ex42.c.html
src/ts/examples/tutorials/ex1f.F.html

Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages