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

Decorator class that wraps any ModelEvaluator object and lumps parameters together using a linear basis matrix. More...

#include <Thyra_DefaultLumpedParameterModelEvaluator.hpp>

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

Related Functions

(Note that these are not member functions.)

template<class Scalar >
RCP
< DefaultLumpedParameterModelEvaluator
< Scalar > > 
defaultLumpedParameterModelEvaluator (const RCP< ModelEvaluator< Scalar > > &thyraModel)
 Non-member constructor. More...
 

Constructors/initializers/accessors/utilities.

 DefaultLumpedParameterModelEvaluator ()
 
void initialize (const RCP< ModelEvaluator< Scalar > > &thyraModel)
 
void uninitialize (RCP< ModelEvaluator< Scalar > > *thyraModel)
 

Public functions overridden from Teuchos::Describable.

std::string description () const
 

Overridden from ParameterListAcceptor

void setParameterList (RCP< Teuchos::ParameterList > const &paramList)
 
RCP< Teuchos::ParameterListgetNonconstParameterList ()
 
RCP< Teuchos::ParameterListunsetParameterList ()
 
RCP< const Teuchos::ParameterListgetParameterList () const
 
RCP< const Teuchos::ParameterListgetValidParameters () const
 

Public functions overridden from ModelEvaulator.

RCP< const VectorSpaceBase
< Scalar > > 
get_p_space (int l) const
 
RCP< const Array< std::string > > get_p_names (int l) const
 
ModelEvaluatorBase::InArgs
< Scalar > 
getNominalValues () const
 
ModelEvaluatorBase::InArgs
< Scalar > 
getLowerBounds () const
 
ModelEvaluatorBase::InArgs
< Scalar > 
getUpperBounds () const
 
void reportFinalPoint (const ModelEvaluatorBase::InArgs< Scalar > &finalPoint, const bool wasSolved)
 

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::ModelEvaluatorDelegatorBase< Scalar >
 ModelEvaluatorDelegatorBase ()
 Constructs to uninitialized. More...
 
 ModelEvaluatorDelegatorBase (const RCP< ModelEvaluator< Scalar > > &model)
 Calls initialize(). More...
 
 ModelEvaluatorDelegatorBase (const RCP< const ModelEvaluator< Scalar > > &model)
 Calls initialize(). More...
 
void initialize (const RCP< ModelEvaluator< Scalar > > &model)
 Initialize given a non-const model evaluator. More...
 
void initialize (const RCP< const ModelEvaluator< Scalar > > &model)
 Initialize given a const model evaluator. More...
 
void uninitialize ()
 Uninitialize. More...
 
virtual bool isUnderlyingModelConst () const
 
virtual RCP< ModelEvaluator
< Scalar > > 
getNonconstUnderlyingModel ()
 
virtual RCP< const
ModelEvaluator< Scalar > > 
getUnderlyingModel () const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_x_space () const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_f_space () const
 
RCP< const VectorSpaceBase
< Scalar > > 
get_f_multiplier_space () const
 
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
 
RCP< const VectorSpaceBase
< Scalar > > 
get_g_multiplier_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< LinearOpWithSolveBase
< Scalar > > 
create_W () 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)
 
- 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
 
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::ModelEvaluatorDelegatorBase< Scalar >
void setLocalVerbosityLevelValidatedParameter (ParameterList *paramList) const
 Set a valid parameter for reading the local verbosity level. More...
 
Teuchos::EVerbosityLevel readLocalVerbosityLevelValidatedParameter (ParameterList &paramList) const
 Read the local verbosity level parameter. More...
 
- 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::DefaultLumpedParameterModelEvaluator< Scalar >

Decorator class that wraps any ModelEvaluator object and lumps parameters together using a linear basis matrix.

Introduction

The purpose of this Decorator class is to provide a linear basis reduction, or "lumping", for a set of parameters. Let p_orig be one of the parameter subvectors of the underlying model *getUnderlyingModel(). This class provides an affine model for the reduced parameters:

  p_orig = B * p + p_orig_base

where B is some MultiVectorBase object with linearly independent columns where B.range()->isCompatible(*p_orig.space())==true and B.domain()->dim() <= B.range()->dim(). The basis matrix B can be set by the client or can be generated automatically in this object for any number of columns. The vector p_orig_base is generally selected to be the nominal values of the original parameters p_orig in which case p == 0 gives the original nominal values and p != 0 gives a perturbation from the nominal values. This is the default option but a different base vector can be set. Note that storing B as a column-wise multi-vector allows the space *p_orig.space() to be any arbitrary vector space, including a parallel distributed vector space in an SPMD program. It is this flexibility that truly makes this decorator class an ANA class.

The reduced parameter subvector p is what is exposed to an ANA client dnd the original parameter subvector p_orig is what is passed to the original underlying model. Note that given p, p_orig is uniquely determined but opposite is not true in general (see Section Mapping between fully and reduced parameters).

Derivatives

The reduced basis B of course affects the definition of all of the derivatives with respect to the given parameter subvector. The relationship between the derivative of any function h with respect to p_orig and p is:

  d(h)/d(p) = d(h)/d(p_orig) * B

When d(h)/d(p) is only needed as a linear operator, both d(h)/d(p_orig) and B would just need to be supplied as general linear operators and then the multiplied linear operator d(h)/d(p_orig)*B could be represented implicitly using a Thyra::DefaultMultipliedLinearOp subclass object.

When d(h)/d(p) is only needed as a column-wise multi-vector, then d(h)/d(p_orig) is only needed as a general linear operator and B must be a column-wise multi-vector.

Lastly, when the row-wise transpose multi-vector form of d(h)/d(p)^T is requested, it can be computed as:

  d(h)/d(p)^T = B^T * d(h)/d(p_orig)^T

which requires that d(h)/d(p_orig)^T be supplied in row-wise transpose multi-vector form but B could actually just be a general linear operator. This would allow for huge spaces for p_orig and p which may be of some use for adjoint-based sensitivity methods.

Since the current implementation in this class requires B be a multi-vector, both forms d(h)/d(p) and d(h)/d(p)^T can be supported. In the future, it would be possible to support the implementation of B as a general linear operator which would mean that very large full and reduced parameter spaces could be supported. However, focusing on small reduced parameter spaces is the very motivation for this decorator subclass and therefore requiring that B be represented as a multi-vector is not a serious limitation.

Mapping between fully and reduced parameters

In cases where p_orig is given and p must be determined, in general there is no solution for p that will satisfy p_orig=B*p+p_orig_base.

To support the (overdetermined) mapping from p_orig to p, we choose p to solve the classic linear least squares problem:

  min   0.5 * (B*p + p_orig_base - p_orig)^T * (B*p + p_orig_base - p_orig)

This well known linear least squares problem has the solution:

  p = inv(B^T * B) * (B^T * (-p_orig_base+p_oirg))

This approach has the unfortunate side effect that we can not completely represent an arbitrary vector p_orig with a reduced p but this is the nature of this approximation for both good and bad.

The one case in which a unique reduced p can be determined is when p_orig was computed from the afffine relationship given p and therefore we expect a zero residual meaning that (p_orig_base-p_oirg) lies in the null-space of B. This is handy when one wants to write out p in the form of p_orig and then reconstruct p again.

Parameter bounds

As mentioned above, one can simply select p_orig_base to be the nominal values of the original full parameters p_orig and therefore allow p==0 to give the exact nominal parameter values which is critical for the initial guess for many ANAs. The problem of handling the bounds on the parameters is more difficult. This approximation really requires that the simple bounds on p_orig be replaced with the general linear inequality constraints:

  p_orig_l - p_orig_base <= B*p <= p_orig_u - p_orig_base

Modeling and enforcing these linear inequality constraints correctly would require that B and p_orig_base be exposed to the client so that the client can build these extra linear inequality constraints into the ANA algorithm. Certain optimization algorithms could then enforce these general inequalities and therefore guarantee that the original parameter bounds are never violated. This is easy to do in both active-set and interior-point methods when the size of the space p_orig.space()->dim() is not too large. When p_orig.space()->dim() is large, handling these inequality constraints can be very difficult to deal with.

An alternative simpler approach would be to try to find lower and upper bounds p_l and p_u that tried to match the original bounds as close as possible while not restricting the feasible region p_l <= p <= p_u too much. For example, one could solve the inequality-constrained least squares problems:

  min     0.5 * (B*p+p_orig_base-p_orig_l)^T * (B*p+p_orig_base-p_orig_l)
  s.t.    B*p+p_orig_base >= p_orig_l

and

  min     0.5 * (B*p+p_orig_base-p_orig_u)^T * (B*p+p_orig_base-p_orig_u)
  s.t.    B*p+p_orig_base <= p_orig_u

but in general it would not be possible to guarantee that the solution p_l and p_u satisfies p_l <= p <= p_u.

Because of the challenges of dealing with bounds on the parameters, this subclass currently throws an exception if it is given an underlying model with finite bounds on the parameters. However, a parameter-list option can be set that will cause the bounds to be ignored and it would be the client's responsibility to deal with the implications of this choice.

Definition at line 240 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.

Constructor & Destructor Documentation

Member Function Documentation

template<class Scalar >
void Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::initialize ( const RCP< ModelEvaluator< Scalar > > &  thyraModel)
template<class Scalar >
void Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::uninitialize ( RCP< ModelEvaluator< Scalar > > *  thyraModel)
template<class Scalar >
std::string Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::description ( ) const
virtual

Reimplemented from Teuchos::Describable.

Definition at line 567 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.

template<class Scalar >
void Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::setParameterList ( RCP< Teuchos::ParameterList > const &  paramList)
virtual
template<class Scalar >
RCP< Teuchos::ParameterList > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::getNonconstParameterList ( )
virtual
template<class Scalar >
RCP< Teuchos::ParameterList > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::unsetParameterList ( )
virtual
template<class Scalar >
RCP< const Teuchos::ParameterList > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::getParameterList ( ) const
virtual
template<class Scalar >
RCP< const Teuchos::ParameterList > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::getValidParameters ( ) const
virtual
template<class Scalar >
RCP< const VectorSpaceBase< Scalar > > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::get_p_space ( int  l) const
virtual
template<class Scalar >
RCP< const Array< std::string > > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::get_p_names ( int  l) const
virtual
template<class Scalar >
ModelEvaluatorBase::InArgs< Scalar > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::getNominalValues ( ) const
virtual
template<class Scalar >
ModelEvaluatorBase::InArgs< Scalar > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::getLowerBounds ( ) const
virtual
template<class Scalar >
ModelEvaluatorBase::InArgs< Scalar > Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::getUpperBounds ( ) const
virtual
template<class Scalar >
void Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::reportFinalPoint ( const ModelEvaluatorBase::InArgs< Scalar > &  finalPoint,
const bool  wasSolved 
)
virtual

Friends And Related Function Documentation

template<class Scalar >
RCP< DefaultLumpedParameterModelEvaluator< Scalar > > defaultLumpedParameterModelEvaluator ( const RCP< ModelEvaluator< Scalar > > &  thyraModel)
related

Non-member constructor.

Definition at line 433 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.


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