EpetraExt
Development
|
These scaling functions implement scaling of input variables and output functions and their derivatives. More...
Classes | |
class | EpetraExt::InArgsGetterSetter_x_dot |
Class that gets and sets x_dot in an InArgs object. More... | |
class | EpetraExt::InArgsGetterSetter_x_dotdot |
Class that gets and sets x_dotdot in an InArgs object. More... | |
class | EpetraExt::InArgsGetterSetter_x |
Class that gets and sets x in an InArgs object. More... | |
class | EpetraExt::InArgsGetterSetter_p |
Class that gets and sets p(l) in an InArgs object. More... | |
class | EpetraExt::OutArgsGetterSetter_f |
Class that gets and sets f in an OutArgs object. More... | |
class | EpetraExt::OutArgsGetterSetter_g |
Class that gets and sets g(j) in an OutArgs object. More... | |
Functions | |
void | EpetraExt::gatherModelNominalValues (const ModelEvaluator &model, ModelEvaluator::InArgs *nominalValues) |
Gather the nominal values from a model evaluator. More... | |
void | EpetraExt::gatherModelBounds (const ModelEvaluator &model, ModelEvaluator::InArgs *lowerBounds, ModelEvaluator::InArgs *upperBounds) |
Gather the lower and upper bounds from a model evaluator. More... | |
void | EpetraExt::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. More... | |
void | EpetraExt::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. More... | |
void | EpetraExt::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. More... | |
void | EpetraExt::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. More... | |
Teuchos::RCP< const Epetra_Vector > | EpetraExt::createInverseModelScalingVector (Teuchos::RCP< const Epetra_Vector > const &scalingVector) |
Create an inverse scaling vector. More... | |
void | EpetraExt::scaleModelVarsGivenInverseScaling (const Epetra_Vector &origVars, const Epetra_Vector &invVarScaling, Epetra_Vector *scaledVars) |
Scale a vector given its inverse scaling vector. More... | |
void | EpetraExt::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. More... | |
void | EpetraExt::unscaleModelVarsGivenInverseScaling (const Epetra_Vector &origVars, const Epetra_Vector &invVarScaling, Epetra_Vector *scaledVars) |
Unscale a vector given its inverse scaling vector. More... | |
void | EpetraExt::scaleModelFuncGivenForwardScaling (const Epetra_Vector &fwdFuncScaling, Epetra_Vector *funcs) |
Scale (in place) an output function vector given its forward scaling vector. More... | |
void | EpetraExt::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 given its function and variable scaling. More... | |
void | EpetraExt::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 scaling. More... | |
These scaling functions implement scaling of input variables and output functions and their derivatives.
The scaling vectors are stored in EpetraExt::ModelEvaluator::InArgs
and EpetraExt::ModelEvaluator::OutArgs
objects in order to enhance maintainability and to avoid programming errors. This will result in some wasted space but it should not be excessive if used carefully.
First, consider scaling of the state function. Reguardless of how the state function scaling is computed, it will be represented as a positive vector s_f
that defines a diagonal matrix S_f = diag(s_f)
that transforms the state function:
f(...) = S_f * f_hat(...)
where f_hat(...)
is the original unscaled state function as computed by the underlying EpetraExt::ModelEvaluator
object and f(...)
is the scaled state function.
Next, consider the scaling of the state varaibles. The scaling for the state variables is defined by a positive vector s_x>/tt> defines a diagonal scaling matrix
S_x = diag(s_x)
that transforms the variables as:
x = S_x * x_hat
where
x_hat
is the original unscaled state variable vector as defined by the underlying EpetraExt::ModelEvaluator
object and x
is the scaled state varaible vector. Note that when the scaled variables x
are passed into evalModel
that they must be unscaled as:
x_hat = inv(S_x) * x
where
inv(S_x)
is the inverse of the diagonals of S_x
which is stored as a positive vector inv_s_x
. Since unscaling the variables as shown above is more common than scaling the original variables, the scaling vector will be stored as inv_s_x
and not as s_x
.
Note how these scalings affect the state function:
f( x_dot, x, ... ) = S_f * f_hat( inv(S_x)*x_dot, inv(S_x)*x, ... )
which has the state/state Jacobian:
W = alpha * d(f)/d(x_dot) + beta * d(f)/d(x) = S_f * ( alpha * d(f_hat)/d(x_hat) + beta * d(f_hat)/d(x) ) * inv(S_x)
Currently, these functions do 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_hat W = S_f * W_hat * inv(S_x) DfDp(l) = S_f * DfDp_hat(l), for l=0...Np-1 g(j) = g_hat(j), for j=0...Ng-1 DgDx_dot(j) = DgDx_dot_hat(j) * inv(S_x), for j=0...Ng-1 DgDx(j) = DgDx_hat(j) * inv(S_x), for j=0...Ng-1 DgDp(j,l) = DgDp_hat(j,l), for j=0...Ng-1, l=0...Np-1
ToDo: Describe how scaling of the state function
S_f
affects the Hessian-vector products an how you just need to scale the Lagrange mutipliers as:
u^T * f(...) = u^T * (S_f * f_hat(...)) = u_f^T * f_hat(...)
where
u_f = S_f * u
.
ToDo: Also describe how scaling of the state varaibles
S_x
affects Hessian-vector products and other related quantities.
These scaling tools must be updated whenever the
InArgs
or OutArgs
classes are augmented. However, not every use case with the model evaluator requires scaling so scaling with respect to some inputs and some outputs may never be needed and therefore never need to be seen by these tools. However, there is some danger in ignoring inputs and outputs in these scaling tools since some objects may be silently unscaled and could cause hard to track down bugs.
ToDo: Finish documentation!
void EpetraExt::gatherModelNominalValues | ( | const ModelEvaluator & | model, |
ModelEvaluator::InArgs * | nominalValues | ||
) |
Gather the nominal values from a model evaluator.
ToDo: Finish documentation!
ToDo: Perhaps refactor the EpetraExt::ModelEvaluator interface to return nominal values as a single InArgs object?
Definition at line 311 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::gatherModelBounds | ( | const ModelEvaluator & | model, |
ModelEvaluator::InArgs * | lowerBounds, | ||
ModelEvaluator::InArgs * | upperBounds | ||
) |
Gather the lower and upper bounds from a model evaluator.
ToDo: Finish documentation!
ToDo: Perhaps refactor the EpetraExt::ModelEvaluator interface to return lower and upper bounds as two different InArgs objects?
Definition at line 349 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::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.
origVars | [in] The orginal unscaled variables. The input data pointed to in this object will not be changed by this function call. |
varScalings | [in] The variable scaling vectors. These scaling vectors must be stored as the inverse scaling vector, such as inv_s_x as described in Scaling Tools for EpetraExt::ModelEvaluator.. |
scaledVars | [in/out] On output, contains copies of the scaled varaibles. On first call, *scaledVars may be empty. Any storage that does not exist will be created. On subsequent calls the storage will be reused. Warning! const casting will be used to allow the modification of the vector objects pointed to by *scaledVars so don't bank on those vectors not being modified. Any vectors pointed to by *scaledVars is fair game to be modified in this function. <it>Precondition:</it> scaledVars!=0 . |
out | [out] If out != 0 then output will be sent to *out . |
verbLevel | [in] Determines the verbosity level for output sent to *out . |
The assumpition, of course, is that the InArgs objects origVars
, varScalings
, and *scaledVars
will all have been created by the same EpetraExt::ModelEvaluator::createOutArgs()
function call.
Definition at line 385 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::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.
ToDo: Finish documentation!
Definition at line 474 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::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.
scaledVars | [in] The scaled varaibles. The input data pointed to in this object will not be modified by this function call. |
varScalings | [in] The variable scaling vectors. These scaling vectors must be stored as the inverse scaling vector, such as inv_s_x as described in Scaling Tools for EpetraExt::ModelEvaluator.. |
origVars | [in/out] On output, contains copies of the unscaled varaibles. On first call, *origVars may be empty. Any storage that does not exist will be created. On subsequent calls the storage will be reused. Warning! const casting will be used to allow the modification of the vector objects pointed to by *origVars so don't bank on those vectors not being modified. Any vectors pointed to by *scaledVars is fair game to be modified in this function. <it>Precondition:</it> origVars!=0 . |
out | [out] If out != 0 then output will be sent to *out . |
verbLevel | [in] Determines the verbosity level for output sent to *out . |
Definition at line 528 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::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.
origFuncs | [in/out] On input, contains the unscaled functions and derivative objects. On output, many to most of the objects pointed to by this object will be scaled in place. See details below. |
varScalings | [in] The variable scaling vectors. These scaling vectors must be stored as the inverse scaling vector, such as inv_s_x as described in Scaling Tools for EpetraExt::ModelEvaluator.. |
funcScalings | [in] The function scaling vectors. These scaling vectors must be stored as the forward scaling vector, such as s_f as described in Scaling Tools for EpetraExt::ModelEvaluator.. |
scaledFuncs | [out] On output, points to in-place scaled functions and derivatives. No new storage is created in this object. Any functions or derivative objects in origFuncs that could not be scaled will not be presented in this object. An output object may not be scaled if the object does not allow scaling. For example, if a derivative object is defined only as an Epetra_Operator object and can not be dynamic cased to a Epetra_RowMatrix , then no explicit scaling will be possible. <it>Precondition:</it> scaledFuncs!=0 . |
allFuncsWhereScaled | [out] On output, determines if all of the functions and derivatives in origFuncs where successfully scaled. If *allFuncsWhereScaled==true on output, then all of the functions in origFuncs where scaled and are represented in *scaledFuncs . If *allFuncsWhereScaled==false on output, then at least one of the functions or derivative objects present in origFuncs was not scaled and is not present in *scaledFuncs . It is up to the client to search *scaledFuncs and compare to origFuncs to see what is missing; Sorry :-(. |
out | [out] If out != 0 then output will be sent to *out . |
verbLevel | [in] Determines the verbosity level for output sent to *out . |
In general, any output objects that are Epetra_MultiVector
(or Epetra_Vector
) or dynamic castable to Epetra_RowMatrix
can be explicitly scaled by this function. Objects that are simply Epetra_Operator
objects can not and will not be scaled by this function and the client is on thier own.
ToDo: Consider using some composite Epetra_Operator classes to create implicitly scaled Epetra_Operator objects and put in an option for doing this or not.
Definition at line 633 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
Teuchos::RCP< const Epetra_Vector > EpetraExt::createInverseModelScalingVector | ( | Teuchos::RCP< const Epetra_Vector > const & | scalingVector | ) |
Create an inverse scaling vector.
This function may actually result in scalingVector
being saved for use later embedded within returnValue
.
Definition at line 797 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::scaleModelVarsGivenInverseScaling | ( | const Epetra_Vector & | origVars, |
const Epetra_Vector & | invVarScaling, | ||
Epetra_Vector * | scaledVars | ||
) |
Scale a vector given its inverse scaling vector.
origVars | [in] The vector of original unscaled varaibles. |
invVarScaling | [in] The inverse scaling vector. |
scaledVars | [out] On output, will contain the scaled varaibles: scaledVars[i] = origVars[i] / invScaleVector[i] , for i=0...n-1 . |
Preconditions:
scaledVars != 0
origVars.Map().SameAs(invVarScaling.Map()) == true
origVars.Map().SameAs(scaledVars->Map()) == true
This function is used by the scaleModelVars()
function to scale each of the varaible vectors.
Definition at line 814 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::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.
Definition at line 831 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::unscaleModelVarsGivenInverseScaling | ( | const Epetra_Vector & | origVars, |
const Epetra_Vector & | invVarScaling, | ||
Epetra_Vector * | scaledVars | ||
) |
Unscale a vector given its inverse scaling vector.
scaledVars | [in] The vector of original unscaled varaibles. |
invVarScaling | [in] The inverse scaling vector. |
origVars | [out] On output, will contain the unscaled varaibles: origVars[i] = invScaleVector[i] * scaledVars[i] , for i=0...n-1 . |
Preconditions:
origVars != 0
scaledVars.Map().SameAs(invVarScaling.Map()) == true
scaledVars.Map().SameAs(origVars->Map()) == true
This function is used by the function unscaleModelVars()
function to unscale each of the varaible vectors.
Definition at line 844 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::scaleModelFuncGivenForwardScaling | ( | const Epetra_Vector & | fwdFuncScaling, |
Epetra_Vector * | funcs | ||
) |
Scale (in place) an output function vector given its forward scaling vector.
fwdFuncScaling | [in] The forward scaling vector. |
funcs | [in/out] On input, contains the vector of unscaled functions. On output, contains the scaled functions: scaledFuncs[i] *= fwdFuncScaling[i] . |
Preconditions:
This function is used by the scaleModelFuncs()
function to scale each of the otuput function vectors.
Definition at line 855 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::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 given its function and variable scaling.
invVarScaling | [in] If invVarScaling !=0 , then this represents the inverse varaible scaling (e.g. inv_s_x ). If invVarScaling==0 , the identity scaling is assumed. |
fwdFuncScaling | [in] If fwdFuncScaling !=0 , then this represents the forward function scaling (e.g. s_f ). If fwdFuncScaling==0 , the identity scaling is assumed. |
funcDerivOp | [in/out] If scaling could be performed, then on output, this object will be scaled. Otherwise, it will not be scaled (see didScaling ).. |
didScaling | [out] On output *didScaling==true if the scaling could be performed. |
Preconditions:
This function is used by the scaleModelFuncs()
function to scale each of the otuput function first derivative objects.
Definition at line 867 of file EpetraExt_ModelEvaluatorScalingTools.cpp.
void EpetraExt::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 scaling.
origFuncDeriv | [in/out] The vector of original unscaled function derivative. If this object can be scaled, then on output it will be scaled in place and scaledFuncDeriv will also point to the scaled derivative object. |
invVarScaling | [in] If invVarScaling !=0 , then this represents the inverse varaible scaling (e.g. inv_s_x ). If invVarScaling==0 , the identity scaling is assumed. |
fwdFuncScaling | [in] If fwdFuncScaling !=0 , then this represents the forward function scaling (e.g. s_f ). If fwdFuncScaling==0 , the identity scaling is assumed. |
scaledFuncDeriv | [out] If scaling could be performed, then on output, this object will point to the scaled function derivative. If not, then scaledFuncDeriv.isEmpty() == true on output. |
didScaling | [out] On output *didScaling==true if the scaling could be performed. |
Preconditions:
This function is used by the scaleModelFuncs()
function to scale each of the otuput function first derivative objects.
Definition at line 891 of file EpetraExt_ModelEvaluatorScalingTools.cpp.