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

This class wraps any ModelEvaluator object and adds a simple, but fairly general, inverse response function. More...

#include <Thyra_DefaultInverseModelEvaluator.hpp>

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

Public Member Functions

 STANDARD_CONST_COMPOSITION_MEMBERS (VectorBase< Scalar >, observationTarget)
 Observation target vector ot. More...
 
 STANDARD_CONST_COMPOSITION_MEMBERS (VectorBase< Scalar >, parameterBase)
 Parameter base vector pt. More...
 
 STANDARD_CONST_COMPOSITION_MEMBERS (LinearOpBase< Scalar >, observationMatchWeightingOp)
 Observation match weighting operator Q_o. More...
 
 STANDARD_CONST_COMPOSITION_MEMBERS (LinearOpBase< Scalar >, parameterRegularizationWeightingOp)
 Parameter regulization weighting operator Q_p. More...
 
 STANDARD_NONCONST_COMPOSITION_MEMBERS (MultiVectorFileIOBase< Scalar >, observationTargetIO)
 MultiVectorFileIOBase object used to read the observation target vector ot as directed by the parameter list. More...
 
 STANDARD_NONCONST_COMPOSITION_MEMBERS (MultiVectorFileIOBase< Scalar >, parameterBaseIO)
 MultiVectorFileIOBase object used to read the parameter base vector pt as directed by the parameter list. More...
 
- 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...
 

Related Functions

(Note that these are not member functions.)

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

Constructors/initializers/accessors/utilities.

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

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 Teuchos::Describable.

std::string description () const
 

Public functions overridden from ModelEvaulator.

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 }
 
- 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::DefaultInverseModelEvaluator< Scalar >

This class wraps any ModelEvaluator object and adds a simple, but fairly general, inverse response function.

The following response function is added to the end of the supported response functions:

 g_inv(x,p)
   = g_(getUnderlyingModel()->Ng())(x,p,...)
   = observationMultiplier * observationMatch(x,p)
   + parameterMultiplier * parameterRegularization(p)

where observationMatch(x,p) is some scalar-valued function that gives the match of some state observation, parameterRegularization(p) is some scaled valued function that regularizes the parameters, and observationMultiplier and parameterMultiplier are scalar constant multipliers for the state observation and the parameter regularization respectively.

The state observation matching function and the parameter regularization function can be defined in one of two ways.

If a symmetric positive definite linear operator Q_o is defined, then the state observation matching function is given as:

  observationMatch(x,p) = 0.5 * diff_o(x,p)^T * Q_o * diff_o(x,p)

and if Q_o is not defined, then the state observation matching function is given as:

  observationMatch(x,p) = (0.5/no) * diff_o(x,p)^T * diff_o(x,p)

where

diff_o(x,p) = o(x,p) - ot

where ot is the target vector for some observation (see below) and p is one of the parameter subvectors supported by the underlying model.

The observation function o(x,p) can be the state vector itself o(x,p) = x for obs_idx < 0, or can be any of the built-in response functions o(x,p) = g(obs_idx)(x,p) when 0 <= obs_idx < getUnderlyingModel()->Ng().

The parameter regularization function also has one of two definitions.

If a symmetric positive definite linear operator Q_p is defined, then the parameter regularization function is given as:

  parameterRegularization(p) = 0.5 * diff_p(p)^T * Q_p * diff_p(p)

and if Q_p is not defined, then the parameter regularization function is given as:

  parameterRegularization(p) = (0.5/np) * diff_p(p)^T * diff_p(p)

where

diff_p(p) = p - pt

where pt is a nomial parameter vector for which violations are penalized against.

Since this decorator class adds a response function, then this->Ng() == getUnderlyingModel()->Ng() + 1.

Let's consider the derivatives of this inverse function.

The first derivatives are given by:

 d(g_inv)/d(x) = observationMultiplier * d(observationMatch)/d(x)

 d(g_inv)/d(p) = observationMultiplier * d(observationMatch)/d(p)
               + parameterMultiplier * d(parameterRegularization)/d(p)

where the derivatives of observationMatch(x,p) and parameterRegularization(p) are given by:

                             /  diff_o(x,p)^T * Q_o * d(o)/d(x) : Q_o defined
  d(observationMatch)/d(x) = |
                             \  (1/no) * diff_o(x,p)^T * d(o)/d(x) : Q_o not defined


                             /  diff_o(x,p)^T * Q_o * d(o)/d(p) : Q_o defined
  d(observationMatch)/d(p) = |
                             \  (1/no) * diff_o(x,p)^T * d(o)/d(p) : Q_o not defined


                                    /  diff_p(p)^T * Q_p : Q_p defined
  d(parameterRegularization)/d(p) = |
                                    \  (1/np) * diff_p(p)^T : Q_p not defined

Of course when obs_idx < -1 where o(x,p) = x then d(o)/d(x) = I and d(o)/d(p) = 0 which also gives d(observationMatch)/d(p) = 0.

Also, we typically want these derivatives in gradient form which gives:

 d(g_inv)/d(x)^T = observationMultiplier * d(observationMatch)/d(x)^T


 d(g_inv)/d(p)^T = observationMultiplier * d(observationMatch)/d(p)^T
                 + parameterMultiplier * d(parameterRegularization)/d(p)^T


                               /  d(o)/d(x)^T * Q_o * diff_o(x,p) : Q_o defined
  d(observationMatch)/d(x)^T = |
                               \  (1/no) * d(o)/d(x)^T * diff_o(x,p) : Q_o not defined


                               /  d(o)/d(p)^T * Q_o * diff_o(x,p) : Q_o defined
  d(observationMatch)/d(p)^T = |
                               \  (1/no) * d(o)/d(p)^T * diff_o(x,p) : Q_o not defined


                                      /  Q_p * diff_p(p) : Q_p defined
  d(parameterRegularization)/d(p)^T = |
                                      \  (1/np) * diff_p(p) : Q_p not defined

When obs_idx >= 0, this implementation currently requires that (DoDx^T) and (DoDp^T) be computed and returned by the underlying model as multi-vector objects. In the future, we really only need the action of DoDx^T and DoDp^T onto vectors as shown above.

Another feature supported by this class is the ability to tack on parameter regularization to an existing response function. This mode is enabled by setting the parameter "Observation Pass Through" to true. This results in the observation matching term to be defined as:

  observationMatch(x,p) = o(x,p)

and has the derivatives:

  d(observationMatch)/d(x)^T = d(o)/d(x)^T


  d(observationMatch)/d(p)^T = d(o)/d(p)^ T

Everything else about the above discussion.

Note: In this case, of course, the observation response function must have dimension 1.

Definition at line 263 of file Thyra_DefaultInverseModelEvaluator.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Thyra::DefaultInverseModelEvaluator< Scalar >::DefaultInverseModelEvaluator ( )

Definition at line 541 of file Thyra_DefaultInverseModelEvaluator.hpp.

Member Function Documentation

template<class Scalar >
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS ( VectorBase< Scalar >  ,
observationTarget   
)

Observation target vector ot.

template<class Scalar >
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS ( VectorBase< Scalar >  ,
parameterBase   
)

Parameter base vector pt.

template<class Scalar >
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS ( LinearOpBase< Scalar >  ,
observationMatchWeightingOp   
)

Observation match weighting operator Q_o.

template<class Scalar >
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS ( LinearOpBase< Scalar >  ,
parameterRegularizationWeightingOp   
)

Parameter regulization weighting operator Q_p.

template<class Scalar >
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_NONCONST_COMPOSITION_MEMBERS ( MultiVectorFileIOBase< Scalar >  ,
observationTargetIO   
)

MultiVectorFileIOBase object used to read the observation target vector ot as directed by the parameter list.

template<class Scalar >
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_NONCONST_COMPOSITION_MEMBERS ( MultiVectorFileIOBase< Scalar >  ,
parameterBaseIO   
)

MultiVectorFileIOBase object used to read the parameter base vector pt as directed by the parameter list.

template<class Scalar >
void Thyra::DefaultInverseModelEvaluator< Scalar >::initialize ( const RCP< ModelEvaluator< Scalar > > &  thyraModel)

Definition at line 552 of file Thyra_DefaultInverseModelEvaluator.hpp.

template<class Scalar >
void Thyra::DefaultInverseModelEvaluator< Scalar >::uninitialize ( RCP< ModelEvaluator< Scalar > > *  thyraModel)

Definition at line 565 of file Thyra_DefaultInverseModelEvaluator.hpp.

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

Note that observationTargetIO() and parameterBaseIO() must be set before calling this function in order to use the parameter sublist to read in the vectors observationTarget() and parameterBase().

Implements Teuchos::ParameterListAcceptor.

Definition at line 578 of file Thyra_DefaultInverseModelEvaluator.hpp.

template<class Scalar >
RCP< Teuchos::ParameterList > Thyra::DefaultInverseModelEvaluator< Scalar >::getNonconstParameterList ( )
virtual
template<class Scalar >
RCP< Teuchos::ParameterList > Thyra::DefaultInverseModelEvaluator< Scalar >::unsetParameterList ( )
virtual
template<class Scalar >
RCP< const Teuchos::ParameterList > Thyra::DefaultInverseModelEvaluator< Scalar >::getParameterList ( ) const
virtual

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 689 of file Thyra_DefaultInverseModelEvaluator.hpp.

template<class Scalar >
RCP< const Teuchos::ParameterList > Thyra::DefaultInverseModelEvaluator< Scalar >::getValidParameters ( ) const
virtual

Note that observationTargetIO() and parameterBaseIO() must be set before calling this function in order to have the sublists added that will allow the vectors observationTarget() and parameterBase() to be read in latter when the parameter list is set..

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 697 of file Thyra_DefaultInverseModelEvaluator.hpp.

template<class Scalar >
std::string Thyra::DefaultInverseModelEvaluator< Scalar >::description ( ) const
virtual

Reimplemented from Teuchos::Describable.

Definition at line 784 of file Thyra_DefaultInverseModelEvaluator.hpp.

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

Friends And Related Function Documentation

template<class Scalar >
RCP< DefaultInverseModelEvaluator< Scalar > > defaultInverseModelEvaluator ( const RCP< ModelEvaluator< Scalar > > &  thyraModel)
related

Non-member constructor.

Definition at line 441 of file Thyra_DefaultInverseModelEvaluator.hpp.


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