Tempus
Version of the Day
Time Integration
|
Backward Euler time stepper. More...
#include <Tempus_StepperBackwardEuler_decl.hpp>
Public Member Functions | |
StepperBackwardEuler () | |
Default constructor. More... | |
StepperBackwardEuler (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null) | |
Constructor. More... | |
virtual Scalar | getAlpha (const Scalar dt) const |
Return alpha = d(xDot)/dx. More... | |
virtual Scalar | getBeta (const Scalar) const |
Return beta = d(x)/dx. More... | |
virtual void | computePredictor (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) |
Compute predictor given the supplied stepper. More... | |
Basic stepper methods | |
virtual void | setObserver (Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null) |
Set Observer. More... | |
void | setPredictor (std::string predictorName) |
Set the predictor. More... | |
void | setPredictor (Teuchos::RCP< Teuchos::ParameterList >predPL=Teuchos::null) |
Set the predictor to the supplied Parameter sublist. This adds a new predictor Parameter sublist to the Stepper's ParameterList. If the predictor sublist is null, it tests if the predictor is set in the Stepper's ParameterList. More... | |
virtual void | initialize () |
Initialize during construction and after changing input parameters. More... | |
virtual void | setInitialConditions (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) |
Set the initial conditions and make them consistent. More... | |
virtual void | takeStep (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) |
Take the specified timestep, dt, and return true if successful. More... | |
virtual Teuchos::RCP < Tempus::StepperState< Scalar > > | getDefaultStepperState () |
Get a default (initial) StepperState. More... | |
virtual Scalar | getOrder () const |
virtual Scalar | getOrderMin () const |
virtual Scalar | getOrderMax () const |
virtual bool | isExplicit () const |
virtual bool | isImplicit () const |
virtual bool | isExplicitImplicit () const |
virtual bool | isOneStepMethod () const |
virtual bool | isMultiStepMethod () const |
virtual OrderODE | getOrderODE () const |
ParameterList methods | |
void | setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &pl) |
Teuchos::RCP < Teuchos::ParameterList > | getNonconstParameterList () |
Teuchos::RCP < Teuchos::ParameterList > | unsetParameterList () |
Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const |
Teuchos::RCP < Teuchos::ParameterList > | getDefaultParameters () const |
Overridden from Teuchos::Describable | |
virtual std::string | description () const |
virtual void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
Implementation of StepperOptimizationInterface | |
virtual int | stencilLength () const |
Return the number of solution vectors in the time step stencil. More... | |
virtual void | computeStepResidual (Thyra::VectorBase< Scalar > &residual, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index) const |
Compute time step residual. More... | |
virtual void | computeStepJacobian (Thyra::LinearOpBase< Scalar > &jacobian, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index, const int deriv_index) const |
Compute time step Jacobian. More... | |
virtual void | computeStepParamDeriv (Thyra::LinearOpBase< Scalar > &deriv, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index) const |
Compute time step derivative w.r.t. model parameters. More... | |
virtual void | computeStepSolver (Thyra::LinearOpWithSolveBase< Scalar > &jacobian_solver, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index) const |
Compute time step Jacobian solver. More... | |
Public Member Functions inherited from Tempus::StepperImplicit< Scalar > | |
virtual void | setModel (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel) |
virtual void | setNonConstModel (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &appModel) |
virtual Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () |
virtual Teuchos::RCP< const WrapperModelEvaluator< Scalar > > | getWrapperModel () |
virtual void | setSolver (std::string solverName) |
Set solver via ParameterList solver name. More... | |
virtual void | setSolver (Teuchos::RCP< Teuchos::ParameterList > solverPL=Teuchos::null) |
Set solver via solver ParameterList. More... | |
virtual void | setSolver (Teuchos::RCP< Thyra::NonlinearSolverBase< Scalar > > solver) |
Set solver. More... | |
virtual Teuchos::RCP < Thyra::NonlinearSolverBase < Scalar > > | getSolver () const |
Get solver. More... | |
virtual std::string | getStepperType () const |
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, f(x, xDot, t, p), residual. More... | |
virtual void | setInitialGuess (Teuchos::RCP< const Thyra::VectorBase< Scalar > > initial_guess) |
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 |
virtual bool | getEmbedded () const |
virtual void | setUseFSAL (bool a) |
virtual bool | getUseFSAL () const |
virtual void | setICConsistency (std::string s) |
virtual std::string | getICConsistency () const |
virtual void | setICConsistencyCheck (bool c) |
virtual bool | getICConsistencyCheck () const |
virtual void | setStepperXDot (Teuchos::RCP< Thyra::VectorBase< Scalar > > xDot) |
Set xDot for Stepper storage. More... | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getStepperXDot (Teuchos::RCP< SolutionState< Scalar > > state) |
Get xDot from SolutionState or Stepper storage. More... | |
virtual Teuchos::RCP < Thyra::VectorBase< Scalar > > | getStepperXDotDot (Teuchos::RCP< SolutionState< Scalar > > state) |
Get xDotDot from SolutionState or Stepper storage. More... | |
Public Member Functions inherited from Tempus::Stepper< Scalar > | |
virtual void | modelWarning () const |
void | getValidParametersBasic (Teuchos::RCP< Teuchos::ParameterList > pl) const |
virtual void | createSubSteppers (std::vector< Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > >) |
void | validExplicitODE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model) const |
Validate that the model supports explicit ODE evaluation, f(x,t) [=xdot]. More... | |
void | validSecondOrderExplicitODE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model) const |
Validate that the model supports explicit second order ODE evaluation, f(x,xdot,t) [=xdotdot]. More... | |
void | validImplicitODE_DAE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model) const |
Validate ME supports implicit ODE/DAE evaluation, f(xdot,x,t) [= 0]. More... | |
void | validSecondOrderODE_DAE (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model) const |
Validate ME supports 2nd order implicit ODE/DAE evaluation, f(xdotdot,xdot,x,t) [= 0]. More... | |
Teuchos::RCP < Teuchos::ParameterList > | defaultSolverParameters () const |
Public Member Functions inherited from Tempus::StepperOptimizationInterface< Scalar > | |
StepperOptimizationInterface () | |
virtual | ~StepperOptimizationInterface () |
Private Member Functions | |
void | computeStepResidDerivImpl (const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs, const Teuchos::Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x, const Teuchos::Array< Scalar > &t, const Thyra::VectorBase< Scalar > &p, const int param_index, const int deriv_index=0) const |
Implementation of computeStep*() methods. More... | |
Private Attributes | |
Teuchos::RCP< Stepper< Scalar > > | predictorStepper_ |
Teuchos::RCP < StepperBackwardEulerObserver < Scalar > > | stepperBEObserver_ |
Additional Inherited Members | |
Protected Attributes inherited from Tempus::StepperImplicit< Scalar > | |
Teuchos::RCP < Teuchos::ParameterList > | stepperPL_ |
Teuchos::RCP < WrapperModelEvaluator < Scalar > > | wrapperModel_ |
Teuchos::RCP < Thyra::NonlinearSolverBase < Scalar > > | solver_ |
Teuchos::RCP< const Thyra::VectorBase< Scalar > > | initial_guess_ |
Teuchos::RCP< StepperObserver < Scalar > > | stepperObserver_ |
Teuchos::RCP < Thyra::VectorBase< Scalar > > | stepperXDot_ |
Teuchos::RCP < Thyra::VectorBase< Scalar > > | stepperXDotDot_ |
Backward Euler time stepper.
For the implicit ODE system, , the solution, and , is determined using a solver (e.g., a non-linear solver, like NOX).
Algorithm The single-timestep algorithm for Backward Euler is simply,
The First-Step-As-Last (FSAL) principle is not needed with Backward Euler. The default is to set useFSAL=false, however useFSAL=true will also work but have no affect (i.e., no-op).
Definition at line 36 of file Tempus_StepperBackwardEuler_decl.hpp.
Tempus::StepperBackwardEuler< Scalar >::StepperBackwardEuler | ( | ) |
Default constructor.
Definition at line 26 of file Tempus_StepperBackwardEuler_impl.hpp.
Tempus::StepperBackwardEuler< Scalar >::StepperBackwardEuler | ( | const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | appModel, |
Teuchos::RCP< Teuchos::ParameterList > | pList = Teuchos::null |
||
) |
Constructor.
Definition at line 34 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Compute predictor given the supplied stepper.
Definition at line 234 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Compute time step Jacobian.
deriv_index determines which component of x the derivative should be computed with respect to.
Implements Tempus::StepperOptimizationInterface< Scalar >.
Definition at line 394 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Compute time step derivative w.r.t. model parameters.
Implements Tempus::StepperOptimizationInterface< Scalar >.
Definition at line 411 of file Tempus_StepperBackwardEuler_impl.hpp.
|
private |
Implementation of computeStep*() methods.
Definition at line 444 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Compute time step residual.
Implements Tempus::StepperOptimizationInterface< Scalar >.
Definition at line 379 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Compute time step Jacobian solver.
Derivative is always w.r.t. the most current solution vector
Implements Tempus::StepperOptimizationInterface< Scalar >.
Definition at line 428 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Definition at line 277 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Definition at line 269 of file Tempus_StepperBackwardEuler_impl.hpp.
|
inlinevirtual |
Return alpha = d(xDot)/dx.
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 93 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Return beta = d(x)/dx.
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 95 of file Tempus_StepperBackwardEuler_decl.hpp.
|
virtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 334 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Get a default (initial) StepperState.
Provide a StepperState to the SolutionState. This Stepper does not have any special state data, so just provide the base class StepperState with the Stepper description. This can be checked to ensure that the input StepperState can be used by this Stepper.
Implements Tempus::Stepper< Scalar >.
Definition at line 260 of file Tempus_StepperBackwardEuler_impl.hpp.
Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperBackwardEuler< Scalar >::getNonconstParameterList | ( | ) |
Definition at line 353 of file Tempus_StepperBackwardEuler_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 78 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 80 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 79 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 89 of file Tempus_StepperBackwardEuler_decl.hpp.
Teuchos::RCP< const Teuchos::ParameterList > Tempus::StepperBackwardEuler< Scalar >::getValidParameters | ( | ) | const |
Definition at line 318 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Initialize during construction and after changing input parameters.
Implements Tempus::Stepper< Scalar >.
Definition at line 132 of file Tempus_StepperBackwardEuler_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 82 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 84 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 83 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 87 of file Tempus_StepperBackwardEuler_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 86 of file Tempus_StepperBackwardEuler_decl.hpp.
|
virtual |
Set the initial conditions and make them consistent.
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 147 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Set Observer.
Implements Tempus::Stepper< Scalar >.
Definition at line 111 of file Tempus_StepperBackwardEuler_impl.hpp.
void Tempus::StepperBackwardEuler< Scalar >::setParameterList | ( | const Teuchos::RCP< Teuchos::ParameterList > & | pl | ) |
Definition at line 287 of file Tempus_StepperBackwardEuler_impl.hpp.
void Tempus::StepperBackwardEuler< Scalar >::setPredictor | ( | std::string | predictorName | ) |
Set the predictor.
Set the predictor to a pre-defined predictor in the ParameterList. The predictor is set to predictorName sublist in the Stepper's ParameterList. The predictorName sublist should already be defined in the Stepper's ParameterList. Otherwise it will fail.
Definition at line 56 of file Tempus_StepperBackwardEuler_impl.hpp.
void Tempus::StepperBackwardEuler< Scalar >::setPredictor | ( | Teuchos::RCP< Teuchos::ParameterList > | predPL = Teuchos::null | ) |
Set the predictor to the supplied Parameter sublist. This adds a new predictor Parameter sublist to the Stepper's ParameterList. If the predictor sublist is null, it tests if the predictor is set in the Stepper's ParameterList.
Definition at line 75 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Return the number of solution vectors in the time step stencil.
Implements Tempus::StepperOptimizationInterface< Scalar >.
Definition at line 371 of file Tempus_StepperBackwardEuler_impl.hpp.
|
virtual |
Take the specified timestep, dt, and return true if successful.
Implements Tempus::Stepper< Scalar >.
Definition at line 172 of file Tempus_StepperBackwardEuler_impl.hpp.
Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperBackwardEuler< Scalar >::unsetParameterList | ( | ) |
Definition at line 361 of file Tempus_StepperBackwardEuler_impl.hpp.
|
private |
Definition at line 160 of file Tempus_StepperBackwardEuler_decl.hpp.
|
private |
Definition at line 161 of file Tempus_StepperBackwardEuler_decl.hpp.