Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
Tempus_Test::SinCosModel< Scalar > Class Template Reference

Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation

\[ \mathbf{\ddot{x}}=-\mathbf{x} \]

with a few enhancements. We start with the exact solution to the differential equation

\begin{eqnarray*} x_{0}(t) & = & a+b*\sin((f/L)*t+\phi)\\ x_{1}(t) & = & b*(f/L)*\cos((f/L)*t+\phi) \end{eqnarray*}

then the form of the model is

\begin{eqnarray*} \frac{d}{dt}x_{0}(t) & = & x_{1}(t)\\ \frac{d}{dt}x_{1}(t) & = & \left(\frac{f}{L}\right)^{2}(a-x_{0}(t)) \end{eqnarray*}

where the default parameter values are $a=0$, $f=1$, and $L=1$, and the initial conditions

\begin{eqnarray*} x_{0}(t_{0}=0) & = & \gamma_{0}[=0]\\ x_{1}(t_{0}=0) & = & \gamma_{1}[=1] \end{eqnarray*}

determine the remaining coefficients

\begin{eqnarray*} \phi & = & \arctan(((f/L)/\gamma_{1})*(\gamma_{0}-a))-(f/L)*t_{0}[=0]\\ b & = & \gamma_{1}/((f/L)*cos((f/L)*t_{0}+\phi))[=1] \end{eqnarray*}

. More...

#include <SinCosModel_decl.hpp>

Inheritance diagram for Tempus_Test::SinCosModel< Scalar >:
Thyra::StateFuncModelEvaluatorBase< Scalar > Teuchos::ParameterListAcceptorDefaultBase Thyra::ModelEvaluatorDefaultBase< Scalar > Teuchos::ParameterListAcceptor Thyra::ModelEvaluator< class > Thyra::ModelEvaluatorBase Teuchos::Describable Teuchos::LabeledObject Tempus_Test::SinCosModelAdjoint< Scalar >

Public Member Functions

 SinCosModel (Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
getExactSolution (double t) const
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
getExactSensSolution (int j, double t) const
 
- Public Member Functions inherited from Thyra::StateFuncModelEvaluatorBase< Scalar >
RCP< const VectorSpaceBase
< Scalar > > 
get_p_space (int l) const
 
RCP< const Teuchos::Array
< std::string > > 
get_p_names (int l) const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_g_space (int j) const
 
Teuchos::ArrayView< const
std::string > 
get_g_names (int j) const
 
ModelEvaluatorBase::InArgs
< Scalar > 
getNominalValues () const
 
ModelEvaluatorBase::InArgs
< Scalar > 
getLowerBounds () const
 
ModelEvaluatorBase::InArgs
< Scalar > 
getUpperBounds () const
 
RCP< LinearOpBase< Scalar > > create_W_op () const
 
RCP< PreconditionerBase< Scalar > > create_W_prec () const
 
RCP< const
LinearOpWithSolveFactoryBase
< Scalar > > 
get_W_factory () const
 
void reportFinalPoint (const ModelEvaluatorBase::InArgs< Scalar > &finalPoint, const bool wasSolved)
 
- Public Member Functions inherited from Thyra::ModelEvaluator< class >
virtual int Np () const =0
 
virtual int Ng () const =0
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
get_f_multiplier_space () const =0
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
get_g_multiplier_space (int j) const =0
 
virtual RCP< LinearOpBase
< Scalar > > 
create_DfDp_op (int l) const =0
 
virtual RCP< LinearOpBase
< Scalar > > 
create_DgDx_dot_op (int j) const =0
 
virtual RCP< LinearOpBase
< Scalar > > 
create_DgDx_op (int j) const =0
 
virtual RCP< LinearOpBase
< Scalar > > 
create_DgDp_op (int j, int l) const =0
 
virtual
ModelEvaluatorBase::OutArgs
< Scalar > 
createOutArgs () const =0
 
virtual void evalModel (const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const =0
 
 ModelEvaluatorBase ()
 
std::string toString (ModelEvaluatorBase::EInArgsMembers)
 
std::string toString (ModelEvaluatorBase::EOutArgsMembers)
 
std::string toString (ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation)
 
ModelEvaluatorBase::EDerivativeMultiVectorOrientation getOtherDerivativeMultiVectorOrientation (ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation)
 
- Public Member Functions inherited from Teuchos::Describable
virtual std::string description () const
 
virtual void describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
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::ParameterListAcceptorDefaultBase
RCP< ParameterListgetNonconstParameterList ()
 
RCP< ParameterListunsetParameterList ()
 
RCP< const ParameterListgetParameterList () const
 

Protected Attributes

int dim_
 Number of state unknowns (2) More...
 
int Np_
 Number of parameter vectors (1) More...
 
int np_
 Number of parameters in this vector (2) More...
 
int Ng_
 Number of observation functions (1) More...
 
int ng_
 Number of elements in this observation function (1) More...
 
bool haveIC_
 false => no nominal values are provided (default=true) More...
 
bool acceptModelParams_
 Changes inArgs to require parameters. More...
 
bool useDfDpAsTangent_
 Treat DfDp OutArg as tangent (df/dx*dx/dp+df/dp) More...
 
bool isInitialized_
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
inArgs_
 
Thyra::ModelEvaluatorBase::OutArgs
< Scalar > 
outArgs_
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
nominalValues_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
x_space_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
f_space_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
p_space_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
g_space_
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
DxDp_space_
 
Scalar a_
 Model parameter. More...
 
Scalar f_
 Model parameter. More...
 
Scalar L_
 Model parameter. More...
 
Scalar phi_
 Parameter determined from the IC. More...
 
Scalar b_
 Parameter determined from the IC. More...
 
Scalar t0_ic_
 Time value where the initial condition is specified. More...
 
Scalar x0_ic_
 Initial condition for x0. More...
 
Scalar x1_ic_
 Initial condition for x1. More...
 

Private Member Functions

void setupInOutArgs_ () const
 
void calculateCoeffFromIC_ ()
 

Public functions overridden from ModelEvaluator.

Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_x_space () const
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_f_space () const
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
getNominalValues () const
 
Teuchos::RCP
< Thyra::LinearOpWithSolveBase
< Scalar > > 
create_W () const
 
Teuchos::RCP
< Thyra::LinearOpBase< Scalar > > 
create_W_op () const
 
Teuchos::RCP< const
Thyra::LinearOpWithSolveFactoryBase
< Scalar > > 
get_W_factory () const
 
Thyra::ModelEvaluatorBase::InArgs
< Scalar > 
createInArgs () const
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_p_space (int l) const
 
Teuchos::RCP< const
Teuchos::Array< std::string > > 
get_p_names (int l) const
 
Teuchos::RCP< const
Thyra::VectorSpaceBase< Scalar > > 
get_g_space (int j) const
 

Public functions overridden from ParameterListAcceptor.

void setParameterList (Teuchos::RCP< Teuchos::ParameterList > const &paramList)
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 

Private functions overridden from ModelEvaluatorDefaultBase.

Thyra::ModelEvaluatorBase::OutArgs
< Scalar > 
createOutArgsImpl () const
 
void evalModelImpl (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs_bar, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs_bar) const
 

Additional Inherited Members

- Public Types inherited from Thyra::ModelEvaluator< class >
enum  EInArgsMembers
 
enum  EInArgs_p_mp
 
enum  EEvalType
 
enum  EDerivativeMultiVectorOrientation
 
enum  EDerivativeLinearOp
 
enum  EDerivativeLinearity
 
enum  ERankStatus
 
enum  EOutArgsMembers
 
enum  EOutArgsDfDp
 
enum  EOutArgsDgDx_dot
 
enum  EOutArgsDgDx
 
enum  EOutArgsDgDp
 
enum  EOutArgsDfDp_mp
 
enum  EOutArgs_g_mp
 
enum  EOutArgsDgDx_dot_mp
 
enum  EOutArgsDgDx_mp
 
enum  EOutArgsDgDp_mp
 
typedef Teuchos::ScalarTraits
< Scalar >::magnitudeType 
ScalarMag
 
- Static Public Attributes inherited from Thyra::ModelEvaluator< class >
static const int NUM_E_IN_ARGS_MEMBERS
 
static const int NUM_E_OUT_ARGS_MEMBERS
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 
- Protected Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
void setMyParamList (const RCP< ParameterList > &paramList)
 
RCP< ParameterListgetMyNonconstParamList ()
 
RCP< const ParameterListgetMyParamList () const
 

Detailed Description

template<class Scalar>
class Tempus_Test::SinCosModel< Scalar >

Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation

\[ \mathbf{\ddot{x}}=-\mathbf{x} \]

with a few enhancements. We start with the exact solution to the differential equation

\begin{eqnarray*} x_{0}(t) & = & a+b*\sin((f/L)*t+\phi)\\ x_{1}(t) & = & b*(f/L)*\cos((f/L)*t+\phi) \end{eqnarray*}

then the form of the model is

\begin{eqnarray*} \frac{d}{dt}x_{0}(t) & = & x_{1}(t)\\ \frac{d}{dt}x_{1}(t) & = & \left(\frac{f}{L}\right)^{2}(a-x_{0}(t)) \end{eqnarray*}

where the default parameter values are $a=0$, $f=1$, and $L=1$, and the initial conditions

\begin{eqnarray*} x_{0}(t_{0}=0) & = & \gamma_{0}[=0]\\ x_{1}(t_{0}=0) & = & \gamma_{1}[=1] \end{eqnarray*}

determine the remaining coefficients

\begin{eqnarray*} \phi & = & \arctan(((f/L)/\gamma_{1})*(\gamma_{0}-a))-(f/L)*t_{0}[=0]\\ b & = & \gamma_{1}/((f/L)*cos((f/L)*t_{0}+\phi))[=1] \end{eqnarray*}

.

Therefore this model has three model parameters and two initial conditions which effect the exact solution as above.

\[ \mathbf{p}=(a,f,L) \]

\[ \dot{\mathbf{x}}=\mathbf{F}(\mathbf{x},t,\mathbf{p}) \]

where

\begin{eqnarray*} F_{0} & = & x_{1}\\ F_{1} & = & \left(\frac{f}{L}\right)^{2}(a-x_{0}) \end{eqnarray*}

The exact sensitivities, $\mathbf{s}=\partial\mathbf{x}/\partial\mathbf{p}$, for the problem are specified as

\[ \mathbf{s}(t)=\left[\begin{array}{cc} 1 & 0\\ \left(\frac{b}{L}\right)t\,\cos\left(\left(\frac{f}{L}\right)t+\phi\right) & \left(\frac{b}{L}\right)\cos\left(\left(\frac{f}{L}\right)t+\phi\right)-\frac{b\, f\, t}{L^{2}}\sin\left(\left(\frac{f}{L}\right)t+\phi\right)\\ -\frac{b\, f\, t}{L^{2}}\cos\left(\left(\frac{f}{L}\right)t+\phi\right) & -\frac{b\, f}{L^{2}}\cos\left(\left(\frac{f}{L}\right)t+\phi\right)+\frac{b\, f^{2}\, t}{L^{3}}\sin\left(\left(\frac{f}{L}\right)t+\phi\right) \end{array}\right] \]

and for the default initial conditions, $\phi=0$ and $b=1$

\[ \mathbf{s}(t=0)=\left[\begin{array}{cc} 1 & 0\\ 0 & \frac{b}{L}\\ 0 & -\frac{f}{L^{2}} \end{array}\right] \]

The time differentiated sensitivities, $\dot{\mathbf{s}}=\partial\mathbf{s}/\partial t=\partial/\partial t(\partial\mathbf{x}/\partial\mathbf{p})=\partial/\partial\mathbf{p}(\partial\mathbf{x}/\partial t)$ are

\[ \dot{\mathbf{s}}(t)=\left[\begin{array}{cc} 0 & 0\\ \left(\frac{b}{L}\right)\cos\left(\left(\frac{f}{L}\right)t+\phi\right)-\frac{b\, f\, t}{L^{2}}\sin\left(\left(\frac{f}{L}\right)t+\phi\right) & -\frac{2b\, f}{L^{2}}\sin\left(\left(\frac{f}{L}\right)t+\phi\right)\left(\frac{b}{L}\right)-\frac{b\, f^{2}\, t}{L^{3}}\cos\left(\left(\frac{f}{L}\right)t+\phi\right)\\ -\frac{b\, f}{L^{2}}\cos\left(\left(\frac{f}{L}\right)t+\phi\right)+\frac{b\, f^{2}\, t}{L^{3}}\sin\left(\left(\frac{f}{L}\right)t+\phi\right) & \frac{2b\, f^{2}}{L^{3}}\sin\left(\left(\frac{f}{L}\right)t+\phi\right)+\frac{b\, f^{3}\, t}{L^{4}}\cos\left(\left(\frac{f}{L}\right)t+\phi\right) \end{array}\right] \]

Definition at line 108 of file SinCosModel_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Tempus_Test::SinCosModel< Scalar >::SinCosModel ( Teuchos::RCP< Teuchos::ParameterList pList = Teuchos::null)

Definition at line 29 of file SinCosModel_impl.hpp.

Member Function Documentation

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > Tempus_Test::SinCosModel< Scalar >::getExactSolution ( double  t) const

Definition at line 61 of file SinCosModel_impl.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > Tempus_Test::SinCosModel< Scalar >::getExactSensSolution ( int  j,
double  t 
) const

Definition at line 97 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > Tempus_Test::SinCosModel< Scalar >::get_x_space ( ) const
virtual

Implements Thyra::ModelEvaluator< class >.

Definition at line 160 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > Tempus_Test::SinCosModel< Scalar >::get_f_space ( ) const
virtual

Implements Thyra::ModelEvaluator< class >.

Definition at line 167 of file SinCosModel_impl.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > Tempus_Test::SinCosModel< Scalar >::getNominalValues ( ) const
virtual

Implements Thyra::ModelEvaluator< class >.

Definition at line 174 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< Thyra::LinearOpWithSolveBase< Scalar > > Tempus_Test::SinCosModel< Scalar >::create_W ( ) const
virtual

Implements Thyra::ModelEvaluator< class >.

Reimplemented in Tempus_Test::SinCosModelAdjoint< Scalar >.

Definition at line 183 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > Tempus_Test::SinCosModel< Scalar >::create_W_op ( ) const
virtual

Implements Thyra::ModelEvaluator< class >.

Reimplemented in Tempus_Test::SinCosModelAdjoint< Scalar >.

Definition at line 226 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > Tempus_Test::SinCosModel< Scalar >::get_W_factory ( ) const
virtual

Implements Thyra::ModelEvaluator< class >.

Definition at line 236 of file SinCosModel_impl.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > Tempus_Test::SinCosModel< Scalar >::createInArgs ( ) const
virtual

Implements Thyra::ModelEvaluator< class >.

Reimplemented in Tempus_Test::SinCosModelAdjoint< Scalar >.

Definition at line 255 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > Tempus_Test::SinCosModel< Scalar >::get_p_space ( int  l) const
virtual

Implements Thyra::ModelEvaluator< class >.

Definition at line 436 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Teuchos::Array< std::string > > Tempus_Test::SinCosModel< Scalar >::get_p_names ( int  l) const
virtual

Implements Thyra::ModelEvaluator< class >.

Definition at line 451 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > Tempus_Test::SinCosModel< Scalar >::get_g_space ( int  j) const
virtual

Implements Thyra::ModelEvaluator< class >.

Definition at line 473 of file SinCosModel_impl.hpp.

template<class Scalar >
void Tempus_Test::SinCosModel< Scalar >::setParameterList ( Teuchos::RCP< Teuchos::ParameterList > const &  paramList)
virtual

Implements Teuchos::ParameterListAcceptor.

Definition at line 555 of file SinCosModel_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Teuchos::ParameterList > Tempus_Test::SinCosModel< Scalar >::getValidParameters ( ) const
virtual

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 587 of file SinCosModel_impl.hpp.

template<class Scalar >
void Tempus_Test::SinCosModel< Scalar >::setupInOutArgs_ ( ) const
private

Definition at line 482 of file SinCosModel_impl.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::OutArgs< Scalar > Tempus_Test::SinCosModel< Scalar >::createOutArgsImpl ( ) const
private

Definition at line 266 of file SinCosModel_impl.hpp.

template<class Scalar >
void Tempus_Test::SinCosModel< Scalar >::evalModelImpl ( const Thyra::ModelEvaluatorBase::InArgs< Scalar > &  inArgs_bar,
const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &  outArgs_bar 
) const
private

Definition at line 273 of file SinCosModel_impl.hpp.

template<class Scalar >
void Tempus_Test::SinCosModel< Scalar >::calculateCoeffFromIC_ ( )
private

Definition at line 611 of file SinCosModel_impl.hpp.

Member Data Documentation

template<class Scalar >
int Tempus_Test::SinCosModel< Scalar >::dim_
protected

Number of state unknowns (2)

Definition at line 159 of file SinCosModel_decl.hpp.

template<class Scalar >
int Tempus_Test::SinCosModel< Scalar >::Np_
protected

Number of parameter vectors (1)

Definition at line 160 of file SinCosModel_decl.hpp.

template<class Scalar >
int Tempus_Test::SinCosModel< Scalar >::np_
protected

Number of parameters in this vector (2)

Definition at line 161 of file SinCosModel_decl.hpp.

template<class Scalar >
int Tempus_Test::SinCosModel< Scalar >::Ng_
protected

Number of observation functions (1)

Definition at line 162 of file SinCosModel_decl.hpp.

template<class Scalar >
int Tempus_Test::SinCosModel< Scalar >::ng_
protected

Number of elements in this observation function (1)

Definition at line 163 of file SinCosModel_decl.hpp.

template<class Scalar >
bool Tempus_Test::SinCosModel< Scalar >::haveIC_
protected

false => no nominal values are provided (default=true)

Definition at line 164 of file SinCosModel_decl.hpp.

template<class Scalar >
bool Tempus_Test::SinCosModel< Scalar >::acceptModelParams_
protected

Changes inArgs to require parameters.

Definition at line 165 of file SinCosModel_decl.hpp.

template<class Scalar >
bool Tempus_Test::SinCosModel< Scalar >::useDfDpAsTangent_
protected

Treat DfDp OutArg as tangent (df/dx*dx/dp+df/dp)

Definition at line 166 of file SinCosModel_decl.hpp.

template<class Scalar >
bool Tempus_Test::SinCosModel< Scalar >::isInitialized_
mutableprotected

Definition at line 167 of file SinCosModel_decl.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs<Scalar> Tempus_Test::SinCosModel< Scalar >::inArgs_
mutableprotected

Definition at line 168 of file SinCosModel_decl.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::OutArgs<Scalar> Tempus_Test::SinCosModel< Scalar >::outArgs_
mutableprotected

Definition at line 169 of file SinCosModel_decl.hpp.

template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs<Scalar> Tempus_Test::SinCosModel< Scalar >::nominalValues_
mutableprotected

Definition at line 170 of file SinCosModel_decl.hpp.

template<class Scalar >
Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > Tempus_Test::SinCosModel< Scalar >::x_space_
protected

Definition at line 171 of file SinCosModel_decl.hpp.

template<class Scalar >
Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > Tempus_Test::SinCosModel< Scalar >::f_space_
protected

Definition at line 172 of file SinCosModel_decl.hpp.

template<class Scalar >
Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > Tempus_Test::SinCosModel< Scalar >::p_space_
protected

Definition at line 173 of file SinCosModel_decl.hpp.

template<class Scalar >
Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > Tempus_Test::SinCosModel< Scalar >::g_space_
protected

Definition at line 174 of file SinCosModel_decl.hpp.

template<class Scalar >
Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > Tempus_Test::SinCosModel< Scalar >::DxDp_space_
protected

Definition at line 175 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::a_
protected

Model parameter.

Definition at line 179 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::f_
protected

Model parameter.

Definition at line 180 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::L_
protected

Model parameter.

Definition at line 181 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::phi_
protected

Parameter determined from the IC.

Definition at line 182 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::b_
protected

Parameter determined from the IC.

Definition at line 183 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::t0_ic_
protected

Time value where the initial condition is specified.

Definition at line 184 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::x0_ic_
protected

Initial condition for x0.

Definition at line 185 of file SinCosModel_decl.hpp.

template<class Scalar >
Scalar Tempus_Test::SinCosModel< Scalar >::x1_ic_
protected

Initial condition for x1.

Definition at line 186 of file SinCosModel_decl.hpp.


The documentation for this class was generated from the following files: