|
Tempus Version of the Day
Time Integration
|
Thyra Base interface for implicit time steppers. More...
#include <Tempus_StepperImplicit_decl.hpp>
Overridden from Teuchos::Describable | |
| Teuchos::RCP< WrapperModelEvaluator< Scalar > > | wrapperModel_ |
| Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > | solver_ |
| Teuchos::RCP< const Thyra::VectorBase< Scalar > > | initialGuess_ |
| bool | zeroInitialGuess_ |
| std::string | solverName_ |
| virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override |
| virtual bool | isValidSetup (Teuchos::FancyOStream &out) const override |
| virtual Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const override |
| Teuchos::RCP< Teuchos::ParameterList > | getValidParametersBasicImplicit () const |
| void | setStepperImplicitValues (Teuchos::RCP< Teuchos::ParameterList > pl) |
| Set StepperImplicit member data from the ParameterList. | |
| void | setStepperSolverValues (Teuchos::RCP< Teuchos::ParameterList > pl) |
| Set solver from ParameterList. | |
| void | setSolverName (std::string i) |
| Set the Solver Name. | |
| std::string | getSolverName () const |
| Get the Solver Name. | |
Basic implicit stepper methods | |
| virtual void | setModel (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel) override |
| Set the model. | |
| virtual Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () const override |
| virtual Teuchos::RCP< const WrapperModelEvaluator< Scalar > > | getWrapperModel () |
| virtual void | setDefaultSolver () |
| virtual void | setSolver (Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > solver) override |
| Set solver. | |
| virtual Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > | getSolver () const override |
| Get solver. | |
| virtual void | setInitialConditions (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) override |
| Set the initial conditions and make them consistent. | |
| virtual Scalar | getAlpha (const Scalar dt) const =0 |
| Return alpha = d(xDot)/dx. | |
| virtual Scalar | getBeta (const Scalar dt) const =0 |
| Return beta = d(x)/dx. | |
| const Thyra::SolveStatus< Scalar > | solveImplicitODE (const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xDot, const Scalar time, const Teuchos::RCP< ImplicitODEParameters< Scalar > > &p, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &y=Teuchos::null, const int index=-1) |
| Solve implicit ODE, f(x, xDot, t, p) = 0. | |
| void | evaluateImplicitODE (Teuchos::RCP< Thyra::VectorBase< Scalar > > &f, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xDot, const Scalar time, const Teuchos::RCP< ImplicitODEParameters< Scalar > > &p) |
| Evaluate implicit ODE residual, f(x, xDot, t, p). | |
| virtual void | setInitialGuess (Teuchos::RCP< const Thyra::VectorBase< Scalar > > initialGuess) override |
| Pass initial guess to Newton solver (only relevant for implicit solvers) | |
| virtual void | setZeroInitialGuess (bool zIG) |
| virtual bool | getZeroInitialGuess () const |
| virtual Scalar | getInitTimeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &) const override |
Additional Inherited Members | |
Public Member Functions inherited from Tempus::Stepper< Scalar > | |
| virtual std::string | description () const |
| void | setStepperValues (const Teuchos::RCP< Teuchos::ParameterList > pl) |
| Set Stepper member data from ParameterList. | |
| Teuchos::RCP< Teuchos::ParameterList > | getValidParametersBasic () const |
| Add basic parameters to Steppers ParameterList. | |
| virtual void | initialize () |
| Initialize after construction and changing input parameters. | |
| virtual bool | isInitialized () |
| True if stepper's member data is initialized. | |
| virtual void | checkInitialized () |
| Check initialization, and error out on failure. | |
| virtual void | takeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)=0 |
| Take the specified timestep, dt, and return true if successful. | |
| virtual Teuchos::RCP< Tempus::StepperState< Scalar > > | getDefaultStepperState ()=0 |
| virtual Scalar | getOrder () const =0 |
| virtual Scalar | getOrderMin () const =0 |
| virtual Scalar | getOrderMax () const =0 |
| virtual bool | isExplicit () const =0 |
| virtual bool | isImplicit () const =0 |
| virtual bool | isExplicitImplicit () const =0 |
| virtual bool | isOneStepMethod () const =0 |
| virtual bool | isMultiStepMethod () const =0 |
| void | setStepperName (std::string s) |
| Set the stepper name. | |
| std::string | getStepperName () const |
| Get the stepper name. | |
| std::string | getStepperType () const |
| Get the stepper type. The stepper type is used as an identifier for the stepper, and can only be set by the derived Stepper class. | |
| virtual void | setUseFSAL (bool a) |
| void | setUseFSALTrueOnly (bool a) |
| void | setUseFSALFalseOnly (bool a) |
| bool | getUseFSAL () const |
| void | setICConsistency (std::string s) |
| std::string | getICConsistency () const |
| void | setICConsistencyCheck (bool c) |
| bool | getICConsistencyCheck () const |
| virtual OrderODE | getOrderODE () const =0 |
| virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperX () |
| Get Stepper x. | |
| virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperXDot () |
| Get Stepper xDot. | |
| virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperXDotDot () |
| Get Stepper xDotDot. | |
| virtual Teuchos::RCP< Thyra::VectorBase< Scalar > > | getStepperXDotDot (Teuchos::RCP< SolutionState< Scalar > > state) |
| Get xDotDot from SolutionState or Stepper storage. | |
Protected Member Functions inherited from Tempus::Stepper< Scalar > | |
| virtual void | setStepperX (Teuchos::RCP< Thyra::VectorBase< Scalar > > x) |
| Set x for Stepper storage. | |
| virtual void | setStepperXDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDot) |
| Set xDot for Stepper storage. | |
| virtual void | setStepperXDotDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDotDot) |
| Set x for Stepper storage. | |
| void | setStepperType (std::string s) |
| Set the stepper type. | |
Protected Attributes inherited from Tempus::Stepper< Scalar > | |
| bool | useFSAL_ = false |
| Use First-Same-As-Last (FSAL) principle. | |
| bool | isInitialized_ |
| True if stepper's member data is initialized. | |
Thyra Base interface for implicit time steppers.
For first-order ODEs, we can write the implicit ODE as
![\[
\mathcal{F}(\dot{x}_n,x_n,t_n) = 0
\]](form_207.png)
where 





![\[
\dot{x}_n(x_n) = \frac{x_n - x_{n-1}}{\Delta t}
\]](form_210.png)
Defining the Iteration Matrix
Often we use Newton's method or one of its variations to solve for 
![\[
\left[
\frac{\partial}{\partial x_n}
\left(
\mathcal{F}(\dot{x}_n,x_n,t_n)
\right)
\right] \Delta x_n^\nu = - \mathcal{F}(\dot{x}_n^\nu,x_n^\nu,t_n)
\]](form_211.png)
where 

![\[
\left[
\frac{\partial \dot{x}_n(x_n) }{\partial x_n}
\frac{\partial}{\partial \dot{x}_n}
\left(
\mathcal{F}(\dot{x}_n,x_n,t_n)
\right)
+
\frac{\partial x_n}{\partial x_n}
\frac{\partial}{\partial x_n}
\left(
\mathcal{F}(\dot{x}_n,x_n,t_n)
\right)
\right] \Delta x_n^\nu = - \mathcal{F}(\dot{x}_n^\nu,x_n^\nu,t_n)
\]](form_214.png)
Defining the iteration matrix, 
![\[
W \Delta x_n^\nu = - \mathcal{F}(\dot{x}_n^\nu,x_n^\nu,t_n)
\]](form_215.png)
using 
![\[
W = \alpha \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n}
+ \beta \frac{\partial \mathcal{F}_n}{\partial x_n}
\]](form_217.png)
and
![\[
W = \alpha M + \beta J
\]](form_218.png)
where
![\[
\alpha \equiv \frac{\partial \dot{x}_n(x_n) }{\partial x_n},
\quad \quad
\beta \equiv \frac{\partial x_n}{\partial x_n} = 1,
\quad \quad
M = \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n},
\quad \quad
J = \frac{\partial \mathcal{F}_n}{\partial x_n}
\]](form_219.png)
and 

Note that sometimes it is helpful to set 





As a concrete example, the time derivative for Backward Euler is
![\[
\dot{x}_n(x_n) = \frac{x_n - x_{n-1}}{\Delta t}
\]](form_210.png)
thus
![\[
\alpha \equiv \frac{\partial \dot{x}_n(x_n) }{\partial x_n}
= \frac{1}{\Delta t}
\quad \quad
\beta \equiv \frac{\partial x_n}{\partial x_n} = 1
\]](form_225.png)
and the iteration matrix for Backward Euler is
![\[
W = \frac{1}{\Delta t} \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n}
+ \frac{\partial \mathcal{F}_n}{\partial x_n}
\]](form_226.png)
Dangers of multiplying through by 

Starting with a simple implicit ODE and multiplying through by 
![\[
\mathcal{F}_n = \Delta t \dot{x}_n - \Delta t \bar{f}(x_n,t_n) = 0
\]](form_227.png)
For the Backward Euler stepper, we recall from above that
![\[
\dot{x}_n(x_n) = \frac{x_n - x_{n-1}}{\Delta t}
\quad\quad
\alpha \equiv \frac{\partial \dot{x}_n(x_n) }{\partial x_n}
= \frac{1}{\Delta t}
\quad \quad
\beta \equiv \frac{\partial x_n}{\partial x_n} = 1
\]](form_228.png)
and we can find for our simple implicit ODE, 
![\[
M = \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n} = \Delta t,
\quad \quad
J = \frac{\partial \mathcal{F}_n}{\partial x_n}
= -\Delta t \frac{\partial \bar{f}_n}{\partial x_n}
\]](form_230.png)
Thus this iteration matrix, 
![\[
W^\ast = \alpha \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n}
+ \beta \frac{\partial \mathcal{F}_n}{\partial x_n}
= \frac{1}{\Delta t} \Delta t
+ (1) \left( - \Delta t \frac{\partial \bar{f}_n}{\partial x_n}
\right)
\]](form_232.png)
or simply
![\[
W^\ast = 1 - \Delta t \frac{\partial \bar{f}_n}{\partial x_n}
\]](form_233.png)
Note that 










Dangers of explicitly including time-derivative definition. If we explicitly include the time-derivative defintion for Backward Euler, we find for our simple implicit ODE,
![\[
\mathcal{F}_n = \frac{x_n - x_{n-1}}{\Delta t} - \bar{f}(x_n,t_n) = 0
\]](form_240.png)
that the iteration matrix is
![\[
W^{\ast\ast} =
\alpha \frac{\partial \mathcal{F}_n}{\partial \dot{x}_n}
+ \beta \frac{\partial \mathcal{F}_n}{\partial x_n}
= \frac{1}{\Delta t} (0)
+ (1) \left(\frac{1}{\Delta t}
- \frac{\partial \bar{f}_n}{\partial x_n}
\right)
\]](form_241.png)
or simply
![\[
W^{\ast\ast} =
\frac{1}{\Delta t} - \frac{\partial \bar{f}_n}{\partial x_n}
\]](form_242.png)
which is the same as 





Definition at line 197 of file Tempus_StepperImplicit_decl.hpp.
|
overridevirtual |
Set the model.
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperBDF2< Scalar >, Tempus::StepperHHTAlpha< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperDIRK< Scalar >.
Definition at line 19 of file Tempus_StepperImplicit_impl.hpp.
|
inlineoverridevirtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperIMEX_RK< Scalar >, and Tempus::StepperIMEX_RK_Partition< Scalar >.
Definition at line 206 of file Tempus_StepperImplicit_decl.hpp.
|
inlinevirtual |
Definition at line 214 of file Tempus_StepperImplicit_decl.hpp.
|
virtual |
Definition at line 224 of file Tempus_StepperImplicit_impl.hpp.
|
overridevirtual |
Set solver.
Reimplemented from Tempus::Stepper< Scalar >.
Definition at line 238 of file Tempus_StepperImplicit_impl.hpp.
|
inlineoverridevirtual |
Get solver.
Reimplemented from Tempus::Stepper< Scalar >.
Definition at line 225 of file Tempus_StepperImplicit_decl.hpp.
|
overridevirtual |
Set the initial conditions and make them consistent.
Implements Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperHHTAlpha< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperDIRK< Scalar >.
Definition at line 34 of file Tempus_StepperImplicit_impl.hpp.
|
pure virtual |
Return alpha = d(xDot)/dx.
Implemented in Tempus::StepperBDF2< Scalar >, Tempus::StepperHHTAlpha< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperDIRK< Scalar >.
|
pure virtual |
Return beta = d(x)/dx.
Implemented in Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperHHTAlpha< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperDIRK< Scalar >.
| const Thyra::SolveStatus< Scalar > Tempus::StepperImplicit< Scalar >::solveImplicitODE | ( | const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | x, |
| const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | xDot, | ||
| const Scalar | time, | ||
| const Teuchos::RCP< ImplicitODEParameters< Scalar > > & | p, | ||
| const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | y = Teuchos::null, |
||
| const int | index = -1 |
||
| ) |
Solve implicit ODE, f(x, xDot, t, p) = 0.
Definition at line 252 of file Tempus_StepperImplicit_impl.hpp.
| void Tempus::StepperImplicit< Scalar >::evaluateImplicitODE | ( | Teuchos::RCP< Thyra::VectorBase< Scalar > > & | f, |
| const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | x, | ||
| const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | xDot, | ||
| const Scalar | time, | ||
| const Teuchos::RCP< ImplicitODEParameters< Scalar > > & | p | ||
| ) |
Evaluate implicit ODE residual, f(x, xDot, t, p).
Definition at line 282 of file Tempus_StepperImplicit_impl.hpp.
|
inlineoverridevirtual |
Pass initial guess to Newton solver (only relevant for implicit solvers)
Implements Tempus::Stepper< Scalar >.
Definition at line 257 of file Tempus_StepperImplicit_decl.hpp.
|
inlinevirtual |
Set parameter so that the initial guess is set to zero (=True) or use last timestep (=False).
Definition at line 266 of file Tempus_StepperImplicit_decl.hpp.
|
inlinevirtual |
Definition at line 272 of file Tempus_StepperImplicit_decl.hpp.
|
inlineoverridevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 274 of file Tempus_StepperImplicit_decl.hpp.
|
overridevirtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperBDF2< Scalar >, Tempus::StepperHHTAlpha< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperDIRK< Scalar >.
Definition at line 307 of file Tempus_StepperImplicit_impl.hpp.
|
overridevirtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperBDF2< Scalar >, Tempus::StepperHHTAlpha< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperDIRK< Scalar >.
Definition at line 320 of file Tempus_StepperImplicit_impl.hpp.
|
overridevirtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperBDF2< Scalar >, Tempus::StepperHHTAlpha< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperSDIRK_2Stage2ndOrder< Scalar >, Tempus::StepperSDIRK_2Stage3rdOrder< Scalar >, Tempus::StepperDIRK_1StageTheta< Scalar >, Tempus::StepperEDIRK_2StageTheta< Scalar >, Tempus::StepperDIRK_General< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperDIRK< Scalar >.
Definition at line 352 of file Tempus_StepperImplicit_impl.hpp.
| Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperImplicit< Scalar >::getValidParametersBasicImplicit | ( | ) | const |
Definition at line 359 of file Tempus_StepperImplicit_impl.hpp.
| void Tempus::StepperImplicit< Scalar >::setStepperImplicitValues | ( | Teuchos::RCP< Teuchos::ParameterList > | pl | ) |
Set StepperImplicit member data from the ParameterList.
Definition at line 373 of file Tempus_StepperImplicit_impl.hpp.
| void Tempus::StepperImplicit< Scalar >::setStepperSolverValues | ( | Teuchos::RCP< Teuchos::ParameterList > | pl | ) |
Set solver from ParameterList.
Definition at line 386 of file Tempus_StepperImplicit_impl.hpp.
|
inline |
Set the Solver Name.
Definition at line 303 of file Tempus_StepperImplicit_decl.hpp.
|
inline |
Get the Solver Name.
Definition at line 305 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 308 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 309 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 310 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 311 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 312 of file Tempus_StepperImplicit_decl.hpp.