Thyra
Version of the Day
|
This class wraps any ModelEvaluator object and adds a simple, but fairly general, inverse response function. More...
#include <Thyra_DefaultInverseModelEvaluator.hpp>
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... | |
Related Functions inherited from Thyra::ModelEvaluatorDefaultBase< Scalar > | |
template<class Scalar > | |
RCP < ModelEvaluatorBase::InArgs < Scalar > > | clone (const ModelEvaluatorBase::InArgs< Scalar > &inArgs) |
Create a clone of an InArgs object. More... | |
template<class Scalar > | |
ModelEvaluatorBase::Derivative < Scalar > | derivativeGradient (const RCP< MultiVectorBase< Scalar > > &grad) |
template<class Scalar > | |
ModelEvaluatorBase::DerivativeMultiVector < Scalar > | create_DfDp_mv (const ModelEvaluator< Scalar > &model, int l, ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
template<class Scalar > | |
ModelEvaluatorBase::DerivativeMultiVector < Scalar > | create_DgDx_dot_mv (const ModelEvaluator< Scalar > &model, int j, ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
template<class Scalar > | |
ModelEvaluatorBase::DerivativeMultiVector < Scalar > | create_DgDx_mv (const ModelEvaluator< Scalar > &model, int j, ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
template<class Scalar > | |
ModelEvaluatorBase::DerivativeMultiVector < Scalar > | create_DgDp_mv (const ModelEvaluator< Scalar > &model, int j, int l, ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
template<class Scalar > | |
ModelEvaluatorBase::DerivativeMultiVector < Scalar > | get_dmv (const ModelEvaluatorBase::Derivative< Scalar > &deriv, const std::string &derivName) |
template<class Scalar > | |
RCP< MultiVectorBase< Scalar > > | get_mv (const ModelEvaluatorBase::Derivative< Scalar > &deriv, const std::string &derivName, ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
template<class Scalar > | |
void | assertDerivSpaces (const std::string &modelEvalDescription, const ModelEvaluatorBase::Derivative< Scalar > &deriv, const std::string &deriv_name, const VectorSpaceBase< Scalar > &fnc_space, const std::string &fnc_space_name, const VectorSpaceBase< Scalar > &var_space, const std::string &var_space_name) |
Assert that that Thyra objects imbedded in a Derivative object matches its function and variable spaces. More... | |
template<class Scalar > | |
void | assertInArgsOutArgsSetup (const std::string &modelEvalDescription, const ModelEvaluatorBase::InArgs< Scalar > &inArgs, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs) |
Assert that an InArgs and OutArgs object are setup consistently. More... | |
template<class Scalar > | |
void | assertInArgsEvalObjects (const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::InArgs< Scalar > &inArgs) |
Assert that the objects in an InArgs object match a given model. More... | |
template<class Scalar > | |
void | assertOutArgsEvalObjects (const ModelEvaluator< Scalar > &model, const ModelEvaluatorBase::OutArgs< Scalar > &outArgs, const ModelEvaluatorBase::InArgs< Scalar > *inArgs=0) |
Assert that the objects in an OutArgs object match a given model. More... | |
Related Functions inherited from Thyra::ModelEvaluatorBase | |
std::string | toString (ModelEvaluatorBase::EInArgsMembers) |
std::string | toString (ModelEvaluatorBase::EOutArgsMembers) |
std::string | toString (ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
ModelEvaluatorBase::EDerivativeMultiVectorOrientation | getOtherDerivativeMultiVectorOrientation (ModelEvaluatorBase::EDerivativeMultiVectorOrientation orientation) |
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 ¶mList) |
RCP< Teuchos::ParameterList > | getNonconstParameterList () |
RCP< Teuchos::ParameterList > | unsetParameterList () |
RCP< const Teuchos::ParameterList > | getParameterList () const |
RCP< const Teuchos::ParameterList > | getValidParameters () 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 |
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 231 of file Thyra_DefaultInverseModelEvaluator.hpp.
Thyra::DefaultInverseModelEvaluator< Scalar >::DefaultInverseModelEvaluator | ( | ) |
Definition at line 509 of file Thyra_DefaultInverseModelEvaluator.hpp.
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS | ( | VectorBase< Scalar > | , |
observationTarget | |||
) |
Observation target vector ot
.
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS | ( | VectorBase< Scalar > | , |
parameterBase | |||
) |
Parameter base vector pt
.
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS | ( | LinearOpBase< Scalar > | , |
observationMatchWeightingOp | |||
) |
Observation match weighting operator Q_o
.
Thyra::DefaultInverseModelEvaluator< Scalar >::STANDARD_CONST_COMPOSITION_MEMBERS | ( | LinearOpBase< Scalar > | , |
parameterRegularizationWeightingOp | |||
) |
Parameter regulization weighting operator Q_p
.
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.
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.
void Thyra::DefaultInverseModelEvaluator< Scalar >::initialize | ( | const RCP< ModelEvaluator< Scalar > > & | thyraModel | ) |
Definition at line 520 of file Thyra_DefaultInverseModelEvaluator.hpp.
void Thyra::DefaultInverseModelEvaluator< Scalar >::uninitialize | ( | RCP< ModelEvaluator< Scalar > > * | thyraModel | ) |
Definition at line 533 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
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 546 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 639 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 647 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
virtual |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 657 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
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 665 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 752 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 716 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 728 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 740 of file Thyra_DefaultInverseModelEvaluator.hpp.
|
related |
Non-member constructor.
Definition at line 409 of file Thyra_DefaultInverseModelEvaluator.hpp.