Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Thyra::DefaultMultiPeriodModelEvaluator< Scalar > Class Template Reference

Composite subclass that takes a single ModelEvaluator object and represents it as a single aggregate multi-period ModelEvalator object. More...

#include <Thyra_DefaultMultiPeriodModelEvaluator.hpp>

Inheritance diagram for Thyra::DefaultMultiPeriodModelEvaluator< Scalar >:
Inheritance graph
[legend]

Constructors/Intializers/Accessors

 DefaultMultiPeriodModelEvaluator ()
 
 DefaultMultiPeriodModelEvaluator (const int N, const Array< RCP< ModelEvaluator< Scalar > > > &periodModels, const Array< int > &z_indexes, const Array< Array< RCP< const VectorBase< Scalar > > > > &z, const int g_index, const Array< Scalar > g_weights, const RCP< const ProductVectorSpaceBase< Scalar > > &x_bar_space=Teuchos::null, const RCP< const ProductVectorSpaceBase< Scalar > > &f_bar_space=Teuchos::null, const RCP< LinearOpWithSolveFactoryBase< Scalar > > &W_bar_factory=Teuchos::null)
 Calls intialize(...). More...
 
void initialize (const int N, const Array< RCP< ModelEvaluator< Scalar > > > &periodModels, const Array< int > &z_indexes, const Array< Array< RCP< const VectorBase< Scalar > > > > &z, const int g_index, const Array< Scalar > g_weights, const RCP< const ProductVectorSpaceBase< Scalar > > &x_bar_space=Teuchos::null, const RCP< const ProductVectorSpaceBase< Scalar > > &f_bar_space=Teuchos::null, const RCP< LinearOpWithSolveFactoryBase< Scalar > > &W_bar_factory=Teuchos::null)
 Initialize. More...
 
void reset_z (const Array< Array< RCP< const VectorBase< Scalar > > > > &z)
 Reset z. More...
 

Public functions overridden from ModelEvaulator.

int Np () const
 
int Ng () const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_x_space () const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_f_space () const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_p_space (int l) const
 
RCP< const Array< std::string > > get_p_names (int l) const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_g_space (int j) const
 
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
 
ModelEvaluatorBase::InArgs
< Scalar > 
createInArgs () const
 
void reportFinalPoint (const ModelEvaluatorBase::InArgs< Scalar > &finalPoint, const bool wasSolved)
 Ignores the final point. More...
 

Additional Inherited Members

- Public Types inherited from Thyra::ModelEvaluator< Scalar >
typedef Teuchos::ScalarTraits
< Scalar >::magnitudeType 
ScalarMag
 
- Public Types inherited from Thyra::ModelEvaluatorBase
enum  EInArgsMembers {
  IN_ARG_x_dot_dot, IN_ARG_x_dot, IN_ARG_x, IN_ARG_x_dot_poly,
  IN_ARG_x_poly, IN_ARG_x_dot_mp, IN_ARG_x_mp, IN_ARG_t,
  IN_ARG_alpha, IN_ARG_beta, IN_ARG_W_x_dot_dot_coeff, IN_ARG_step_size,
  IN_ARG_stage_number
}
 
enum  EInArgs_p_mp { IN_ARG_p_mp }
 
enum  EEvalType { , EVAL_TYPE_APPROX_DERIV, EVAL_TYPE_VERY_APPROX_DERIV }
 The type of an evaluation. More...
 
enum  EDerivativeMultiVectorOrientation { DERIV_MV_JACOBIAN_FORM, DERIV_MV_GRADIENT_FORM, DERIV_MV_BY_COL = DERIV_MV_JACOBIAN_FORM, DERIV_TRANS_MV_BY_ROW = DERIV_MV_GRADIENT_FORM }
 
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_mp, OUT_ARG_W_mp,
  OUT_ARG_W_op, OUT_ARG_W_prec, OUT_ARG_f_poly
}
 
enum  EOutArgsDfDp { OUT_ARG_DfDp }
 
enum  EOutArgsDgDx_dot { OUT_ARG_DgDx_dot }
 
enum  EOutArgsDgDx { OUT_ARG_DgDx }
 
enum  EOutArgsDgDp { OUT_ARG_DgDp }
 
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_mp { OUT_ARG_DgDx_mp }
 
enum  EOutArgsDgDp_mp { OUT_ARG_DgDp_mp }
 
- Public Member Functions inherited from Thyra::ModelEvaluatorDefaultBase< Scalar >
int Np () const
 
int Ng () const
 
RCP< LinearOpBase< Scalar > > create_DfDp_op (int l) const
 
RCP< LinearOpBase< Scalar > > create_DgDx_dot_op (int j) const
 
RCP< LinearOpBase< Scalar > > create_DgDx_op (int j) const
 
RCP< LinearOpBase< Scalar > > create_DgDp_op (int j, int l) const
 
RCP< LinearOpWithSolveBase
< Scalar > > 
create_W () const
 
ModelEvaluatorBase::OutArgs
< Scalar > 
createOutArgs () const
 
void evalModel (const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
get_f_multiplier_space () const
 
virtual RCP< const
VectorSpaceBase< Scalar > > 
get_g_multiplier_space (int j) const
 
- Public Member Functions inherited from Thyra::ModelEvaluator< Scalar >
- Public Member Functions inherited from Thyra::ModelEvaluatorBase
 ModelEvaluatorBase ()
 constructor More...
 
- Static Public Attributes inherited from Thyra::ModelEvaluatorBase
static const int NUM_E_IN_ARGS_MEMBERS =13
 
static const int NUM_E_OUT_ARGS_MEMBERS =7
 
- Protected Member Functions inherited from Thyra::ModelEvaluatorDefaultBase< Scalar >
 ModelEvaluatorDefaultBase ()
 
void initializeDefaultBase ()
 Function called by subclasses to fully initialize this object on any important change. More...
 
void resetDefaultBase ()
 Sets the the DefaultBase to an uninitialized state, forcing lazy initialization when needed. More...
 

Detailed Description

template<class Scalar>
class Thyra::DefaultMultiPeriodModelEvaluator< Scalar >

Composite subclass that takes a single ModelEvaluator object and represents it as a single aggregate multi-period ModelEvalator object.

Mathematically, this subclass is used to form a multi-period (or multi-point) problem of the form:

  (x_bar,{p_l}\p_period) -> f_bar

  (x_bar,{p_l}\p_period) -> g_bar

where

  x_bar = [ x[0]; x[1]; ...; x[N-1]]

  {p_l}\p_period = { p_(0), p_(1), ..., p_(z_index-1), p_(z_index+1), ..., p_(Np-1) }

  f_bar(...) = [ f(x[0],{p_l},z[0]); f(x[1],{p_l},z[1]); ..., f(x[N-1],{p_l},z[N-1]) ]

  g_bar(...) = sum( g_weights[i]*g[g_index](x[0],{p_l},z[i]), i = 0,...,N-1 )

Above, the notation {p_l} is meant to represent the set of all parameter subvectors in each of the constituent models minus the parameter subvector used to define period data.

This gives the first derivative objects:

                     [ W[0], 0,     ...  0      ]
                     [ 0,    W[i],  ...  0      ]
  W_bar = DfDx_bar = [ .,    .,     ...  .      ]
                     [ 0,    0,     ...  W[N-1] ]

 
                [ DfDp[0][l]   ]
                [ DfDp[1][l]   ]
  DfDp_bar[l] = [ ...          ]
                [ DfDp[N=1][l] ]


                 [ g_wieght[0] * DgDx_dot[0]     ]
                 [ g_wieght[1] * DgDx_dot[1]     ]
  DgDx_dot_bar = [ ...                           ]
                 [ g_wieght[N-1] * DgDx_dot[N-1] ]


             [ g_wieght[0] * DgDx[0]     ]
             [ g_wieght[1] * DgDx[1]     ]
  DgDx_bar = [ ...                       ]
             [ g_wieght[N-1] * DgDx[N-1] ]


  DgDp_bar = sum( g_weights[i]*DgDp[i][g_index], i = 0,...,N-1 )

This class could be made much more general but for now ???

Definition at line 131 of file Thyra_DefaultMultiPeriodModelEvaluator.hpp.

Constructor & Destructor Documentation

Definition at line 368 of file Thyra_DefaultMultiPeriodModelEvaluator.hpp.

template<class Scalar >
Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::DefaultMultiPeriodModelEvaluator ( const int  N,
const Array< RCP< ModelEvaluator< Scalar > > > &  periodModels,
const Array< int > &  z_indexes,
const Array< Array< RCP< const VectorBase< Scalar > > > > &  z,
const int  g_index,
const Array< Scalar >  g_weights,
const RCP< const ProductVectorSpaceBase< Scalar > > &  x_bar_space = Teuchos::null,
const RCP< const ProductVectorSpaceBase< Scalar > > &  f_bar_space = Teuchos::null,
const RCP< LinearOpWithSolveFactoryBase< Scalar > > &  W_bar_factory = Teuchos::null 
)

Calls intialize(...).

Definition at line 374 of file Thyra_DefaultMultiPeriodModelEvaluator.hpp.

Member Function Documentation

template<class Scalar >
void Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::initialize ( const int  N,
const Array< RCP< ModelEvaluator< Scalar > > > &  periodModels,
const Array< int > &  z_indexes,
const Array< Array< RCP< const VectorBase< Scalar > > > > &  z,
const int  g_index,
const Array< Scalar >  g_weights,
const RCP< const ProductVectorSpaceBase< Scalar > > &  x_bar_space = Teuchos::null,
const RCP< const ProductVectorSpaceBase< Scalar > > &  f_bar_space = Teuchos::null,
const RCP< LinearOpWithSolveFactoryBase< Scalar > > &  W_bar_factory = Teuchos::null 
)

Initialize.

Parameters
N[in] The number of periods.
periodModels[in] Array (length N) of the period models. For now, each of the period models at periodModels[i] must have an identical structure for every input and output. The reason that different models are passed in is so that different behaviors of the output functions and different nominal values for each period's x and bounds can be specified. The nominal values from periodModels[0] will be used to set the nominal values for the non-peroid parameters.
z_indexes[in] Array of sorted zero-based indexes of the model's parameter subvector that represents z[i]. Note, this array must be sorted from smallest to largets and can not have any duplicate entires!
z[in] Array (length N) of the period defining auxiliary paramter subvectors. Each entry z[i] is an array of subvectors where z[i][k] is the subvector p(z_indexes[k]) in the underlying perild model. Therefore, the array z[i] must be ordered according to z_indexes. Note that z[i][k] is allowed to be null in which case the underlying default value for this parameter will be used.
g_index[in] The index of the response function that will be used for the period objective function.
g_weights[in] Array (length N) of the g_weights for the auxiliary response functions in the summation for g_bar(...) shown above.
x_bar_space[in] The product vector space that represents the space for x_bar as defined above. If x_bar_space.get()==NULL then a default version of this product space will be created internally.
f_bar_space[in] The product vector space that represents the space for f_bar as defined above. If f_bar_space.get()==NULL then a default version of this product space will be created internally.
W_bar_factory[in] Factory object that is used to create the block LOWS object that will be used to represent the block diagonal object W_bar. If is_null(W_bar_factory)==true on input then a DefaultBlockedTriangularLinearOpWithSolveFactory object will be created and used internally.

Preconditions:

  • N > 0
  • periodModels.size()==N.
  • z_indexes.size() > 0 and z_indexes[k] >= 0 and z_indexes[k] are sorted low to high and are unique.
  • z.size()==N and z[i].size()==z_indexes.size(), for i=0...N-1

Definition at line 395 of file Thyra_DefaultMultiPeriodModelEvaluator.hpp.

template<class Scalar >
void Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::reset_z ( const Array< Array< RCP< const VectorBase< Scalar > > > > &  z)

Reset z.

Parameters
z[in] See initialize().

Preconditions:

Definition at line 482 of file Thyra_DefaultMultiPeriodModelEvaluator.hpp.

template<class Scalar >
int Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::Np ( ) const
virtual
template<class Scalar >
int Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::Ng ( ) const
virtual
template<class Scalar >
RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::get_x_space ( ) const
virtual
template<class Scalar >
RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::get_f_space ( ) const
virtual
template<class Scalar >
RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::get_p_space ( int  l) const
virtual
template<class Scalar >
RCP< const Array< std::string > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::get_p_names ( int  l) const
virtual
template<class Scalar >
RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::get_g_space ( int  j) const
virtual
template<class Scalar >
ArrayView< const std::string > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::get_g_names ( int  j) const
virtual
template<class Scalar >
ModelEvaluatorBase::InArgs< Scalar > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::getNominalValues ( ) const
virtual
template<class Scalar >
ModelEvaluatorBase::InArgs< Scalar > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::getLowerBounds ( ) const
virtual
template<class Scalar >
ModelEvaluatorBase::InArgs< Scalar > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::getUpperBounds ( ) const
virtual
template<class Scalar >
RCP< LinearOpBase< Scalar > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::create_W_op ( ) const
virtual
template<class Scalar >
RCP< PreconditionerBase< Scalar > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::create_W_prec ( ) const
virtual
template<class Scalar >
RCP< const LinearOpWithSolveFactoryBase< Scalar > > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::get_W_factory ( ) const
virtual
template<class Scalar >
ModelEvaluatorBase::InArgs< Scalar > Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::createInArgs ( ) const
virtual
template<class Scalar >
void Thyra::DefaultMultiPeriodModelEvaluator< Scalar >::reportFinalPoint ( const ModelEvaluatorBase::InArgs< Scalar > &  finalPoint,
const bool  wasSolved 
)
virtual

Ignores the final point.

Implements Thyra::ModelEvaluator< Scalar >.

Definition at line 645 of file Thyra_DefaultMultiPeriodModelEvaluator.hpp.


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