Tempus
Version of the Day
Time Integration
|
Thyra Base interface for implicit time steppers. More...
#include <Tempus_StepperImplicit_decl.hpp>
Public Member Functions | |
virtual bool | isValidSetup (Teuchos::FancyOStream &out) const |
virtual Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const |
Teuchos::RCP < Teuchos::ParameterList > | getValidParametersBasicImplicit () const |
void | setStepperImplicitValues (Teuchos::RCP< Teuchos::ParameterList > pl) |
Set StepperImplicit member data from the ParameterList. More... | |
void | setStepperSolverValues (Teuchos::RCP< Teuchos::ParameterList > pl) |
Set solver from ParameterList. More... | |
void | setSolverName (std::string i) |
Set the Solver Name. More... | |
std::string | getSolverName () const |
Get the Solver Name. More... | |
Public Member Functions inherited from Tempus::Stepper< Scalar > | |
void | setStepperValues (const Teuchos::RCP< Teuchos::ParameterList > pl) |
Set Stepper member data from ParameterList. More... | |
Teuchos::RCP < Teuchos::ParameterList > | getValidParametersBasic () const |
Add basic parameters to Steppers ParameterList. More... | |
virtual void | setNonConstModel (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &) |
virtual void | initialize () |
Initialize after construction and changing input parameters. More... | |
virtual bool | isInitialized () |
True if stepper's member data is initialized. More... | |
virtual void | checkInitialized () |
Check initialization, and error out on failure. More... | |
virtual void | takeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)=0 |
Take the specified timestep, dt, and return true if successful. More... | |
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. More... | |
std::string | getStepperName () const |
Get the stepper name. More... | |
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. More... | |
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. More... | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getStepperXDot () |
Get Stepper xDot. More... | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getStepperXDotDot () |
Get Stepper xDotDot. More... | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getStepperXDotDot (Teuchos::RCP< SolutionState< Scalar > > state) |
Get xDotDot from SolutionState or Stepper storage. More... | |
virtual std::string | description () const |
Public Member Functions inherited from Teuchos::Describable | |
void | describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
virtual | ~Describable () |
LabeledObject () | |
virtual | ~LabeledObject () |
virtual void | setObjectLabel (const std::string &objectLabel) |
virtual std::string | getObjectLabel () const |
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
Public Member Functions inherited from Teuchos::VerboseObject< Stepper< Scalar > > | |
VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObject & | setVerbLevel (const EVerbosityLevel verbLevel) const |
virtual const VerboseObject & | setOverridingVerbLevel (const EVerbosityLevel verbLevel) const |
virtual EVerbosityLevel | getVerbLevel () const |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT RCP< const ParameterList > | getValidVerboseObjectSublist () |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | setupVerboseObjectSublist (ParameterList *paramList) |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel) |
void | readVerboseObjectSublist (ParameterList *paramList, VerboseObject< ObjectType > *verboseObject) |
Public Member Functions inherited from Teuchos::VerboseObjectBase | |
virtual | ~VerboseObjectBase () |
VerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObjectBase & | setOStream (const RCP< FancyOStream > &oStream) const |
virtual const VerboseObjectBase & | setOverridingOStream (const RCP< FancyOStream > &oStream) const |
virtual VerboseObjectBase & | setLinePrefix (const std::string &linePrefix) |
virtual RCP< FancyOStream > | getOStream () const |
virtual RCP< FancyOStream > | getOverridingOStream () const |
virtual std::string | getLinePrefix () const |
virtual OSTab | getOSTab (const int tabs=1, const std::string &linePrefix="") const |
Protected Attributes | |
Teuchos::RCP < WrapperModelEvaluator < Scalar > > | wrapperModel_ |
Teuchos::RCP < Thyra::NonlinearSolverBase < Scalar > > | solver_ |
Teuchos::RCP< const Thyra::VectorBase< Scalar > > | initialGuess_ |
bool | zeroInitialGuess_ |
std::string | solverName_ |
Protected Attributes inherited from Tempus::Stepper< Scalar > | |
bool | useFSAL_ = false |
Use First-Same-As-Last (FSAL) principle. More... | |
bool | isInitialized_ = false |
True if stepper's member data is initialized. More... | |
Basic implicit stepper methods | |
virtual void | setModel (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel) |
virtual Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () |
virtual Teuchos::RCP< const WrapperModelEvaluator< Scalar > > | getWrapperModel () |
virtual void | setDefaultSolver () |
virtual void | setSolver (Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > solver) |
Set solver. More... | |
virtual Teuchos::RCP < Thyra::NonlinearSolverBase < Scalar > > | getSolver () const |
Get solver. More... | |
virtual void | setInitialConditions (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) |
Set the initial conditions and make them consistent. More... | |
virtual Scalar | getAlpha (const Scalar dt) const =0 |
Return alpha = d(xDot)/dx. More... | |
virtual Scalar | getBeta (const Scalar dt) const =0 |
Return beta = d(x)/dx. More... | |
const Thyra::SolveStatus< Scalar > | solveImplicitODE (const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x) |
Solve problem using x in-place. (Needs to be deprecated!) More... | |
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) |
Solve implicit ODE, f(x, xDot, t, p) = 0. More... | |
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). More... | |
virtual void | setInitialGuess (Teuchos::RCP< const Thyra::VectorBase< Scalar > > initialGuess) |
Pass initial guess to Newton solver (only relevant for implicit solvers) More... | |
virtual void | setZeroInitialGuess (bool zIG) |
Set parameter so that the initial guess is set to zero (=True) or use last timestep (=False). More... | |
virtual bool | getZeroInitialGuess () const |
virtual Scalar | getInitTimeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &) const |
Overridden from Teuchos::Describable | |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
Additional Inherited Members | |
Static Public Member Functions inherited from Teuchos::VerboseObject< Stepper< Scalar > > | |
static void | setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel) |
static EVerbosityLevel | getDefaultVerbLevel () |
Static Public Member Functions inherited from Teuchos::VerboseObjectBase | |
static void | setDefaultOStream (const RCP< FancyOStream > &defaultOStream) |
static RCP< FancyOStream > | getDefaultOStream () |
Static Public Attributes inherited from Teuchos::Describable | |
static const EVerbosityLevel | verbLevel_default |
Protected Member Functions inherited from Tempus::Stepper< Scalar > | |
virtual void | setStepperX (Teuchos::RCP< Thyra::VectorBase< Scalar > > x) |
Set x for Stepper storage. More... | |
virtual void | setStepperXDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDot) |
Set xDot for Stepper storage. More... | |
virtual void | setStepperXDotDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDotDot) |
Set x for Stepper storage. More... | |
void | setStepperType (std::string s) |
Set the stepper type. More... | |
Protected Member Functions inherited from Teuchos::VerboseObject< Stepper< Scalar > > | |
void | initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) |
Protected Member Functions inherited from Teuchos::VerboseObjectBase | |
void | initializeVerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) |
virtual void | informUpdatedVerbosityState () const |
Thyra Base interface for implicit time steppers.
For first-order ODEs, we can write the implicit ODE as
where is the solution vector, is the time derivative, is the time and indicates the time level. Note that is different for each time stepper and is a function of other solution states, e.g., for Backward Euler,
Defining the Iteration Matrix
Often we use Newton's method or one of its variations to solve for , such as
where and is the iteration index. Using the chain rule for a function with multiple variables, we can write
Defining the iteration matrix, , we have
using , where
and
where
and is the mass matrix and is the Jacobian.
Note that sometimes it is helpful to set and to obtain the Jacobian, , from the iteration matrix (i.e., ModelEvaluator), or set and to obtain the mass matrix, , from the iteration matrix (i.e., the ModelEvaluator).
As a concrete example, the time derivative for Backward Euler is
thus
and the iteration matrix for Backward Euler is
Dangers of multiplying through by . In some time-integration schemes, the application might want to multiply the governing equations by the time-step size, , for scaling or other reasons. Here we illustrate what that means and the complications that follow from it.
Starting with a simple implicit ODE and multiplying through by , we have
For the Backward Euler stepper, we recall from above that
and we can find for our simple implicit ODE, ,
Thus this iteration matrix, , is
or simply
Note that is not the same as from above (i.e., ). But we should not infer from this is that or , as those definitions are unchanged (i.e., and ). However, the mass matrix, , the Jacobian, and the residual, , all need to include in their evaluations (i.e., be included in the ModelEvaluator return values for these terms).
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,
that the iteration matrix is
or simply
which is the same as from above for Backward Euler, but again we should not infer that or . However the major drawback is the mass matrix, , the Jacobian, , and the residual, (i.e., the ModelEvaluator) are explicitly written only for Backward Euler. The application would need to write separate ModelEvaluators for each stepper, thus destroying the ability to re-use the ModelEvaluator with any stepper.
Definition at line 228 of file Tempus_StepperImplicit_decl.hpp.
|
virtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, and Tempus::StepperHHTAlpha< Scalar >.
Definition at line 20 of file Tempus_StepperImplicit_impl.hpp.
|
inlinevirtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperIMEX_RK_Partition< Scalar >, and Tempus::StepperIMEX_RK< Scalar >.
Definition at line 237 of file Tempus_StepperImplicit_decl.hpp.
|
inlinevirtual |
Definition at line 245 of file Tempus_StepperImplicit_decl.hpp.
|
virtual |
Definition at line 212 of file Tempus_StepperImplicit_impl.hpp.
|
virtual |
Set solver.
Reimplemented from Tempus::Stepper< Scalar >.
Definition at line 227 of file Tempus_StepperImplicit_impl.hpp.
|
inlinevirtual |
Get solver.
Reimplemented from Tempus::Stepper< Scalar >.
Definition at line 253 of file Tempus_StepperImplicit_decl.hpp.
|
virtual |
Set the initial conditions and make them consistent.
Implements Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperDIRK< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, and Tempus::StepperHHTAlpha< Scalar >.
Definition at line 36 of file Tempus_StepperImplicit_impl.hpp.
|
pure virtual |
Return alpha = d(xDot)/dx.
Implemented in Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperDIRK< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperHHTAlpha< Scalar >.
|
pure virtual |
Return beta = d(x)/dx.
Implemented in Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperDIRK< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperTrapezoidal< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperHHTAlpha< Scalar >.
const Thyra::SolveStatus< Scalar > Tempus::StepperImplicit< Scalar >::solveImplicitODE | ( | const Teuchos::RCP< Thyra::VectorBase< Scalar > > & | x | ) |
Solve problem using x in-place. (Needs to be deprecated!)
Definition at line 242 of file Tempus_StepperImplicit_impl.hpp.
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 | ||
) |
Solve implicit ODE, f(x, xDot, t, p) = 0.
Definition at line 256 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 302 of file Tempus_StepperImplicit_impl.hpp.
|
inlinevirtual |
Pass initial guess to Newton solver (only relevant for implicit solvers)
Implements Tempus::Stepper< Scalar >.
Definition at line 286 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 294 of file Tempus_StepperImplicit_decl.hpp.
|
inlinevirtual |
Definition at line 299 of file Tempus_StepperImplicit_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 301 of file Tempus_StepperImplicit_decl.hpp.
|
virtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperDIRK< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperBackwardEuler< Scalar >, Tempus::StepperTrapezoidal< Scalar >, and Tempus::StepperHHTAlpha< Scalar >.
Definition at line 328 of file Tempus_StepperImplicit_impl.hpp.
|
virtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperDIRK< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperBackwardEuler< Scalar >, Tempus::StepperTrapezoidal< Scalar >, and Tempus::StepperHHTAlpha< Scalar >.
Definition at line 341 of file Tempus_StepperImplicit_impl.hpp.
|
virtual |
Reimplemented from Tempus::Stepper< Scalar >.
Reimplemented in Tempus::StepperDIRK_General< Scalar >, Tempus::StepperEDIRK_2StageTheta< Scalar >, Tempus::StepperDIRK_1StageTheta< Scalar >, Tempus::StepperSDIRK_2Stage3rdOrder< Scalar >, Tempus::StepperSDIRK_2Stage2ndOrder< Scalar >, Tempus::StepperIMEX_RK_Partition< Scalar >, Tempus::StepperIMEX_RK< Scalar >, Tempus::StepperDIRK< Scalar >, Tempus::StepperBDF2< Scalar >, Tempus::StepperNewmarkImplicitAForm< Scalar >, Tempus::StepperNewmarkImplicitDForm< Scalar >, Tempus::StepperBackwardEuler< Scalar >, and Tempus::StepperHHTAlpha< Scalar >.
Definition at line 373 of file Tempus_StepperImplicit_impl.hpp.
Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperImplicit< Scalar >::getValidParametersBasicImplicit | ( | ) | const |
Definition at line 381 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 397 of file Tempus_StepperImplicit_impl.hpp.
void Tempus::StepperImplicit< Scalar >::setStepperSolverValues | ( | Teuchos::RCP< Teuchos::ParameterList > | pl | ) |
Set solver from ParameterList.
Definition at line 412 of file Tempus_StepperImplicit_impl.hpp.
|
inline |
Set the Solver Name.
Definition at line 325 of file Tempus_StepperImplicit_decl.hpp.
|
inline |
Get the Solver Name.
Definition at line 327 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 331 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 332 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 333 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 334 of file Tempus_StepperImplicit_decl.hpp.
|
protected |
Definition at line 335 of file Tempus_StepperImplicit_decl.hpp.