42#ifndef EPETRA_EXT_MODEL_EVALUATOR_SCALING_TOOLS_H
43#define EPETRA_EXT_MODEL_EVALUATOR_SCALING_TOOLS_H
47#include "Teuchos_Utils.hpp"
188 const ModelEvaluator &model,
189 ModelEvaluator::InArgs *nominalValues
201 const ModelEvaluator &model,
202 ModelEvaluator::InArgs *lowerBounds,
203 ModelEvaluator::InArgs *upperBounds
242 const ModelEvaluator::InArgs &origVars,
243 const ModelEvaluator::InArgs &varScalings,
244 ModelEvaluator::InArgs *scaledVars,
245 Teuchos::FancyOStream *out = 0,
246 Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_LOW
254 const ModelEvaluator::InArgs &origLowerBounds,
255 const ModelEvaluator::InArgs &origUpperBounds,
257 const ModelEvaluator::InArgs &varScalings,
258 ModelEvaluator::InArgs *scaledLowerBounds,
259 ModelEvaluator::InArgs *scaledUpperBounds,
260 Teuchos::FancyOStream *out = 0,
261 Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_LOW
296 const ModelEvaluator::InArgs &scaledVars,
297 const ModelEvaluator::InArgs &varScalings,
298 ModelEvaluator::InArgs *origVars,
299 Teuchos::FancyOStream *out = 0,
300 Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_LOW
362 const ModelEvaluator::OutArgs &origFuncs,
363 const ModelEvaluator::InArgs &varScalings,
364 const ModelEvaluator::OutArgs &funcScalings,
365 ModelEvaluator::OutArgs *scaledFuncs,
366 bool *allFuncsWhereScaled,
367 Teuchos::FancyOStream *out = 0,
368 Teuchos::EVerbosityLevel verbLevel = Teuchos::VERB_LOW
377Teuchos::RCP<const Epetra_Vector>
379 Teuchos::RCP<const Epetra_Vector>
const& scalingVector
548 const ModelEvaluator::Derivative &origFuncDeriv,
551 ModelEvaluator::Derivative *scaledFuncDeriv,
560 std::string
getName()
const {
return "x_dot"; }
562 Teuchos::RCP<const Epetra_Vector>
569 const Teuchos::RCP<const Epetra_Vector> &x_dot,
574 TEUCHOS_TEST_FOR_EXCEPT(!inArgs);
586 std::string
getName()
const {
return "x_dotdot"; }
588 Teuchos::RCP<const Epetra_Vector>
595 const Teuchos::RCP<const Epetra_Vector> &x_dotdot,
600 TEUCHOS_TEST_FOR_EXCEPT(!inArgs);
614 Teuchos::RCP<const Epetra_Vector>
617 return inArgs.
get_x();
621 const Teuchos::RCP<const Epetra_Vector> &x,
626 TEUCHOS_TEST_FOR_EXCEPT(!inArgs);
641 {
return "p["+Teuchos::Utils::toString(l_)+
"]"; }
643 Teuchos::RCP<const Epetra_Vector>
646 return inArgs.
get_p(l_);
650 const Teuchos::RCP<const Epetra_Vector> &p_l,
655 TEUCHOS_TEST_FOR_EXCEPT(!inArgs);
657 inArgs->
set_p(l_,p_l);
673 Teuchos::RCP<Epetra_Vector>
676 return outArgs.
get_f();
680 const Teuchos::RCP<Epetra_Vector> &f,
685 TEUCHOS_TEST_FOR_EXCEPT(!outArgs);
699 Teuchos::RCP<Epetra_Vector>
702 return outArgs.
get_g(j_);
706 const Teuchos::RCP<Epetra_Vector> &g_j,
711 TEUCHOS_TEST_FOR_EXCEPT(!outArgs);
713 outArgs->
set_g(j_,g_j);
Class that gets and sets p(l) in an InArgs object.
void setVector(const Teuchos::RCP< const Epetra_Vector > &p_l, ModelEvaluator::InArgs *inArgs) const
std::string getName() const
Teuchos::RCP< const Epetra_Vector > getVector(const ModelEvaluator::InArgs &inArgs) const
InArgsGetterSetter_p(int l)
Class that gets and sets x_dot in an InArgs object.
std::string getName() const
Teuchos::RCP< const Epetra_Vector > getVector(const ModelEvaluator::InArgs &inArgs) const
void setVector(const Teuchos::RCP< const Epetra_Vector > &x_dot, ModelEvaluator::InArgs *inArgs) const
Class that gets and sets x_dotdot in an InArgs object.
std::string getName() const
void setVector(const Teuchos::RCP< const Epetra_Vector > &x_dotdot, ModelEvaluator::InArgs *inArgs) const
Teuchos::RCP< const Epetra_Vector > getVector(const ModelEvaluator::InArgs &inArgs) const
Class that gets and sets x in an InArgs object.
std::string getName() const
Teuchos::RCP< const Epetra_Vector > getVector(const ModelEvaluator::InArgs &inArgs) const
void setVector(const Teuchos::RCP< const Epetra_Vector > &x, ModelEvaluator::InArgs *inArgs) const
void set_x_dot(const Teuchos::RCP< const Epetra_Vector > &x_dot)
Teuchos::RCP< const Epetra_Vector > get_x_dotdot() const
Teuchos::RCP< const Epetra_Vector > get_p(int l) const
void set_x_dotdot(const Teuchos::RCP< const Epetra_Vector > &x_dotdot)
void set_x(const Teuchos::RCP< const Epetra_Vector > &x)
Teuchos::RCP< const Epetra_Vector > get_x() const
Set solution vector Taylor polynomial.
Teuchos::RCP< const Epetra_Vector > get_x_dot() const
void set_p(int l, const Teuchos::RCP< const Epetra_Vector > &p_l)
void set_g(int j, const Evaluation< Epetra_Vector > &g_j)
Set g(j) where 0 <= j && j < this->Ng().
Evaluation< Epetra_Vector > get_g(int j) const
Get g(j) where 0 <= j && j < this->Ng().
Evaluation< Epetra_Vector > get_f() const
void set_f(const Evaluation< Epetra_Vector > &f)
Class that gets and sets f in an OutArgs object.
Teuchos::RCP< Epetra_Vector > getVector(const ModelEvaluator::OutArgs &outArgs) const
void setVector(const Teuchos::RCP< Epetra_Vector > &f, ModelEvaluator::OutArgs *outArgs) const
Class that gets and sets g(j) in an OutArgs object.
void setVector(const Teuchos::RCP< Epetra_Vector > &g_j, ModelEvaluator::OutArgs *outArgs) const
OutArgsGetterSetter_g(int j)
Teuchos::RCP< Epetra_Vector > getVector(const ModelEvaluator::OutArgs &outArgs) const
EpetraExt::BlockCrsMatrix: A class for constructing a distributed block matrix.
void scaleModelFuncFirstDeriv(const ModelEvaluator::Derivative &origFuncDeriv, const Epetra_Vector *invVarScaling, const Epetra_Vector *fwdFuncScaling, ModelEvaluator::Derivative *scaledFuncDeriv, bool *didScaling)
Scale (in place) an output first-order function derivative object given its function and variable sca...
void scaleModelFuncs(const ModelEvaluator::OutArgs &origFuncs, const ModelEvaluator::InArgs &varScalings, const ModelEvaluator::OutArgs &funcScalings, ModelEvaluator::OutArgs *scaledFuncs, bool *allFuncsWhereScaled, Teuchos::FancyOStream *out=0, Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_LOW)
Scale the output functions and their derivative objects.
void scaleModelVars(const ModelEvaluator::InArgs &origVars, const ModelEvaluator::InArgs &varScalings, ModelEvaluator::InArgs *scaledVars, Teuchos::FancyOStream *out=0, Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_LOW)
Scale the original unscaled variables into the scaled variables.
void scaleModelFuncFirstDerivOp(const Epetra_Vector *invVarScaling, const Epetra_Vector *fwdFuncScaling, Epetra_Operator *funcDerivOp, bool *didScaling)
Scale (in place) an output first-order function derivative object represented as an Epetra_Operator g...
void scaleModelVarsGivenInverseScaling(const Epetra_Vector &origVars, const Epetra_Vector &invVarScaling, Epetra_Vector *scaledVars)
Scale a vector given its inverse scaling vector.
void unscaleModelVars(const ModelEvaluator::InArgs &scaledVars, const ModelEvaluator::InArgs &varScalings, ModelEvaluator::InArgs *origVars, Teuchos::FancyOStream *out=0, Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_LOW)
Unscale the scaled variables.
void gatherModelNominalValues(const ModelEvaluator &model, ModelEvaluator::InArgs *nominalValues)
Gather the nominal values from a model evaluator.
void scaleModelBounds(const ModelEvaluator::InArgs &origLowerBounds, const ModelEvaluator::InArgs &origUpperBounds, const double infBnd, const ModelEvaluator::InArgs &varScalings, ModelEvaluator::InArgs *scaledLowerBounds, ModelEvaluator::InArgs *scaledUpperBounds, Teuchos::FancyOStream *out=0, Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_LOW)
Scale the lower and upper model variable bounds.
Teuchos::RCP< const Epetra_Vector > createInverseModelScalingVector(Teuchos::RCP< const Epetra_Vector > const &scalingVector)
Create an inverse scaling vector.
void scaleModelFuncGivenForwardScaling(const Epetra_Vector &fwdFuncScaling, Epetra_Vector *funcs)
Scale (in place) an output function vector given its forward scaling vector.
void gatherModelBounds(const ModelEvaluator &model, ModelEvaluator::InArgs *lowerBounds, ModelEvaluator::InArgs *upperBounds)
Gather the lower and upper bounds from a model evaluator.
void scaleModelVarBoundsGivenInverseScaling(const Epetra_Vector &origLowerBounds, const Epetra_Vector &origUpperBounds, const double infBnd, const Epetra_Vector &invVarScaling, Epetra_Vector *scaledLowerBounds, Epetra_Vector *scaledUpperBounds)
Scale the model variable bounds.
void unscaleModelVarsGivenInverseScaling(const Epetra_Vector &origVars, const Epetra_Vector &invVarScaling, Epetra_Vector *scaledVars)
Unscale a vector given its inverse scaling vector.