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

Simple parallel response-only ModelEvaluator. More...

#include <Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_decl.hpp>

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

Constructors/Initializers/Accessors.

 DiagonalQuadraticResponseOnlyModelEvaluator (const int localDim, const RCP< const Teuchos::Comm< Ordinal > > &comm=Teuchos::null)
 
void setSolutionVector (const RCP< const VectorBase< Scalar > > &ps)
 Set the solution vector ps . More...
 
const RCP< const VectorBase
< Scalar > > 
getSolutionVector () const
 Get the solution vector ps . More...
 
void setDiagonalVector (const RCP< const VectorBase< Scalar > > &diag)
 Set the diagonal vector diag. More...
 
void setDiagonalBarVector (const RCP< const VectorBase< Scalar > > &diag_bar)
 Set the diagonal vector diag_bar. More...
 
void setNonlinearTermFactor (const Scalar &nonlinearTermFactor)
 Set nonlinear term factory. More...
 
void setScalarOffset (const Scalar &g_offset)
 Set offset scalar g_offset . More...
 

Public functions overridden from ModelEvaulator.

int Np () const
 
int Ng () const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_p_space (int l) const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_g_space (int j) const
 
ModelEvaluatorBase::InArgs
< Scalar > 
createInArgs () const
 

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::ResponseOnlyModelEvaluatorBase< Scalar >
RCP< const VectorSpaceBase
< Scalar > > 
get_x_space () const
 Throws exception. More...
 
RCP< const Teuchos::Array
< std::string > > 
get_p_names (int l) const
 Returns null. More...
 
Teuchos::ArrayView< const
std::string > 
get_g_names (int j) const
 Returns null. More...
 
RCP< const VectorSpaceBase
< Scalar > > 
get_f_space () const
 Throws exception. More...
 
ModelEvaluatorBase::InArgs
< Scalar > 
getNominalValues () const
 Returns this->createInArgs(). More...
 
ModelEvaluatorBase::InArgs
< Scalar > 
getLowerBounds () const
 Returns this->createInArgs(). More...
 
ModelEvaluatorBase::InArgs
< Scalar > 
getUpperBounds () const
 Returns this->createInArgs(). More...
 
RCP< LinearOpWithSolveBase
< Scalar > > 
create_W () const
 Thorws exception. More...
 
RCP< LinearOpBase< Scalar > > create_W_op () const
 Thorws exception. More...
 
RCP< PreconditionerBase< Scalar > > create_W_prec () const
 Thorws exception. More...
 
RCP< const
LinearOpWithSolveFactoryBase
< Scalar > > 
get_W_factory () const
 Thorws exception. More...
 
void reportFinalPoint (const ModelEvaluatorBase::InArgs< Scalar > &finalPoint, const bool wasSolved)
 Does nothing and ignores input. More...
 
- Public Member Functions inherited from Thyra::ModelEvaluatorDefaultBase< Scalar >
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
 
ModelEvaluatorBase::OutArgs
< Scalar > 
createOutArgs () const
 
void evalModel (const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs) 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::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >

Simple parallel response-only ModelEvaluator.

The representation of the model in coefficient form is:

  g(p) = 0.5 * sum( diag[i] * (p[i] - ps[i])^2, i=0...n-1 )
         + 0.5 * nonlinearTermFactor * sum( (p[i] - ps[i])^3, i=0...n-1 )
         + g_offset

In vector coefficient form it becomes:

  g(p) = 0.5 * (p-ps)^T * D * (p-ps) + cubicTerm(p) + g_offset

where D = diag(diag) and cubitTerm(p) is given above.

This test model implementation also supports the definition of diagonal scalar product implementation. The impact of this is the definition of new abstract Thyra Euclidean vectors of the form:

  p_bar = E_p * p

  g_grad_bar = E_p * g_grad

In this notation, we say that the Thyra vectors themselves p_bar are vectors in a Euclidean space (from the standpoint of the ANA clients) and that p are the coefficient vectors in the new non-Euclidean space S_p with scalar product:

  y_bar^T * x_bar = y^T * E_p^T * E_p * x
                  = y^T * S_p * x
                  = <y, x>_p

The ideas is that by providing this new vector space basis E_p and inner product <.,.>_p, we change the function that the ANA sees.

This testing class allows a client to specify the desirable diagonal D_bar matrix by passing in the vector diag_bar. From diag_bar and diag, the diagonal s_diag for S_p = diag(s_diag) as:

  s_diag[i] = diag_bar[i] / diag[i]

The inner product is therefore:

  scalarProd(y_bar, x_bar) = sum( s_diag[i] * y[i] * x[i], i=0...n-1 )

Definition at line 126 of file Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::DiagonalQuadraticResponseOnlyModelEvaluator ( const int  localDim,
const RCP< const Teuchos::Comm< Ordinal > > &  comm = Teuchos::null 
)

Member Function Documentation

template<class Scalar >
void Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::setSolutionVector ( const RCP< const VectorBase< Scalar > > &  ps)

Set the solution vector ps .

Definition at line 117 of file Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_def.hpp.

template<class Scalar >
const RCP< const Thyra::VectorBase< Scalar > > Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::getSolutionVector ( ) const

Get the solution vector ps .

Definition at line 126 of file Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_def.hpp.

template<class Scalar >
void Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::setDiagonalVector ( const RCP< const VectorBase< Scalar > > &  diag)

Set the diagonal vector diag.

Definition at line 133 of file Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_def.hpp.

template<class Scalar >
void Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::setDiagonalBarVector ( const RCP< const VectorBase< Scalar > > &  diag_bar)

Set the diagonal vector diag_bar.

NOTE: You must call setDiagonalVector(diag) first in order to set the objective diagonal.

Definition at line 142 of file Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_def.hpp.

template<class Scalar >
void Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::setNonlinearTermFactor ( const Scalar &  nonlinearTermFactor)

Set nonlinear term factory.

Definition at line 174 of file Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_def.hpp.

template<class Scalar >
void Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::setScalarOffset ( const Scalar &  g_offset)

Set offset scalar g_offset .

Definition at line 182 of file Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_def.hpp.

template<class Scalar >
int Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::Np ( ) const
virtual
template<class Scalar >
int Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::Ng ( ) const
virtual
template<class Scalar >
RCP< const Thyra::VectorSpaceBase< Scalar > > Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::get_p_space ( int  l) const
virtual
template<class Scalar >
RCP< const Thyra::VectorSpaceBase< Scalar > > Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::get_g_space ( int  j) const
virtual
template<class Scalar >
Thyra::ModelEvaluatorBase::InArgs< Scalar > Thyra::DiagonalQuadraticResponseOnlyModelEvaluator< Scalar >::createInArgs ( ) const
virtual

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