Tempus
Version of the Day
Time Integration
|
Newmark time stepper in acceleration form (a-form). More...
#include <Tempus_StepperNewmarkImplicitAForm_decl.hpp>
Public Member Functions | |
StepperNewmarkImplicitAForm () | |
Default constructor. More... | |
StepperNewmarkImplicitAForm (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null) | |
Constructor. More... | |
virtual Scalar | getW_xDotDot_coeff (const Scalar) const |
Return W_xDotxDot_coeff = d(xDotDot)/d(xDotDot). More... | |
virtual Scalar | getAlpha (const Scalar dt) const |
Return alpha = d(xDot)/d(xDotDot). More... | |
virtual Scalar | getBeta (const Scalar dt) const |
Return beta = d(x)/d(xDotDot). More... | |
void | predictVelocity (Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const |
void | predictDisplacement (Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const |
void | correctVelocity (Thyra::VectorBase< Scalar > &v, const Thyra::VectorBase< Scalar > &vPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const |
void | correctDisplacement (Thyra::VectorBase< Scalar > &d, const Thyra::VectorBase< Scalar > &dPred, const Thyra::VectorBase< Scalar > &a, const Scalar dt) const |
Basic stepper methods | |
virtual void | setModel (const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel) |
virtual void | setObserver (Teuchos::RCP< StepperObserver< Scalar > >=Teuchos::null) |
Set Observer. 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 | 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 | |
Scalar | beta_ |
Scalar | gamma_ |
Teuchos::RCP < Teuchos::FancyOStream > | out_ |
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_ |
Newmark time stepper in acceleration form (a-form).
Here, we implement the Newmark scheme in predictor/corrector form; see equations (34)-(35) in: A. Mota, W. Klug, M. Ortiz, "Finite element simulation of firearm injury to the human cranium", Computational Mechanics 31(1) 115-121 (2003).
Newmark has two parameters: and , both of which need to be in the range . Newmark can be an explicit or implicit method, depending on the value of the parameter. If , the method is explicit. Regardless of whether the method is implicit or explicit, a linear solve is required. This linear solve can be optimized, however, for the explicit case by lumping the mass matrix. This optimization can be invoked by running "Newmark Explicit d-Form" Stepper through the Piro::TempusSolver class.
Newmark is second order accurate if ; otherwise it is first order accurate. Some additional properties about the Newmark scheme can be found here.
The governing equation solved by this stepper is
For the A-form (i.e., solving for the acceleration, ), we have the following implicit ODE
where and .
Algorithm The algorithm for the Newmark implicit A-form with predictors and correctors is
The First-Step-As-Last (FSAL) principle is part of the Newmark implicit A-Form as the acceleration from the previous time step is used for the predictors. The default is to set useFSAL=true, and useFSAL=false will be ignored.
Definition at line 73 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
Tempus::StepperNewmarkImplicitAForm< Scalar >::StepperNewmarkImplicitAForm | ( | ) |
Default constructor.
Definition at line 90 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
Tempus::StepperNewmarkImplicitAForm< Scalar >::StepperNewmarkImplicitAForm | ( | const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > & | appModel, |
Teuchos::RCP< Teuchos::ParameterList > | pList = Teuchos::null |
||
) |
Constructor.
Definition at line 103 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
void Tempus::StepperNewmarkImplicitAForm< Scalar >::correctDisplacement | ( | Thyra::VectorBase< Scalar > & | d, |
const Thyra::VectorBase< Scalar > & | dPred, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 76 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
void Tempus::StepperNewmarkImplicitAForm< Scalar >::correctVelocity | ( | Thyra::VectorBase< Scalar > & | v, |
const Thyra::VectorBase< Scalar > & | vPred, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 62 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
virtual |
Definition at line 428 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
virtual |
Definition at line 417 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
inlinevirtual |
Return alpha = d(xDot)/d(xDotDot).
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 134 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Return beta = d(x)/d(xDotDot).
Implements Tempus::StepperImplicit< Scalar >.
Definition at line 136 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
virtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 551 of file Tempus_StepperNewmarkImplicitAForm_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 405 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperNewmarkImplicitAForm< Scalar >::getNonconstParameterList | ( | ) |
Definition at line 573 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 114 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 119 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 118 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 128 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
Teuchos::RCP< const Teuchos::ParameterList > Tempus::StepperNewmarkImplicitAForm< Scalar >::getValidParameters | ( | ) | const |
Definition at line 532 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
inlinevirtual |
Return W_xDotxDot_coeff = d(xDotDot)/d(xDotDot).
Definition at line 132 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
virtual |
Initialize during construction and after changing input parameters.
Implements Tempus::Stepper< Scalar >.
Definition at line 140 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 121 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 123 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 122 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 126 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
inlinevirtual |
Implements Tempus::Stepper< Scalar >.
Definition at line 125 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
void Tempus::StepperNewmarkImplicitAForm< Scalar >::predictDisplacement | ( | Thyra::VectorBase< Scalar > & | dPred, |
const Thyra::VectorBase< Scalar > & | d, | ||
const Thyra::VectorBase< Scalar > & | v, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 42 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
void Tempus::StepperNewmarkImplicitAForm< Scalar >::predictVelocity | ( | Thyra::VectorBase< Scalar > & | vPred, |
const Thyra::VectorBase< Scalar > & | v, | ||
const Thyra::VectorBase< Scalar > & | a, | ||
const Scalar | dt | ||
) | const |
Definition at line 28 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
virtual |
Set the initial conditions and make them consistent.
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 156 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
virtual |
Reimplemented from Tempus::StepperImplicit< Scalar >.
Definition at line 125 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
inlinevirtual |
Set Observer.
Implements Tempus::Stepper< Scalar >.
Definition at line 97 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
void Tempus::StepperNewmarkImplicitAForm< Scalar >::setParameterList | ( | const Teuchos::RCP< Teuchos::ParameterList > & | pl | ) |
Definition at line 441 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
virtual |
Take the specified timestep, dt, and return true if successful.
Implements Tempus::Stepper< Scalar >.
Definition at line 324 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
Teuchos::RCP< Teuchos::ParameterList > Tempus::StepperNewmarkImplicitAForm< Scalar >::unsetParameterList | ( | ) |
Definition at line 584 of file Tempus_StepperNewmarkImplicitAForm_impl.hpp.
|
private |
Definition at line 176 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
private |
Definition at line 177 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.
|
private |
Definition at line 179 of file Tempus_StepperNewmarkImplicitAForm_decl.hpp.