Tempus
Version of the Day
Time Integration
|
BDF2 (Backward-Difference-Formula-2) time stepper. More...
#include <Tempus_StepperBDF2_decl.hpp>
Public Member Functions | |
StepperBDF2 () | |
Default constructor. More... | |
StepperBDF2 (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 | getAlpha (const Scalar dt, const Scalar dtOld) const |
virtual Scalar | getBeta (const Scalar) const |
Return beta = d(x)/dx. More... | |
virtual void | computeStartUp (const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) |
Compute the first time step given the supplied startup stepper. More... | |
Basic stepper methods | |
virtual void | setObserver (Teuchos::RCP< StepperObserver< Scalar > > obs=Teuchos::null) |
Set Observer. More... | |
void | setStartUpStepper (std::string startupStepperName) |
Set the stepper to use in first step. More... | |
void | setStartUpStepper (Teuchos::RCP< Teuchos::ParameterList >startUpStepperPL=Teuchos::null) |
Set the start up stepper to the supplied Parameter sublist. 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 |
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 |
Private Attributes | |
Teuchos::RCP< Stepper< Scalar > > | startUpStepper_ |
Teuchos::RCP < StepperBDF2Observer< Scalar > > | stepperBDF2Observer_ |
Scalar | order_ |
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_ |
BDF2 (Backward-Difference-Formula-2) time stepper.
For the implicit ODE system, , the solution, and , is determined using a solver (e.g., a non-linear solver, like NOX). This stepper allows for a variable time-step, . It is a 2-step method.
Algorithm
The First-Step-As-Last (FSAL) principle is not needed BDF2. The default is to set useFSAL=false, however useFSAL=true will also work but have no affect (i.e., no-op).
Definition at line 50 of file Tempus_StepperBDF2_decl.hpp.
Tempus::StepperBDF2< Scalar >::StepperBDF2 | ( | ) |
Default constructor.
Definition at line 26 of file Tempus_StepperBDF2_impl.hpp.
Tempus::StepperBDF2< Scalar >::StepperBDF2 | ( | const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | appModel, |
Teuchos::RCP< Teuchos::ParameterList > | pList = Teuchos::null |
||
) |
Constructor.
Definition at line 34 of file Tempus_StepperBDF2_impl.hpp.
|
virtual |
Compute the first time step given the supplied startup stepper.
Definition at line 265 of file Tempus_StepperBDF2_impl.hpp.
|
virtual |
Definition at line 305 of file Tempus_StepperBDF2_impl.hpp.
|
virtual |
Definition at line 297 of file Tempus_StepperBDF2_impl.hpp.
|
inlinevirtual |
Return alpha = d(xDot)/dx.
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 105 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Definition at line 106 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Return beta = d(x)/dx.
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 109 of file Tempus_StepperBDF2_decl.hpp.
|
virtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 362 of file Tempus_StepperBDF2_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 288 of file Tempus_StepperBDF2_impl.hpp.
Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperBDF2< Scalar >::getNonconstParameterList | ( | ) |
Definition at line 381 of file Tempus_StepperBDF2_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 90 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 92 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 91 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 101 of file Tempus_StepperBDF2_decl.hpp.
Teuchos::RCP< const Teuchos::ParameterList > Tempus::StepperBDF2< Scalar >::getValidParameters | ( | ) | const |
Definition at line 345 of file Tempus_StepperBDF2_impl.hpp.
|
virtual |
Initialize during construction and after changing input parameters.
Implements Tempus::Stepper< Scalar >.
Definition at line 152 of file Tempus_StepperBDF2_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 94 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 96 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 95 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 99 of file Tempus_StepperBDF2_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 98 of file Tempus_StepperBDF2_decl.hpp.
|
virtual |
Set the initial conditions and make them consistent.
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 168 of file Tempus_StepperBDF2_impl.hpp.
|
virtual |
Set Observer.
Implements Tempus::Stepper< Scalar >.
Definition at line 131 of file Tempus_StepperBDF2_impl.hpp.
void Tempus::StepperBDF2< Scalar >::setParameterList | ( | const Teuchos::RCP< Teuchos::ParameterList > & | pl | ) |
Definition at line 315 of file Tempus_StepperBDF2_impl.hpp.
void Tempus::StepperBDF2< Scalar >::setStartUpStepper | ( | std::string | startupStepperName | ) |
Set the stepper to use in first step.
Set the startup stepper to a pre-defined stepper in the ParameterList.
The startup stepper is set to stepperName sublist in the Stepper's ParameterList. The stepperName sublist should already be defined in the Stepper's ParameterList. Otherwise it will fail.
Definition at line 57 of file Tempus_StepperBDF2_impl.hpp.
void Tempus::StepperBDF2< Scalar >::setStartUpStepper | ( | Teuchos::RCP< Teuchos::ParameterList > | startupStepperPL = Teuchos::null | ) |
Set the start up stepper to the supplied Parameter sublist.
This adds a new start up stepper Parameter sublist to the Stepper's ParameterList. If the start up stepper sublist is null, it tests if the stepper sublist is set in the Stepper's ParameterList.
Definition at line 78 of file Tempus_StepperBDF2_impl.hpp.
|
virtual |
Take the specified timestep, dt, and return true if successful.
Implements Tempus::Stepper< Scalar >.
Definition at line 193 of file Tempus_StepperBDF2_impl.hpp.
Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperBDF2< Scalar >::unsetParameterList | ( | ) |
Definition at line 389 of file Tempus_StepperBDF2_impl.hpp.
|
private |
Definition at line 135 of file Tempus_StepperBDF2_decl.hpp.
|
private |
Definition at line 133 of file Tempus_StepperBDF2_decl.hpp.
|
private |
Definition at line 134 of file Tempus_StepperBDF2_decl.hpp.