Stepper class for theta integration scheme common in SNL thermal/fluids codes.
More...
#include <Rythmos_ThetaStepper_decl.hpp>
|
| typedef Teuchos::ScalarTraits< Scalar >::magnitudeType | ScalarMag |
| |
|
(Note that these are not member functions.)
|
| template<class Scalar > |
| RCP< ThetaStepper< Scalar > > | thetaStepper (const RCP< Thyra::ModelEvaluator< Scalar > > &model, const RCP< Thyra::NonlinearSolverBase< Scalar > > &solver, RCP< Teuchos::ParameterList > ¶meterList) |
| | Nonmember constructor. More...
|
| |
| template<class Scalar > |
| bool | isInitialized (const StepperBase< Scalar > &stepper) |
| |
| template<class Scalar > |
| bool | isInitialized (const StepperBase< Scalar > &stepper) |
| |
|
| void | setSolver (const RCP< Thyra::NonlinearSolverBase< Scalar > > &solver) |
| |
| RCP< Thyra::NonlinearSolverBase< Scalar > > | getNonconstSolver () |
| |
| RCP< const Thyra::NonlinearSolverBase< Scalar > > | getSolver () const |
| |
|
| RCP< const Thyra::VectorSpaceBase< Scalar > > | get_x_space () const |
| |
| void | addPoints (const Array< Scalar > &time_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &x_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &xdot_vec) |
| |
| TimeRange< Scalar > | getTimeRange () const |
| |
| void | getPoints (const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, Array< ScalarMag > *accuracy_vec) const |
| |
| void | getNodes (Array< Scalar > *time_vec) const |
| |
| void | removeNodes (Array< Scalar > &time_vec) |
| |
| int | getOrder () const |
| |
|
| void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
| |
template<class Scalar>
class Rythmos::ThetaStepper< Scalar >
Stepper class for theta integration scheme common in SNL thermal/fluids codes.
Definition at line 79 of file Rythmos_ThetaStepper_decl.hpp.
◆ ScalarMag
◆ ThetaStepper()
◆ isImplicit()
Return if this stepper is an implicit stepper.
The default implemntation returns false and therefore, by default, a stepper is considered to be an excplicit stepper.
Reimplemented from Rythmos::StepperBase< Scalar >.
◆ setInterpolator()
◆ getNonconstInterpolator()
◆ getInterpolator()
◆ unSetInterpolator()
◆ setSolver()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::setSolver |
( |
const RCP< Thyra::NonlinearSolverBase< Scalar > > & |
solver | ) |
|
◆ getNonconstSolver()
◆ getSolver()
◆ supportsCloning()
◆ cloneStepperAlgorithm()
Creates copies of all internal data (including the parameter list) except the model which is assumed to stateless.
If a shallow copy of the model is not appropirate for some reasone, then the client can simply reset the model using returnVal->setModel().
Reimplemented from Rythmos::StepperBase< Scalar >.
◆ setModel()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::setModel |
( |
const RCP< const Thyra::ModelEvaluator< Scalar > > & |
model | ) |
|
|
virtual |
◆ setNonconstModel()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::setNonconstModel |
( |
const RCP< Thyra::ModelEvaluator< Scalar > > & |
model | ) |
|
|
virtual |
◆ getModel()
◆ getNonconstModel()
◆ setInitialCondition()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::setInitialCondition |
( |
const Thyra::ModelEvaluatorBase::InArgs< Scalar > & |
initialCondition | ) |
|
|
virtual |
◆ getInitialCondition()
◆ takeStep()
◆ getStepStatus()
◆ get_x_space()
◆ addPoints()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::addPoints |
( |
const Array< Scalar > & |
time_vec, |
|
|
const Array< RCP< const Thyra::VectorBase< Scalar > > > & |
x_vec, |
|
|
const Array< RCP< const Thyra::VectorBase< Scalar > > > & |
xdot_vec |
|
) |
| |
◆ getTimeRange()
◆ getPoints()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::getPoints |
( |
const Array< Scalar > & |
time_vec, |
|
|
Array< RCP< const Thyra::VectorBase< Scalar > > > * |
x_vec, |
|
|
Array< RCP< const Thyra::VectorBase< Scalar > > > * |
xdot_vec, |
|
|
Array< ScalarMag > * |
accuracy_vec |
|
) |
| const |
◆ getNodes()
◆ removeNodes()
◆ getOrder()
◆ setParameterList()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::setParameterList |
( |
RCP< Teuchos::ParameterList > const & |
paramList | ) |
|
◆ getNonconstParameterList()
◆ unsetParameterList()
◆ getValidParameters()
◆ describe()
template<class Scalar >
| void Rythmos::ThetaStepper< Scalar >::describe |
( |
Teuchos::FancyOStream & |
out, |
|
|
const Teuchos::EVerbosityLevel |
verbLevel |
|
) |
| const |
◆ thetaStepper()
template<class Scalar >
| RCP< ThetaStepper< Scalar > > thetaStepper |
( |
const RCP< Thyra::ModelEvaluator< Scalar > > & |
model, |
|
|
const RCP< Thyra::NonlinearSolverBase< Scalar > > & |
solver, |
|
|
RCP< Teuchos::ParameterList > & |
parameterList |
|
) |
| |
|
related |
The documentation for this class was generated from the following file: