Thyra
Version of the Day
|
Decorator class that wraps any ModelEvaluator object and lumps parameters together using a linear basis matrix. More...
#include <Thyra_DefaultLumpedParameterModelEvaluator.hpp>
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... | |
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. | |
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 ¶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 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) |
Decorator class that wraps any ModelEvaluator object and lumps parameters together using a linear basis matrix.
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).
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.
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.
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 208 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::DefaultLumpedParameterModelEvaluator | ( | ) |
Definition at line 495 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
void Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::initialize | ( | const RCP< ModelEvaluator< Scalar > > & | thyraModel | ) |
Definition at line 509 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
void Thyra::DefaultLumpedParameterModelEvaluator< Scalar >::uninitialize | ( | RCP< ModelEvaluator< Scalar > > * | thyraModel | ) |
Definition at line 520 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 535 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 555 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 604 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 612 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 622 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 630 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 676 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 687 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 698 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 707 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 716 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
virtual |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 724 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.
|
related |
Non-member constructor.
Definition at line 401 of file Thyra_DefaultLumpedParameterModelEvaluator.hpp.