EpetraExt::ModelEvaluator object and wraps it as a Thyra::ModelEvaluator object.
More...
#include <Thyra_EpetraModelEvaluator.hpp>
Inheritance diagram for Thyra::EpetraModelEvaluator:
Constructors/initializers/accessors/utilities. | |
| EpetraModelEvaluator () | |
| | |
| EpetraModelEvaluator (const RCP< const EpetraExt::ModelEvaluator > &epetraModel, const RCP< LinearOpWithSolveFactoryBase< double > > &W_factory) | |
| | |
| void | initialize (const RCP< const EpetraExt::ModelEvaluator > &epetraModel, const RCP< LinearOpWithSolveFactoryBase< double > > &W_factory) |
| | |
| RCP< const EpetraExt::ModelEvaluator > | getEpetraModel () const |
| | |
| void | setNominalValues (const ModelEvaluatorBase::InArgs< double > &nominalValues) |
| Set the nominal values. | |
| void | setStateVariableScalingVec (const RCP< const Epetra_Vector > &stateVariableScalingVec) |
Set the state variable scaling vector s_x (see above). | |
| RCP< const Epetra_Vector > | getStateVariableInvScalingVec () const |
Get the state variable scaling vector s_x (see above). | |
| RCP< const Epetra_Vector > | getStateVariableScalingVec () const |
Get the inverse state variable scaling vector inv_s_x (see above). | |
| void | setStateFunctionScalingVec (const RCP< const Epetra_Vector > &stateFunctionScalingVec) |
Set the state function scaling vector s_f (see above). | |
| RCP< const Epetra_Vector > | getStateFunctionScalingVec () const |
Get the state function scaling vector s_f (see above). | |
| void | uninitialize (RCP< const EpetraExt::ModelEvaluator > *epetraModel=NULL, RCP< LinearOpWithSolveFactoryBase< double > > *W_factory=NULL) |
| | |
| const ModelEvaluatorBase::InArgs< double > & | getFinalPoint () const |
| | |
| bool | finalPointWasSolved () const |
| | |
Public functions overridden from Teuchos::Describable. | |
| std::string | description () const |
| | |
Overridden from ParameterListAcceptor | |
| void | setParameterList (RCP< Teuchos::ParameterList > const ¶mList) |
| | |
| RCP< Teuchos::ParameterList > | getParameterList () |
| | |
| RCP< Teuchos::ParameterList > | unsetParameterList () |
| | |
| RCP< const Teuchos::ParameterList > | getParameterList () const |
| | |
| RCP< const Teuchos::ParameterList > | getValidParameters () const |
| | |
Public functions overridden from ModelEvaulator. | |
| int | Np () const |
| | |
| int | Ng () const |
| | |
| RCP< const VectorSpaceBase< double > > | get_x_space () const |
| | |
| RCP< const VectorSpaceBase< double > > | get_f_space () const |
| | |
| RCP< const VectorSpaceBase< double > > | get_p_space (int l) const |
| | |
| RCP< const Teuchos::Array< std::string > > | get_p_names (int l) const |
| | |
| RCP< const VectorSpaceBase< double > > | get_g_space (int j) const |
| | |
| ModelEvaluatorBase::InArgs< double > | getNominalValues () const |
| | |
| ModelEvaluatorBase::InArgs< double > | getLowerBounds () const |
| | |
| ModelEvaluatorBase::InArgs< double > | getUpperBounds () const |
| | |
| RCP< LinearOpWithSolveBase< double > > | create_W () const |
| | |
| RCP< LinearOpBase< double > > | create_W_op () const |
| | |
| ModelEvaluatorBase::InArgs< double > | createInArgs () const |
| | |
| void | reportFinalPoint (const ModelEvaluatorBase::InArgs< double > &finalPoint, const bool wasSolved) |
| | |
Public Types | |
| enum | EStateFunctionScaling |
Related Functions | |
| (Note that these are not member functions.) | |
| RCP< EpetraModelEvaluator > | epetraModelEvaluator (const RCP< const EpetraExt::ModelEvaluator > &epetraModel, const RCP< LinearOpWithSolveFactoryBase< double > > &W_factory) |
| | |
| ModelEvaluatorBase::EDerivativeMultiVectorOrientation | convert (const EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation &mvOrientation) |
| | |
| EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation | convert (const ModelEvaluatorBase::EDerivativeMultiVectorOrientation &mvOrientation) |
| | |
| ModelEvaluatorBase::DerivativeProperties | convert (const EpetraExt::ModelEvaluator::DerivativeProperties &derivativeProperties) |
| | |
| ModelEvaluatorBase::DerivativeSupport | convert (const EpetraExt::ModelEvaluator::DerivativeSupport &derivativeSupport) |
| | |
| EpetraExt::ModelEvaluator::Derivative | convert (const ModelEvaluatorBase::Derivative< double > &derivative, const RCP< const Epetra_Map > &fnc_map, const RCP< const Epetra_Map > &var_map) |
| | |
EpetraExt::ModelEvaluator object and wraps it as a Thyra::ModelEvaluator object.
This class takes care of the basic details of wrapping and unwrapping Epetra from Thyra objects. This class is highly configurable and will be maintained and modified in the future as the basic link between the Epetra world and the Thyra world for nonlinear models and nonlinear algorithms.
The scaling for the state function can be set manually using setStateFunctionScalingVec() or can be computed automatically using the parameter "State Function Scaling" (see documentation output from this->getValidParameters()->print(...)) in the input parameter list set by setParameterList(). Reguardless of how the state function scaling is computed, it will compute a positive vector s_f that defines a diagonal matrix S_f = diag(s_f) that transforms the state function:
f(...) = S_f * f_orig(...)
where f_orig(...) is the original state function as computed by the underlying EpetraExt::ModelEvaluator object and f(...) is the state function as computed by evalModel().
The scaling for the state variables must be set manually using Thyra::setStateVariableScalingVec(). The vector that is set s_x>/tt> defines a diagonal scaling matrix S_x = diag(s_x) that transforms the variables as:
x = S_x * x_orig
where x_orig is the original unscaled state variable vector as defined by the underlying EpetraExt::ModelEvaluator object and x is the scaled state varaible vector as returned from getNominalValues() and as accepted by evalModel(). Note that when the scaled variables x are passed into evalModel that they are unscaled as:
x_orig = inv(S_x) * x
where inv(S_x) is the inverse of the diagonals of S_x which is stored as a vector inv_s_x.
Note how these scalings affect the state function:
f(x,...) = S_f * f_orig( inv(S_x)*x...)
which as the state/state Jacobian:
W = d(f)/d(x) = S_f * d(f_orig)/d(x_orig) * inv(S_x)
Currently, this class does not handle scalings of the parameters p(l) or of the auxilary response functions g(j)(...).
The state varaible and state function scaling gives the following scaled quantities:
f = S_f * f_orig
W = S_f * W_orig * inv(S_x)
DfDp(l) = S_f * DfDp_orig(l)
g(j) = g_orig(j)
DgDx_dot(j) = DgDx_dot_orig(j) * inv(S_x)
DgDx(j) = DgDx_orig(j) * inv(S_x)
DgDp(j,l) = DgDp_orig(j,l)
Since the scaling is done explicitly, the client never even sees the orginal scaling and the linear solver (and contained preconditioner) are computed from the scaled W shown above.
ToDo: Describe how scaling affects the Hessian-vector products an how you just need to scale the Lagrange mutipliers as:
u^T * f(...) = u^T * (S_f * f_orig(...)) = u_f^T * f_orig(...)
where u_f = S_f * u.
ToDo: Finish documentation!
Definition at line 160 of file Thyra_EpetraModelEvaluator.hpp.
| Thyra::EpetraModelEvaluator::EpetraModelEvaluator | ( | ) |
| Thyra::EpetraModelEvaluator::EpetraModelEvaluator | ( | const RCP< const EpetraExt::ModelEvaluator > & | epetraModel, | |
| const RCP< LinearOpWithSolveFactoryBase< double > > & | W_factory | |||
| ) |
| void Thyra::EpetraModelEvaluator::initialize | ( | const RCP< const EpetraExt::ModelEvaluator > & | epetraModel, | |
| const RCP< LinearOpWithSolveFactoryBase< double > > & | W_factory | |||
| ) |
| RCP< const EpetraExt::ModelEvaluator > Thyra::EpetraModelEvaluator::getEpetraModel | ( | ) | const |
| void Thyra::EpetraModelEvaluator::setNominalValues | ( | const ModelEvaluatorBase::InArgs< double > & | nominalValues | ) |
Set the nominal values.
Warning, if scaling is being used, these must be according to the scaled values, not the original unscaled values.
Definition at line 186 of file Thyra_EpetraModelEvaluator.cpp.
| void Thyra::EpetraModelEvaluator::setStateVariableScalingVec | ( | const RCP< const Epetra_Vector > & | stateVariableScalingVec | ) |
Set the state variable scaling vector s_x (see above).
This function must be called after intialize() or the constructur in order to set the scaling vector correctly!
ToDo: Move this into an external strategy class object!
Definition at line 195 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const Epetra_Vector > Thyra::EpetraModelEvaluator::getStateVariableInvScalingVec | ( | ) | const |
Get the state variable scaling vector s_x (see above).
Definition at line 217 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const Epetra_Vector > Thyra::EpetraModelEvaluator::getStateVariableScalingVec | ( | ) | const |
Get the inverse state variable scaling vector inv_s_x (see above).
Definition at line 210 of file Thyra_EpetraModelEvaluator.cpp.
| void Thyra::EpetraModelEvaluator::setStateFunctionScalingVec | ( | const RCP< const Epetra_Vector > & | stateFunctionScalingVec | ) |
Set the state function scaling vector s_f (see above).
Definition at line 224 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const Epetra_Vector > Thyra::EpetraModelEvaluator::getStateFunctionScalingVec | ( | ) | const |
Get the state function scaling vector s_f (see above).
Definition at line 233 of file Thyra_EpetraModelEvaluator.cpp.
| void Thyra::EpetraModelEvaluator::uninitialize | ( | RCP< const EpetraExt::ModelEvaluator > * | epetraModel = NULL, |
|
| RCP< LinearOpWithSolveFactoryBase< double > > * | W_factory = NULL | |||
| ) |
| const ModelEvaluatorBase::InArgs< double > & Thyra::EpetraModelEvaluator::getFinalPoint | ( | ) | const |
| bool Thyra::EpetraModelEvaluator::finalPointWasSolved | ( | ) | const |
| std::string Thyra::EpetraModelEvaluator::description | ( | ) | const [virtual] |
Reimplemented from Teuchos::Describable.
Definition at line 271 of file Thyra_EpetraModelEvaluator.cpp.
| void Thyra::EpetraModelEvaluator::setParameterList | ( | RCP< Teuchos::ParameterList > const & | paramList | ) |
| RCP< Teuchos::ParameterList > Thyra::EpetraModelEvaluator::getParameterList | ( | ) | [virtual] |
Implements Teuchos::ParameterListAcceptor.
Definition at line 314 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< Teuchos::ParameterList > Thyra::EpetraModelEvaluator::unsetParameterList | ( | ) | [virtual] |
Implements Teuchos::ParameterListAcceptor.
Definition at line 321 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const Teuchos::ParameterList > Thyra::EpetraModelEvaluator::getParameterList | ( | ) | const [virtual] |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 330 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const Teuchos::ParameterList > Thyra::EpetraModelEvaluator::getValidParameters | ( | ) | const [virtual] |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 337 of file Thyra_EpetraModelEvaluator.cpp.
| int Thyra::EpetraModelEvaluator::Np | ( | ) | const [virtual] |
Reimplemented from Thyra::ModelEvaluatorDefaultBase< double >.
Definition at line 384 of file Thyra_EpetraModelEvaluator.cpp.
| int Thyra::EpetraModelEvaluator::Ng | ( | ) | const [virtual] |
Reimplemented from Thyra::ModelEvaluatorDefaultBase< double >.
Definition at line 390 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_x_space | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 397 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_f_space | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 404 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_p_space | ( | int | l | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 411 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const Teuchos::Array< std::string > > Thyra::EpetraModelEvaluator::get_p_names | ( | int | l | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 421 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_g_space | ( | int | j | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 431 of file Thyra_EpetraModelEvaluator.cpp.
| ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::getNominalValues | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 439 of file Thyra_EpetraModelEvaluator.cpp.
| ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::getLowerBounds | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 447 of file Thyra_EpetraModelEvaluator.cpp.
| ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::getUpperBounds | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 455 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< LinearOpWithSolveBase< double > > Thyra::EpetraModelEvaluator::create_W | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 463 of file Thyra_EpetraModelEvaluator.cpp.
| RCP< LinearOpBase< double > > Thyra::EpetraModelEvaluator::create_W_op | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 478 of file Thyra_EpetraModelEvaluator.cpp.
| ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::createInArgs | ( | ) | const [virtual] |
Implements Thyra::ModelEvaluator< Scalar >.
Definition at line 484 of file Thyra_EpetraModelEvaluator.cpp.
| void Thyra::EpetraModelEvaluator::reportFinalPoint | ( | const ModelEvaluatorBase::InArgs< double > & | finalPoint, | |
| const bool | wasSolved | |||
| ) |
| RCP< EpetraModelEvaluator > epetraModelEvaluator | ( | const RCP< const EpetraExt::ModelEvaluator > & | epetraModel, | |
| const RCP< LinearOpWithSolveFactoryBase< double > > & | W_factory | |||
| ) | [related] |
| ModelEvaluatorBase::EDerivativeMultiVectorOrientation convert | ( | const EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation & | mvOrientation | ) | [related] |
| EpetraExt::ModelEvaluator::EDerivativeMultiVectorOrientation convert | ( | const ModelEvaluatorBase::EDerivativeMultiVectorOrientation & | mvOrientation | ) | [related] |
| ModelEvaluatorBase::DerivativeProperties convert | ( | const EpetraExt::ModelEvaluator::DerivativeProperties & | derivativeProperties | ) | [related] |
| ModelEvaluatorBase::DerivativeSupport convert | ( | const EpetraExt::ModelEvaluator::DerivativeSupport & | derivativeSupport | ) | [related] |
| EpetraExt::ModelEvaluator::Derivative convert | ( | const ModelEvaluatorBase::Derivative< double > & | derivative, | |
| const RCP< const Epetra_Map > & | fnc_map, | |||
| const RCP< const Epetra_Map > & | var_map | |||
| ) | [related] |
1.4.7