MOOCHO (Single Doxygen Collection)
Version of the Day
|
Strategy interface for computing the product of the derivatives of the functions of an NLP along given directions using finite differences. More...
#include <NLPInterfacePack_CalcFiniteDiffProd.hpp>
Public Types | |
enum | EFDMethodOrder { FD_ORDER_ONE, FD_ORDER_TWO, FD_ORDER_TWO_CENTRAL, FD_ORDER_TWO_AUTO, FD_ORDER_FOUR, FD_ORDER_FOUR_CENTRAL, FD_ORDER_FOUR_AUTO } |
enum | EFDStepSelect { FD_STEP_ABSOLUTE, FD_STEP_RELATIVE } |
Public Member Functions | |
STANDARD_MEMBER_COMPOSITION_MEMBERS (EFDMethodOrder, fd_method_order) | |
STANDARD_MEMBER_COMPOSITION_MEMBERS (EFDStepSelect, fd_step_select) | |
STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, fd_step_size) | |
Pick the size of the finite difference step. More... | |
STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, fd_step_size_min) | |
Pick the minimum step size under which the finite difference product will not be computed. More... | |
STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, fd_step_size_f) | |
Set the step size for f(x) More... | |
STANDARD_MEMBER_COMPOSITION_MEMBERS (value_type, fd_step_size_c) | |
Set the step size for c(x) More... | |
CalcFiniteDiffProd (EFDMethodOrder fd_method_order=FD_ORDER_FOUR_AUTO, EFDStepSelect fd_step_select=FD_STEP_ABSOLUTE, value_type fd_step_size=-1.0, value_type fd_step_size_min=-1.0, value_type fd_step_size_f=-1.0, value_type fd_step_size_c=-1.0) | |
virtual | ~CalcFiniteDiffProd () |
virtual bool | calc_deriv_product (const Vector &xo, const Vector *xl, const Vector *xu, const Vector &v, const value_type *fo, const Vector *co, bool check_nan_inf, NLP *nlp, value_type *Gf_prod, VectorMutable *Gc_prod, std::ostream *out, bool trace=false, bool dump_all=false) const |
Compute the directional derivatives by finite differences. More... | |
Strategy interface for computing the product of the derivatives of the functions of an NLP along given directions using finite differences.
Specifically, this interface can be used to compute products with the finite difference approximations of the gradient FDGf'*v
and/or the Jacobian FDGc'*v
and/or the Jacobian FDGh'*v
. Doing so together may take advantage of any shared calculations in the NLP. These products are computed using finite differences.
One of several different finite differencing schemes can be used.
FD_ORDER_ONE
: O(eps)
one sided finite differences (two (one) evaluations of func). grad(f(x),x)'*v \approx (1/eps) * ( -f(xo) + f(xo + eps * u) )
FD_ORDER_TWO
: O(eps^2)
one sided finite differences (three (two) evaluations of func). grad(f(x),x)'*v \approx (1/(2*eps)) * ( -3*f(xo) + 4*f(xo + eps*v) -f(xo + 2*eps*v))
FD_ORDER_TWO_CENTRAL
: O(eps^2)
two sided central differences (two evaluations of func). grad(f(x),x)'*v \approx (1/(2*eps)) * ( - f(xo - eps*v) + f(xo + eps*v) )
FD_ORDER_TWO_AUTO
: Use FD_ORDER_TWO_CENTRAL
when not limited by bounds, otherwise use FD_ORDER_TWO
. FD_ORDER_FOUR
: O(eps^4)
one sided finite differences (five (four) evaluations of func). grad(f(x),x)'*v \approx (1/(12*eps)) * ( -25*f(xo) + 48*f(xo + eps*v) - 36*f(xo + 2*eps*v) + 16*f(xo + 3*eps*v) - 3*f(xo + 4*eps*v))
FD_ORDER_FOUR_CENTRAL
: O(eps^4)
two sided central differences (four evaluations of func). grad(f(x),x)'*v \approx (1/(12*eps)) * ( f(xo - 2*eps*v) - 8*f(xo - eps*v) + 8*f(xo + eps*v) - f(xo + 2*eps*v))
FD_ORDER_FOUR_AUTO
: Use FD_ORDER_FOUR_CENTRAL
when not limited by bounds, otherwise use FD_ORDER_FOUR
. The client can select the step sizes that are used to compute the finite differences. First, the same step sizes eps
for all f(x), c(x)
and h(x)
can be selected using fd_step_size()
with a positive value. If fd_step_size()
returns a negative number, then a default value will be determined internally that is appropriate for the chosen method (see fd_method_order()
). Whatever step size represented by fd_step_size()
(or a default) will be scaled by ||xo||inf + 1.0
if fd_step_select() == FD_STEP_ABSOLUTE
. Using the same step sizes for f(x), c(x) and h(x) an advantage in that the NLP implementation may be able exploit shared computations. However, for some applications, it may be advantagous to use different step lengths for f(x)
, c(x)
and h(x)
. These individual step lengths can be set using fd_step_size_f()
, fd_step_size_c()
and fd_step_size_h()
respectively. These step lengths will also be scaled the same as for fd_step_size()
if fd_step_select() == FD_STEP_ABSOLUTE
. If any of these individual step size functions returns a negative number, then the value for fd_step_size()
will be used.
The finite difference perturbations may be limited by the relaxed variable bounds:
xl - max_var_bounds_viol <= x <= xu + max_var_bounds_viol
if variable bounds are present. If it is not expected that bounds will limit the step length (or there are no bounds) then the central difference methods FD_ORDER_TWO_CENTRAL
and FD_ORDER_FOUR_CENTRAL
should be preferred since they are more accurate. The method FD_ORDER_FOUR_AUTO
will use fourth order central differences (FD_ORDER_FOUR_CENTRAL
) unless the variable bounds limit the minimum step size in which case fourth order one sided differences (FD_ORDER_FOUR
) will be used. With one sided differences, the implementation may be able to take a larger (in magnitude) negative step than a positive step (and visa-versa) in which case the one sided methods would have an advantage. Note that the one situation where one sided differences is guaranteed to be able to take steps away from the bounds is when xl - max_var_bounds_viol + fd_step_size <= xu
and v = eta(j)
(i.e. eta(j)
is the jth column of identity). The situation v = eta(j)
occurs when the client is computing the full finite difference approximations to Gf
, Gc
and/or Gh
on variable at a time.
Definition at line 120 of file NLPInterfacePack_CalcFiniteDiffProd.hpp.
Definition at line 124 of file NLPInterfacePack_CalcFiniteDiffProd.hpp.
Enumerator | |
---|---|
FD_STEP_ABSOLUTE |
Use absolute step size |
FD_STEP_RELATIVE |
Use relative step size |
Definition at line 134 of file NLPInterfacePack_CalcFiniteDiffProd.hpp.
NLPInterfacePack::CalcFiniteDiffProd::CalcFiniteDiffProd | ( | EFDMethodOrder | fd_method_order = FD_ORDER_FOUR_AUTO , |
EFDStepSelect | fd_step_select = FD_STEP_ABSOLUTE , |
||
value_type | fd_step_size = -1.0 , |
||
value_type | fd_step_size_min = -1.0 , |
||
value_type | fd_step_size_f = -1.0 , |
||
value_type | fd_step_size_c = -1.0 |
||
) |
Definition at line 63 of file NLPInterfacePack_CalcFiniteDiffProd.cpp.
|
inlinevirtual |
Definition at line 173 of file NLPInterfacePack_CalcFiniteDiffProd.hpp.
NLPInterfacePack::CalcFiniteDiffProd::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | EFDMethodOrder | , |
fd_method_order | |||
) |
NLPInterfacePack::CalcFiniteDiffProd::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | EFDStepSelect | , |
fd_step_select | |||
) |
NLPInterfacePack::CalcFiniteDiffProd::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | value_type | , |
fd_step_size | |||
) |
Pick the size of the finite difference step.
If fd_step_size < 0
then the implementation will try to select it based on the order of method fd_method_order()
that is selected.
NLPInterfacePack::CalcFiniteDiffProd::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | value_type | , |
fd_step_size_min | |||
) |
Pick the minimum step size under which the finite difference product will not be computed.
If fd_step_size_min == 0
then the finite difference computation will always be performed. If fd_step_size_min < 0
then the minimum step size will be determined internally.
NLPInterfacePack::CalcFiniteDiffProd::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | value_type | , |
fd_step_size_f | |||
) |
Set the step size for f(x)
NLPInterfacePack::CalcFiniteDiffProd::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | value_type | , |
fd_step_size_c | |||
) |
Set the step size for c(x)
|
virtual |
Compute the directional derivatives by finite differences.
The computation may fail if NaN
or Inf
is encountered durring any of the computations in which case a NaNInfException
exception will be thrown. Otherwise the computation should be completed successfully.
The finite difference peturbations may be limited by the relaxed variable bounds
xl - max_var_bounds_viol <= x <= xu + max_var_bounds_viol
if variable bounds are present. If these bounds do limit the finite difference step size then a warning will be printed to *out (if out!=NULL
) and the derivatives may be very inaccurate. If bounds do limit the steps taken then it is advisable to use the one sided finite difference (FD_ORDER_ONE) and this implementation can move away from the bounds.
xo | [in] Base point for the unknown variables to compute derivatives at. |
xl | [in] If != NULL then this is the lower variable bounds. |
xu | [in] If != NULL then this is the upper variable bounds. If xl != NULL then xu != NULL must also be true and visa-versa or a std::invalid_arguement exception will be thrown. |
v | [in] The vector for which to form the products with. |
fo | [in] If fo != NULL then *fo should be set to the value of f(xo) . Not useful for FD_ORDER_TWO_CENTRAL . |
co | [in] If co != NULL then *co should be set to the value of c(xo) . Not useful for FD_ORDER_TWO_CENTRAL . |
check_nan_inf | [in] If true , the the computed values will be checked for nan and inf. |
nlp | [in] Used to compute f(x), c(x) and h(x). The current set references to nlp->get_f() , nlp->get_c() and nlp->get_h() will be preserved on output. The NLP must be initialized before input. |
Gf_prod | [out] If != NULL , the will contain the finite difference computed product Gf'*v on output. If == NULL , then f(x) is not specifically computed. |
Gc_prod | [out] If != NULL , the will contain the finite difference computed product Gc'*v . If == NULL , then c(x) is not specifically computed. |
out | [in/out] If out != NULL then any waring or error messages are output here. |
true if the finite difference computations were performed.
Returns false
if the finite difference computations were not performed because the variables bounds limited the step length too much.ToDo: Discuss options!
Definition at line 79 of file NLPInterfacePack_CalcFiniteDiffProd.cpp.