EpetraExt  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Related Functions | List of all members
EpetraExt::DiagonalTransientModel Class Reference

Simple transient diagonal model for an implicit or explicit ODE. More...

#include <EpetraExt_DiagonalTransientModel.hpp>

Inheritance diagram for EpetraExt::DiagonalTransientModel:
Inheritance graph
[legend]

Public Types

enum  EGammaFit { GAMMA_FIT_LINEAR, GAMMA_FIT_RANDOM }
 
- Public Types inherited from EpetraExt::ModelEvaluator
enum  EInArgsMembers {
  IN_ARG_x_dot, IN_ARG_x, IN_ARG_x_dot_poly, IN_ARG_x_poly,
  IN_ARG_x_dot_sg, IN_ARG_x_sg, IN_ARG_x_dot_mp, IN_ARG_x_mp,
  IN_ARG_t, IN_ARG_alpha, IN_ARG_beta, IN_ARG_step_size,
  IN_ARG_stage_number, IN_ARG_x_dotdot, IN_ARG_x_dotdot_poly, IN_ARG_x_dotdot_sg,
  IN_ARG_x_dotdot_mp, IN_ARG_omega, IN_ARG_sg_basis, IN_ARG_sg_quadrature,
  IN_ARG_sg_expansion
}
 
enum  EInArgs_p_sg { IN_ARG_p_sg }
 
enum  EInArgs_p_mp { IN_ARG_p_mp }
 
enum  EEvalType { EVAL_TYPE_EXACT, EVAL_TYPE_APPROX_DERIV, EVAL_TYPE_VERY_APPROX_DERIV }
 
enum  EDerivativeMultiVectorOrientation { DERIV_MV_BY_COL, DERIV_TRANS_MV_BY_ROW }
 
enum  EDerivativeLinearOp { DERIV_LINEAR_OP }
 
enum  EDerivativeLinearity { DERIV_LINEARITY_UNKNOWN, DERIV_LINEARITY_CONST, DERIV_LINEARITY_NONCONST }
 
enum  ERankStatus { DERIV_RANK_UNKNOWN, DERIV_RANK_FULL, DERIV_RANK_DEFICIENT }
 
enum  EOutArgsMembers {
  OUT_ARG_f, OUT_ARG_W, OUT_ARG_f_poly, OUT_ARG_f_sg,
  OUT_ARG_W_sg, OUT_ARG_f_mp, OUT_ARG_W_mp, OUT_ARG_WPrec
}
 
enum  EOutArgsDfDp { OUT_ARG_DfDp }
 
enum  EOutArgsDgDx_dot { OUT_ARG_DgDx_dot }
 
enum  EOutArgsDgDx_dotdot { OUT_ARG_DgDx_dotdot }
 
enum  EOutArgsDgDx { OUT_ARG_DgDx }
 
enum  EOutArgsDgDp { OUT_ARG_DgDp }
 
enum  EOutArgsDfDp_sg { OUT_ARG_DfDp_sg }
 
enum  EOutArgs_g_sg { OUT_ARG_g_sg }
 
enum  EOutArgsDgDx_dot_sg { OUT_ARG_DgDx_dot_sg }
 
enum  EOutArgsDgDx_dotdot_sg { OUT_ARG_DgDx_dotdot_sg }
 
enum  EOutArgsDgDx_sg { OUT_ARG_DgDx_sg }
 
enum  EOutArgsDgDp_sg { OUT_ARG_DgDp_sg }
 
enum  EOutArgsDfDp_mp { OUT_ARG_DfDp_mp }
 
enum  EOutArgs_g_mp { OUT_ARG_g_mp }
 
enum  EOutArgsDgDx_dot_mp { OUT_ARG_DgDx_dot_mp }
 
enum  EOutArgsDgDx_dotdot_mp { OUT_ARG_DgDx_dotdot_mp }
 
enum  EOutArgsDgDx_mp { OUT_ARG_DgDx_mp }
 
enum  EOutArgsDgDp_mp { OUT_ARG_DgDp_mp }
 
typedef Teuchos::RCP< const
Stokhos::ProductEpetraVector > 
mp_const_vector_t
 
typedef Teuchos::RCP< const
Stokhos::ProductEpetraMultiVector > 
mp_const_multivector_t
 
typedef Teuchos::RCP< const
Stokhos::ProductEpetraOperator > 
mp_const_operator_t
 
typedef Teuchos::RCP
< Stokhos::ProductEpetraVector > 
mp_vector_t
 
typedef Teuchos::RCP
< Stokhos::ProductEpetraMultiVector > 
mp_multivector_t
 
typedef Teuchos::RCP
< Stokhos::ProductEpetraOperator > 
mp_operator_t
 

Related Functions

(Note that these are not member functions.)

Teuchos::RCP
< DiagonalTransientModel
diagonalTransientModel (Teuchos::RCP< Epetra_Comm > const &epetra_comm, Teuchos::RCP< Teuchos::ParameterList > const &paramList=Teuchos::null)
 Nonmember constructor. More...
 

Constructors, Initializers, Misc.

 DiagonalTransientModel (Teuchos::RCP< Epetra_Comm > const &epetra_comm)
 
Teuchos::RCP< const Epetra_Vectorget_gamma () const
 Return the model vector gamma,. More...
 
Teuchos::RCP< const Epetra_VectorgetExactSolution (const double t, const Epetra_Vector *coeff_s_p=0) const
 Return the exact solution as a function of time. More...
 
Teuchos::RCP< const
Epetra_MultiVector
getExactSensSolution (const double t, const Epetra_Vector *coeff_s_p=0) const
 Return the exact sensitivity of x as a function of time. More...
 

Overridden from ParameterListAcceptor

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

Overridden from EpetraExt::ModelEvaluator .

Teuchos::RCP< const Epetra_Mapget_x_map () const
 
Teuchos::RCP< const Epetra_Mapget_f_map () const
 
Teuchos::RCP< const Epetra_Mapget_p_map (int l) const
 . More...
 
Teuchos::RCP< const
Teuchos::Array< std::string > > 
get_p_names (int l) const
 . More...
 
Teuchos::RCP< const Epetra_Mapget_g_map (int j) const
 . More...
 
Teuchos::RCP< const Epetra_Vectorget_x_init () const
 
Teuchos::RCP< const Epetra_Vectorget_x_dot_init () const
 
Teuchos::RCP< const Epetra_Vectorget_p_init (int l) const
 
Teuchos::RCP< Epetra_Operatorcreate_W () const
 
InArgs createInArgs () const
 
OutArgs createOutArgs () const
 
void evalModel (const InArgs &inArgs, const OutArgs &outArgs) const
 

Detailed Description

Simple transient diagonal model for an implicit or explicit ODE.

The explicit ODE form of the model is:

 x_dot(i) = f_hat(x(i), gamma(i), s(i), t), for i = 0...n-1, on t in [0,t_f]

where:

 f_hat(x(i), gamma(i), s(i), t) = gama(i)*x(i) + exp(gamma(i)*t)*sin(s(i),t)

The implicit ODE form of the model i:

 f(i)(x_dot(i), x(i), t) = x_dot(i) - f_hat(x(i), gamma(i), s(i), t),
 
   for i = 0...n-1, on t in [0,t_f]

This is a diagonal problem so it does not make the greatest test problem but it does make it easy to derive tests for as a starter.

The coefficients s can be exposed as model parameters and are called coeff_s_p in the code. The selection of the coefficients is handled through the

ToDo: Finish Documentation!

Definition at line 99 of file EpetraExt_DiagonalTransientModel.hpp.

Member Enumeration Documentation

Enumerator
GAMMA_FIT_LINEAR 
GAMMA_FIT_RANDOM 

Definition at line 179 of file EpetraExt_DiagonalTransientModel.hpp.

Constructor & Destructor Documentation

EpetraExt::DiagonalTransientModel::DiagonalTransientModel ( Teuchos::RCP< Epetra_Comm > const &  epetra_comm)

Definition at line 166 of file EpetraExt_DiagonalTransientModel.cpp.

Member Function Documentation

Teuchos::RCP< const Epetra_Vector > EpetraExt::DiagonalTransientModel::get_gamma ( ) const

Return the model vector gamma,.

Definition at line 185 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Vector > EpetraExt::DiagonalTransientModel::getExactSolution ( const double  t,
const Epetra_Vector coeff_s_p = 0 
) const

Return the exact solution as a function of time.

Definition at line 192 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_MultiVector > EpetraExt::DiagonalTransientModel::getExactSensSolution ( const double  t,
const Epetra_Vector coeff_s_p = 0 
) const

Return the exact sensitivity of x as a function of time.

Definition at line 211 of file EpetraExt_DiagonalTransientModel.cpp.

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

Definition at line 237 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< Teuchos::ParameterList > EpetraExt::DiagonalTransientModel::getNonconstParameterList ( )

Definition at line 260 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< Teuchos::ParameterList > EpetraExt::DiagonalTransientModel::unsetParameterList ( )

Definition at line 267 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Teuchos::ParameterList > EpetraExt::DiagonalTransientModel::getParameterList ( ) const

Definition at line 276 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Teuchos::ParameterList > EpetraExt::DiagonalTransientModel::getValidParameters ( ) const

Definition at line 283 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Map > EpetraExt::DiagonalTransientModel::get_x_map ( ) const
virtual

Implements EpetraExt::ModelEvaluator.

Definition at line 327 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Map > EpetraExt::DiagonalTransientModel::get_f_map ( ) const
virtual

Implements EpetraExt::ModelEvaluator.

Definition at line 334 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Map > EpetraExt::DiagonalTransientModel::get_p_map ( int  l) const
virtual

.

Reimplemented from EpetraExt::ModelEvaluator.

Definition at line 341 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Teuchos::Array< std::string > > EpetraExt::DiagonalTransientModel::get_p_names ( int  l) const
virtual

.

Reimplemented from EpetraExt::ModelEvaluator.

Definition at line 351 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Map > EpetraExt::DiagonalTransientModel::get_g_map ( int  j) const
virtual

.

Reimplemented from EpetraExt::ModelEvaluator.

Definition at line 361 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Vector > EpetraExt::DiagonalTransientModel::get_x_init ( ) const
virtual

Reimplemented from EpetraExt::ModelEvaluator.

Definition at line 371 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Vector > EpetraExt::DiagonalTransientModel::get_x_dot_init ( ) const
virtual

Reimplemented from EpetraExt::ModelEvaluator.

Definition at line 378 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< const Epetra_Vector > EpetraExt::DiagonalTransientModel::get_p_init ( int  l) const
virtual

Reimplemented from EpetraExt::ModelEvaluator.

Definition at line 385 of file EpetraExt_DiagonalTransientModel.cpp.

Teuchos::RCP< Epetra_Operator > EpetraExt::DiagonalTransientModel::create_W ( ) const
virtual

Reimplemented from EpetraExt::ModelEvaluator.

Definition at line 395 of file EpetraExt_DiagonalTransientModel.cpp.

EpetraExt::ModelEvaluator::InArgs EpetraExt::DiagonalTransientModel::createInArgs ( ) const
virtual

Implements EpetraExt::ModelEvaluator.

Definition at line 404 of file EpetraExt_DiagonalTransientModel.cpp.

EpetraExt::ModelEvaluator::OutArgs EpetraExt::DiagonalTransientModel::createOutArgs ( ) const
virtual

Implements EpetraExt::ModelEvaluator.

Definition at line 420 of file EpetraExt_DiagonalTransientModel.cpp.

void EpetraExt::DiagonalTransientModel::evalModel ( const InArgs inArgs,
const OutArgs outArgs 
) const
virtual

Implements EpetraExt::ModelEvaluator.

Definition at line 457 of file EpetraExt_DiagonalTransientModel.cpp.

Friends And Related Function Documentation

Teuchos::RCP< DiagonalTransientModel > diagonalTransientModel ( Teuchos::RCP< Epetra_Comm > const &  epetra_comm,
Teuchos::RCP< Teuchos::ParameterList > const &  paramList = Teuchos::null 
)
related

Nonmember constructor.


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