Thyra Package Browser (Single Doxygen Collection)
Version of the Day
|
Concrete Adapter subclass that takes an EpetraExt::ModelEvaluator
object and wraps it as a Thyra::ModelEvaluator
object.
More...
#include <Thyra_EpetraModelEvaluator.hpp>
Public Types | |
enum | EStateFunctionScaling { STATE_FUNC_SCALING_NONE, STATE_FUNC_SCALING_ROW_SUM } |
Private Types | |
typedef Teuchos::Array< RCP < const Epetra_Map > > | p_map_t |
typedef Teuchos::Array< RCP < const Epetra_Map > > | g_map_t |
typedef std::vector< bool > | p_map_is_local_t |
typedef std::vector< bool > | g_map_is_local_t |
typedef Teuchos::Array< RCP < const VectorSpaceBase < double > > > | p_space_t |
typedef Teuchos::Array< RCP < const VectorSpaceBase < double > > > | g_space_t |
Private Member Functions | |
void | convertInArgsFromEpetraToThyra (const EpetraExt::ModelEvaluator::InArgs &epetraInArgs, ModelEvaluatorBase::InArgs< double > *inArgs) const |
void | convertInArgsFromThyraToEpetra (const ModelEvaluatorBase::InArgs< double > &inArgs, EpetraExt::ModelEvaluator::InArgs *epetraInArgs) const |
void | convertOutArgsFromThyraToEpetra (const ModelEvaluatorBase::OutArgs< double > &outArgs, EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs, RCP< LinearOpBase< double > > *W_op, RCP< EpetraLinearOp > *efwdW, RCP< Epetra_Operator > *eW) const |
void | preEvalScalingSetup (EpetraExt::ModelEvaluator::InArgs *epetraInArgs, EpetraExt::ModelEvaluator::OutArgs *epetraUnscaledOutArgs, const RCP< Teuchos::FancyOStream > &out, const Teuchos::EVerbosityLevel verbLevel) const |
void | postEvalScalingSetup (const EpetraExt::ModelEvaluator::OutArgs &epetraUnscaledOutArgs, const RCP< Teuchos::FancyOStream > &out, const Teuchos::EVerbosityLevel verbLevel) const |
void | finishConvertingOutArgsFromEpetraToThyra (const EpetraExt::ModelEvaluator::OutArgs &epetraOutArgs, RCP< LinearOpBase< double > > &W_op, RCP< EpetraLinearOp > &efwdW, RCP< Epetra_Operator > &eW, const ModelEvaluatorBase::OutArgs< double > &outArgs) const |
void | updateNominalValuesAndBounds () const |
void | updateInArgsOutArgs () const |
RCP< EpetraLinearOp > | create_epetra_W_op () const |
Private Attributes | |
RCP< const EpetraExt::ModelEvaluator > | epetraModel_ |
RCP< Teuchos::ParameterList > | paramList_ |
RCP < LinearOpWithSolveFactoryBase < double > > | W_factory_ |
RCP< const Epetra_Map > | x_map_ |
p_map_t | p_map_ |
g_map_t | g_map_ |
p_map_is_local_t | p_map_is_local_ |
p_map_is_local_t | g_map_is_local_ |
RCP< const Epetra_Map > | f_map_ |
RCP< const VectorSpaceBase < double > > | x_space_ |
p_space_t | p_space_ |
RCP< const VectorSpaceBase < double > > | f_space_ |
g_space_t | g_space_ |
ModelEvaluatorBase::InArgs < double > | nominalValues_ |
ModelEvaluatorBase::InArgs < double > | lowerBounds_ |
ModelEvaluatorBase::InArgs < double > | upperBounds_ |
bool | nominalValuesAndBoundsAreUpdated_ |
ModelEvaluatorBase::InArgs < double > | finalPoint_ |
EStateFunctionScaling | stateFunctionScaling_ |
RCP< const Epetra_Vector > | stateFunctionScalingVec_ |
RCP< const Epetra_Vector > | stateVariableScalingVec_ |
RCP< const Epetra_Vector > | invStateVariableScalingVec_ |
EpetraExt::ModelEvaluator::InArgs | epetraInArgsScaling_ |
EpetraExt::ModelEvaluator::OutArgs | epetraOutArgsScaling_ |
RCP< Epetra_Vector > | x_unscaled_ |
RCP< Epetra_Vector > | x_dot_unscaled_ |
ModelEvaluatorBase::InArgs < double > | prototypeInArgs_ |
ModelEvaluatorBase::OutArgs < double > | prototypeOutArgs_ |
bool | currentInArgsOutArgs_ |
bool | finalPointWasSolved_ |
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) |
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. More... | |
void | setStateVariableScalingVec (const RCP< const Epetra_Vector > &stateVariableScalingVec) |
Set the state variable scaling vector s_x (see above). More... | |
RCP< const Epetra_Vector > | getStateVariableInvScalingVec () const |
Get the state variable scaling vector s_x (see above). More... | |
RCP< const Epetra_Vector > | getStateVariableScalingVec () const |
Get the inverse state variable scaling vector inv_s_x (see above). More... | |
void | setStateFunctionScalingVec (const RCP< const Epetra_Vector > &stateFunctionScalingVec) |
Set the state function scaling vector s_f (see above). More... | |
RCP< const Epetra_Vector > | getStateFunctionScalingVec () const |
Get the state function scaling vector s_f (see above). More... | |
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 > | getNonconstParameterList () |
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 |
Teuchos::ArrayView< const std::string > | get_g_names (int j) const |
ModelEvaluatorBase::InArgs < double > | getNominalValues () const |
ModelEvaluatorBase::InArgs < double > | getLowerBounds () const |
ModelEvaluatorBase::InArgs < double > | getUpperBounds () const |
RCP< LinearOpBase< double > > | create_W_op () const |
RCP< PreconditionerBase< double > > | create_W_prec () const |
Returns null currently. More... | |
RCP< const LinearOpWithSolveFactoryBase < double > > | get_W_factory () const |
ModelEvaluatorBase::InArgs < double > | createInArgs () const |
void | reportFinalPoint (const ModelEvaluatorBase::InArgs< double > &finalPoint, const bool wasSolved) |
Private functions overridden from ModelEvaulatorDefaultBase. | |
RCP< LinearOpBase< double > > | create_DfDp_op_impl (int l) const |
RCP< LinearOpBase< double > > | create_DgDx_dot_op_impl (int j) const |
RCP< LinearOpBase< double > > | create_DgDx_op_impl (int j) const |
RCP< LinearOpBase< double > > | create_DgDp_op_impl (int j, int l) const |
ModelEvaluatorBase::OutArgs < double > | createOutArgsImpl () const |
void | evalModelImpl (const ModelEvaluatorBase::InArgs< double > &inArgs, const ModelEvaluatorBase::OutArgs< double > &outArgs) const |
Additional Inherited Members |
Concrete Adapter subclass that takes an 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.
This class can handle scaling of the state function f(...) and of the state variables x and all of the affected derivatives.
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 175 of file Thyra_EpetraModelEvaluator.hpp.
|
private |
Definition at line 352 of file Thyra_EpetraModelEvaluator.hpp.
|
private |
Definition at line 353 of file Thyra_EpetraModelEvaluator.hpp.
|
private |
Definition at line 354 of file Thyra_EpetraModelEvaluator.hpp.
|
private |
Definition at line 355 of file Thyra_EpetraModelEvaluator.hpp.
|
private |
Definition at line 358 of file Thyra_EpetraModelEvaluator.hpp.
|
private |
Definition at line 360 of file Thyra_EpetraModelEvaluator.hpp.
Enumerator | |
---|---|
STATE_FUNC_SCALING_NONE | |
STATE_FUNC_SCALING_ROW_SUM |
Definition at line 322 of file Thyra_EpetraModelEvaluator.hpp.
Thyra::EpetraModelEvaluator::EpetraModelEvaluator | ( | ) |
Definition at line 123 of file Thyra_EpetraModelEvaluator.cpp.
Thyra::EpetraModelEvaluator::EpetraModelEvaluator | ( | const RCP< const EpetraExt::ModelEvaluator > & | epetraModel, |
const RCP< LinearOpWithSolveFactoryBase< double > > & | W_factory | ||
) |
Definition at line 129 of file Thyra_EpetraModelEvaluator.cpp.
References initialize().
void Thyra::EpetraModelEvaluator::initialize | ( | const RCP< const EpetraExt::ModelEvaluator > & | epetraModel, |
const RCP< LinearOpWithSolveFactoryBase< double > > & | W_factory | ||
) |
Definition at line 140 of file Thyra_EpetraModelEvaluator.cpp.
References Thyra::create_VectorSpace(), currentInArgsOutArgs_, epetraInArgsScaling_, epetraModel_, epetraOutArgsScaling_, f_map_, f_space_, finalPointWasSolved_, g_map_, g_map_is_local_, g_space_, Teuchos::implicit_cast(), is_null(), nominalValuesAndBoundsAreUpdated_, Teuchos::null, p_map_, p_map_is_local_, p_space_, Teuchos::Array< T >::resize(), Teuchos::Array< T >::size(), stateFunctionScalingVec_, TEUCHOS_TEST_FOR_EXCEPTION, W_factory_, x_map_, and x_space_.
Referenced by EpetraModelEvaluator().
RCP< const EpetraExt::ModelEvaluator > Thyra::EpetraModelEvaluator::getEpetraModel | ( | ) | const |
Definition at line 197 of file Thyra_EpetraModelEvaluator.cpp.
References epetraModel_.
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 203 of file Thyra_EpetraModelEvaluator.cpp.
References nominalValues_.
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 212 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::RCP< T >::assert_not_null(), createInArgs(), invStateVariableScalingVec_, nominalValuesAndBoundsAreUpdated_, Teuchos::null, stateVariableScalingVec_, and TEUCHOS_TEST_FOR_EXCEPT.
RCP< const Epetra_Vector > Thyra::EpetraModelEvaluator::getStateVariableInvScalingVec | ( | ) | const |
Get the state variable scaling vector s_x
(see above).
Definition at line 234 of file Thyra_EpetraModelEvaluator.cpp.
References invStateVariableScalingVec_, and updateNominalValuesAndBounds().
RCP< const Epetra_Vector > Thyra::EpetraModelEvaluator::getStateVariableScalingVec | ( | ) | const |
Get the inverse state variable scaling vector inv_s_x
(see above).
Definition at line 227 of file Thyra_EpetraModelEvaluator.cpp.
References stateVariableScalingVec_.
void Thyra::EpetraModelEvaluator::setStateFunctionScalingVec | ( | const RCP< const Epetra_Vector > & | stateFunctionScalingVec | ) |
Set the state function scaling vector s_f
(see above).
Definition at line 241 of file Thyra_EpetraModelEvaluator.cpp.
References stateFunctionScalingVec_.
RCP< const Epetra_Vector > Thyra::EpetraModelEvaluator::getStateFunctionScalingVec | ( | ) | const |
Get the state function scaling vector s_f
(see above).
Definition at line 250 of file Thyra_EpetraModelEvaluator.cpp.
References stateFunctionScalingVec_.
void Thyra::EpetraModelEvaluator::uninitialize | ( | RCP< const EpetraExt::ModelEvaluator > * | epetraModel = NULL , |
RCP< LinearOpWithSolveFactoryBase< double > > * | W_factory = NULL |
||
) |
Definition at line 256 of file Thyra_EpetraModelEvaluator.cpp.
References currentInArgsOutArgs_, epetraModel_, invStateVariableScalingVec_, Teuchos::null, stateFunctionScalingVec_, stateVariableScalingVec_, and W_factory_.
const ModelEvaluatorBase::InArgs< double > & Thyra::EpetraModelEvaluator::getFinalPoint | ( | ) | const |
Definition at line 273 of file Thyra_EpetraModelEvaluator.cpp.
References finalPoint_.
bool Thyra::EpetraModelEvaluator::finalPointWasSolved | ( | ) | const |
Definition at line 279 of file Thyra_EpetraModelEvaluator.cpp.
References finalPointWasSolved_.
std::string Thyra::EpetraModelEvaluator::description | ( | ) | const |
Definition at line 288 of file Thyra_EpetraModelEvaluator.cpp.
References epetraModel_, Teuchos::RCP< T >::get(), and W_factory_.
Referenced by updateInArgsOutArgs().
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 310 of file Thyra_EpetraModelEvaluator.cpp.
References getValidParameters(), is_null(), Teuchos::null, paramList_, stateFunctionScaling_, stateFunctionScalingVec_, TEUCHOS_TEST_FOR_EXCEPT, and Teuchos::ParameterList::validateParameters().
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 331 of file Thyra_EpetraModelEvaluator.cpp.
References paramList_.
|
virtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 338 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::null, and paramList_.
|
virtual |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 347 of file Thyra_EpetraModelEvaluator.cpp.
References paramList_.
|
virtual |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 354 of file Thyra_EpetraModelEvaluator.cpp.
References is_null(), Teuchos::rcp(), rcp(), Teuchos::ParameterList::set(), STATE_FUNC_SCALING_NONE, and STATE_FUNC_SCALING_ROW_SUM.
Referenced by setParameterList().
int Thyra::EpetraModelEvaluator::Np | ( | ) | const |
Definition at line 403 of file Thyra_EpetraModelEvaluator.cpp.
References p_space_, and Teuchos::Array< T >::size().
Referenced by get_p_names(), and get_p_space().
int Thyra::EpetraModelEvaluator::Ng | ( | ) | const |
Definition at line 409 of file Thyra_EpetraModelEvaluator.cpp.
References g_space_, and Teuchos::Array< T >::size().
Referenced by get_g_names(), and get_g_space().
RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_x_space | ( | ) | const |
Definition at line 416 of file Thyra_EpetraModelEvaluator.cpp.
References x_space_.
Referenced by create_epetra_W_op().
RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_f_space | ( | ) | const |
Definition at line 423 of file Thyra_EpetraModelEvaluator.cpp.
References f_space_.
Referenced by create_epetra_W_op().
RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_p_space | ( | int | l | ) | const |
Definition at line 430 of file Thyra_EpetraModelEvaluator.cpp.
References Np(), p_space_, and TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE.
RCP< const Teuchos::Array< std::string > > Thyra::EpetraModelEvaluator::get_p_names | ( | int | l | ) | const |
Definition at line 440 of file Thyra_EpetraModelEvaluator.cpp.
References epetraModel_, Np(), and TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE.
RCP< const VectorSpaceBase< double > > Thyra::EpetraModelEvaluator::get_g_space | ( | int | j | ) | const |
Definition at line 450 of file Thyra_EpetraModelEvaluator.cpp.
References g_space_, Ng(), and TEUCHOS_TEST_FOR_EXCEPT.
Teuchos::ArrayView< const std::string > Thyra::EpetraModelEvaluator::get_g_names | ( | int | j | ) | const |
Definition at line 458 of file Thyra_EpetraModelEvaluator.cpp.
References epetraModel_, Ng(), and TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE.
ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::getNominalValues | ( | ) | const |
Definition at line 468 of file Thyra_EpetraModelEvaluator.cpp.
References nominalValues_, and updateNominalValuesAndBounds().
Referenced by evalModelImpl().
ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::getLowerBounds | ( | ) | const |
Definition at line 476 of file Thyra_EpetraModelEvaluator.cpp.
References lowerBounds_, and updateNominalValuesAndBounds().
ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::getUpperBounds | ( | ) | const |
Definition at line 484 of file Thyra_EpetraModelEvaluator.cpp.
References updateNominalValuesAndBounds(), and upperBounds_.
RCP< LinearOpBase< double > > Thyra::EpetraModelEvaluator::create_W_op | ( | ) | const |
Definition at line 492 of file Thyra_EpetraModelEvaluator.cpp.
References create_epetra_W_op().
RCP< PreconditionerBase< double > > Thyra::EpetraModelEvaluator::create_W_prec | ( | ) | const |
Returns null currently.
Definition at line 499 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::null.
RCP< const LinearOpWithSolveFactoryBase< double > > Thyra::EpetraModelEvaluator::get_W_factory | ( | ) | const |
ModelEvaluatorBase::InArgs< double > Thyra::EpetraModelEvaluator::createInArgs | ( | ) | const |
Definition at line 512 of file Thyra_EpetraModelEvaluator.cpp.
References currentInArgsOutArgs_, prototypeInArgs_, and updateInArgsOutArgs().
Referenced by reportFinalPoint(), setStateVariableScalingVec(), and updateNominalValuesAndBounds().
void Thyra::EpetraModelEvaluator::reportFinalPoint | ( | const ModelEvaluatorBase::InArgs< double > & | finalPoint, |
const bool | wasSolved | ||
) |
Definition at line 520 of file Thyra_EpetraModelEvaluator.cpp.
References createInArgs(), finalPoint_, and finalPointWasSolved_.
|
private |
Definition at line 535 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::null, TEUCHOS_TEST_FOR_EXCEPT, and TEUCHOS_UNREACHABLE_RETURN.
|
private |
Definition at line 543 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::null, TEUCHOS_TEST_FOR_EXCEPT, and TEUCHOS_UNREACHABLE_RETURN.
|
private |
Definition at line 551 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::null, TEUCHOS_TEST_FOR_EXCEPT, and TEUCHOS_UNREACHABLE_RETURN.
|
private |
Definition at line 559 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::null, TEUCHOS_TEST_FOR_EXCEPT, and TEUCHOS_UNREACHABLE_RETURN.
|
private |
Definition at line 567 of file Thyra_EpetraModelEvaluator.cpp.
References currentInArgsOutArgs_, prototypeOutArgs_, and updateInArgsOutArgs().
|
private |
Definition at line 576 of file Thyra_EpetraModelEvaluator.cpp.
References convertInArgsFromThyraToEpetra(), convertOutArgsFromThyraToEpetra(), epetraInArgsScaling_, epetraModel_, epetraOutArgsScaling_, finishConvertingOutArgsFromEpetraToThyra(), getNominalValues(), Teuchos::includesVerbLevel(), includesVerbLevel(), is_null(), Teuchos::null, postEvalScalingSetup(), preEvalScalingSetup(), Teuchos::rcp(), Teuchos::Time::start(), STATE_FUNC_SCALING_NONE, stateFunctionScaling_, stateFunctionScalingVec_, Teuchos::Time::stop(), TEUCHOS_TEST_FOR_EXCEPTION, Teuchos::Time::totalElapsedTime(), updateNominalValuesAndBounds(), Teuchos::VERB_LOW, and W_factory_.
|
private |
Definition at line 804 of file Thyra_EpetraModelEvaluator.cpp.
References Thyra::create_Vector(), Teuchos::implicit_cast(), p_space_, TEUCHOS_TEST_FOR_EXCEPT, and x_space_.
Referenced by updateNominalValuesAndBounds().
|
private |
Definition at line 847 of file Thyra_EpetraModelEvaluator.cpp.
References Teuchos::RCP< T >::get(), Thyra::get_Epetra_Vector(), p_map_, rcp(), Teuchos::rcp(), TEUCHOS_TEST_FOR_EXCEPT, and x_map_.
Referenced by evalModelImpl().
|
private |
Definition at line 949 of file Thyra_EpetraModelEvaluator.cpp.
References convert(), f(), f_map_, g_map_, Teuchos::RCP< T >::get(), Thyra::get_Epetra_Vector(), Teuchos::implicit_cast(), is_null(), nonnull(), p_map_, Teuchos::rcp(), TEUCHOS_ASSERT, and x_map_.
Referenced by evalModelImpl().
|
private |
Definition at line 1144 of file Thyra_EpetraModelEvaluator.cpp.
References epetraModel_, Teuchos::RCP< T >::get(), is_null(), STATE_FUNC_SCALING_ROW_SUM, stateFunctionScaling_, TEUCHOS_ASSERT, and Teuchos::VERB_LOW.
Referenced by evalModelImpl().
|
private |
Definition at line 1207 of file Thyra_EpetraModelEvaluator.cpp.
References epetraModel_, epetraOutArgsScaling_, Teuchos::RCP< T >::get(), Teuchos::includesVerbLevel(), includesVerbLevel(), rcp(), Teuchos::rcp(), STATE_FUNC_SCALING_ROW_SUM, stateFunctionScaling_, stateFunctionScalingVec_, TEUCHOS_TEST_FOR_EXCEPT, and Teuchos::VERB_LOW.
Referenced by evalModelImpl().
|
private |
Definition at line 1274 of file Thyra_EpetraModelEvaluator.cpp.
References nonnull().
Referenced by evalModelImpl().
|
private |
Definition at line 1303 of file Thyra_EpetraModelEvaluator.cpp.
References convertInArgsFromEpetraToThyra(), createInArgs(), epetraInArgsScaling_, epetraModel_, Teuchos::implicit_cast(), invStateVariableScalingVec_, is_null(), lowerBounds_, nominalValues_, nominalValuesAndBoundsAreUpdated_, Teuchos::rcp(), stateVariableScalingVec_, and upperBounds_.
Referenced by evalModelImpl(), getLowerBounds(), getNominalValues(), getStateVariableInvScalingVec(), and getUpperBounds().
|
private |
Definition at line 1376 of file Thyra_EpetraModelEvaluator.cpp.
References convert(), currentInArgsOutArgs_, description(), epetraModel_, prototypeInArgs_, and prototypeOutArgs_.
Referenced by createInArgs(), and createOutArgsImpl().
|
private |
Definition at line 1504 of file Thyra_EpetraModelEvaluator.cpp.
References epetraModel_, get_f_space(), and get_x_space().
Referenced by create_W_op().
|
related |
|
related |
Referenced by Thyra::convert(), convertOutArgsFromThyraToEpetra(), and updateInArgsOutArgs().
|
related |
|
related |
|
related |
|
related |
|
private |
Definition at line 365 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by create_epetra_W_op(), description(), evalModelImpl(), get_g_names(), get_p_names(), getEpetraModel(), initialize(), postEvalScalingSetup(), preEvalScalingSetup(), uninitialize(), updateInArgsOutArgs(), and updateNominalValuesAndBounds().
|
private |
Definition at line 367 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by getNonconstParameterList(), getParameterList(), setParameterList(), and unsetParameterList().
|
private |
Definition at line 369 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by description(), evalModelImpl(), get_W_factory(), initialize(), and uninitialize().
|
private |
Definition at line 371 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by convertInArgsFromThyraToEpetra(), convertOutArgsFromThyraToEpetra(), and initialize().
|
private |
Definition at line 372 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by convertInArgsFromThyraToEpetra(), convertOutArgsFromThyraToEpetra(), and initialize().
|
private |
Definition at line 373 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by convertOutArgsFromThyraToEpetra(), and initialize().
|
private |
Definition at line 374 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by initialize().
|
private |
Definition at line 375 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by initialize().
|
private |
Definition at line 376 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by convertOutArgsFromThyraToEpetra(), and initialize().
|
private |
Definition at line 378 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by convertInArgsFromEpetraToThyra(), get_x_space(), and initialize().
|
private |
Definition at line 379 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by convertInArgsFromEpetraToThyra(), get_p_space(), initialize(), and Np().
|
private |
Definition at line 380 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by get_f_space(), and initialize().
|
private |
Definition at line 381 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by get_g_space(), initialize(), and Ng().
|
mutableprivate |
Definition at line 383 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by getNominalValues(), setNominalValues(), and updateNominalValuesAndBounds().
|
mutableprivate |
Definition at line 384 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by getLowerBounds(), and updateNominalValuesAndBounds().
|
mutableprivate |
Definition at line 385 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by getUpperBounds(), and updateNominalValuesAndBounds().
|
mutableprivate |
Definition at line 386 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by initialize(), setStateVariableScalingVec(), and updateNominalValuesAndBounds().
|
private |
Definition at line 388 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by getFinalPoint(), and reportFinalPoint().
|
private |
Definition at line 390 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by evalModelImpl(), postEvalScalingSetup(), preEvalScalingSetup(), and setParameterList().
|
mutableprivate |
Definition at line 391 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by evalModelImpl(), getStateFunctionScalingVec(), initialize(), postEvalScalingSetup(), setParameterList(), setStateFunctionScalingVec(), and uninitialize().
|
private |
Definition at line 393 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by getStateVariableScalingVec(), setStateVariableScalingVec(), uninitialize(), and updateNominalValuesAndBounds().
|
mutableprivate |
Definition at line 394 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by getStateVariableInvScalingVec(), setStateVariableScalingVec(), uninitialize(), and updateNominalValuesAndBounds().
|
mutableprivate |
Definition at line 395 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by evalModelImpl(), initialize(), and updateNominalValuesAndBounds().
|
mutableprivate |
Definition at line 396 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by evalModelImpl(), initialize(), and postEvalScalingSetup().
|
mutableprivate |
Definition at line 398 of file Thyra_EpetraModelEvaluator.hpp.
|
mutableprivate |
Definition at line 399 of file Thyra_EpetraModelEvaluator.hpp.
|
mutableprivate |
Definition at line 401 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by createInArgs(), and updateInArgsOutArgs().
|
mutableprivate |
Definition at line 402 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by createOutArgsImpl(), and updateInArgsOutArgs().
|
mutableprivate |
Definition at line 403 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by createInArgs(), createOutArgsImpl(), initialize(), uninitialize(), and updateInArgsOutArgs().
|
private |
Definition at line 405 of file Thyra_EpetraModelEvaluator.hpp.
Referenced by finalPointWasSolved(), initialize(), and reportFinalPoint().