|
MOOCHO
Version of the Day
|
Simple scalable serial NLP subclass. More...
#include <NLPInterfacePack_ExampleNLPBanded.hpp>

Constructors / initializers | |
| ExampleNLPBanded (size_type nD, size_type nI, size_type bw=1, size_type mU=0, size_type mI=0, value_type xo=0.1, value_type xDl=-NLP::infinite_bound(), value_type xDu=+NLP::infinite_bound(), value_type xIl=-NLP::infinite_bound(), value_type xIu=+NLP::infinite_bound(), value_type hl=-NLP::infinite_bound(), value_type hu=+NLP::infinite_bound(), bool nlp_selects_basis=false, value_type diag_scal=10.0, value_type diag_vary=1.0, bool sym_basis=false, value_type f_offset=0.0, value_type co=0.0, bool ignore_constraints=false) | |
| Constructor. More... | |
Overridden public members from NLP | |
| void | initialize (bool test_setup) |
| bool | is_initialized () const |
| value_type | max_var_bounds_viol () const |
Overridden from NLPVarReductPerm | |
| bool | nlp_selects_basis () const |
Overridden protected methods from NLPSerialPreprocess | |
| bool | imp_nlp_has_changed () const |
| size_type | imp_n_orig () const |
| size_type | imp_m_orig () const |
| size_type | imp_mI_orig () const |
| const DVectorSlice | imp_xinit_orig () const |
| bool | imp_has_var_bounds () const |
| const DVectorSlice | imp_xl_orig () const |
| const DVectorSlice | imp_xu_orig () const |
| const DVectorSlice | imp_hl_orig () const |
| const DVectorSlice | imp_hu_orig () const |
| void | imp_calc_f_orig (const DVectorSlice &x_full, bool newx, const ZeroOrderInfoSerial &zero_order_info) const |
| void | imp_calc_c_orig (const DVectorSlice &x_full, bool newx, const ZeroOrderInfoSerial &zero_order_info) const |
| void | imp_calc_h_orig (const DVectorSlice &x_full, bool newx, const ZeroOrderInfoSerial &zero_order_info) const |
| void | imp_calc_Gf_orig (const DVectorSlice &x_full, bool newx, const ObjGradInfoSerial &obj_grad_info) const |
| bool | imp_get_next_basis (IVector *var_perm_full, IVector *equ_perm_full, size_type *rank_full, size_type *rank) |
| void | imp_report_orig_final_solution (const DVectorSlice &x_orig, const DVectorSlice *lambda_orig, const DVectorSlice *lambdaI_orig, const DVectorSlice *nu_orig, bool is_optimal) |
Overridden protected methods from NLPSerialPreprocessExplJac | |
| size_type | imp_Gc_nz_orig () const |
| size_type | imp_Gh_nz_orig () const |
| void | imp_calc_Gc_orig (const DVectorSlice &x_full, bool newx, const FirstOrderExplInfo &first_order_expl_info) const |
| void | imp_calc_Gh_orig (const DVectorSlice &x_full, bool newx, const FirstOrderExplInfo &first_order_expl_info) const |
Simple scalable serial NLP subclass.
This example NLP is a scalable problem where the basis of the jacobian of the equality constraints is a banded (band width = bw) symmetric positive definite matrix. Both the number of dependnet and independent variables can be varied.
To setup this NLP, the client specifies:
nD : the number of dependent variables nI : the number of independent variables (nI <= nD) bw : the band width (defined as the number of diagonals, including the center diagonal i==j that have non-zero element) mU : the number of undecomposed dependent constraints mI : the number of general constraints diag_scal : Constant scaling factor for the diagonal elements diag_vary : The scaling factor for the diagonal elements (to produce illconditioning) between the largets and the smallest. sym_basis : True if the basis (selected by the NLP) is symmetric or not. This NLP is defined as:
min f(x) = (1/2) * sum( x(i)^2, for i = 1..n )
s.t.
c(j) = ( ds(j)*x(j) \
- sum( 3/(k)*x(j-k), k=1...klu(j) ) |
- sum( fu/(k)*x(j+k), k=1...kuu(j) ) | for j = 1...nD
) * (x(nD+q(j)) + 1)^2 |
+ co(j) == 0 /
c(nD+jU) = c(jU) + co(nD+jU) == 0 } for jU = 1...mU
hl(jI) <= x(jI) - x(nD+q(jI)) <= hu(jI) } for jI = 1...mI
xl(i) <= x(i) <= xu(i) } for i = 1...n
where:
n = nD + nI
m = nD + mU
mI = mI
ds(j) = diag_scal * ( (diag_vary - 1)/(nD -1) * (j - 1) + 1 )
/ 3 : if sym_basis = true
fu = |
\ 6 : if sym_basis = false
/ 2 : if floor((j-1)/nI) < nD % nI
q(j) = floor((j-1)/nI) + |
\ 1 : if floor((j-1)/nI) >= nD % nI
/ bw-1 : if j - bw >= 0 \
klu(j) = | |
\ j-1 : if j - bw <= 1 |
| for j=1...nD
/ bw-1 : if j + bw-1 <= nD |
kuu(j) = | |
\ nD-j : if j - bw <= 1 /
In the above formuation, the sums are not computed if the upper bounds on k are zero. The term co(j) is an adjustable term that can be used to manipulate the solution. Note that if co(nD+jI) != 0 above, then the undecomposed dependent equality constraints are inconsistent with the decomposed equalities and therefore the NLP is infeasible. An infeasible NLP can also be created by manipulating xl(i), xu(i), hl(jI), hu(jI) and co(j).
For the above NLP, the Jacobian of the decomposed equalities has Jacobian elements:
/ -3/(j-i) * (x(nD+q(j)) + 1)^2 : i - klu(i) <= j < i
|
| ds(j) * (x(nD+q(j)) + 1)^2 : i == j
d(c(j))/d(x(i)) = |
| -3/(i-j) * (x(nD+q(j)) + 1)^2 : i < j <= i + kuu(i)
|
| 2 * (c(j) - co(j)) / (x(nD+q(j)) + 1) : i == nD + q
|
\ 0 : otherwise
, for j = 1...nD, i = 1...nD+nI
The above definition shows that for the independent variables, the Jacobian elements are written in terms of the constraint c(j). This fact is exploited in the computational routines when this->multi_calc() == true.
For nD == 7, nI == 2, bw = 2 with floor(nD/nI) = 3 and nD % nI = 1, the Jacobian Gc' looks like:
1 | x x x |
2 | x x x x |
3 | x x x x |
4 | x x x x |
5 | x x x x |
6 | x x x x |
7 | x x x |
- - - - - - - - -
1 2 3 4 5 6 7 8 9
ToDo: Finish documentation!
Definition at line 159 of file NLPInterfacePack_ExampleNLPBanded.hpp.
| NLPInterfacePack::ExampleNLPBanded::ExampleNLPBanded | ( | size_type | nD, |
| size_type | nI, | ||
| size_type | bw = 1, |
||
| size_type | mU = 0, |
||
| size_type | mI = 0, |
||
| value_type | xo = 0.1, |
||
| value_type | xDl = -NLP::infinite_bound(), |
||
| value_type | xDu = +NLP::infinite_bound(), |
||
| value_type | xIl = -NLP::infinite_bound(), |
||
| value_type | xIu = +NLP::infinite_bound(), |
||
| value_type | hl = -NLP::infinite_bound(), |
||
| value_type | hu = +NLP::infinite_bound(), |
||
| bool | nlp_selects_basis = false, |
||
| value_type | diag_scal = 10.0, |
||
| value_type | diag_vary = 1.0, |
||
| bool | sym_basis = false, |
||
| value_type | f_offset = 0.0, |
||
| value_type | co = 0.0, |
||
| bool | ignore_constraints = false |
||
| ) |
Constructor.
ToDo: Finish documentation!
Definition at line 53 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
virtual |
Reimplemented from NLPInterfacePack::NLPSerialPreprocessExplJac.
Definition at line 140 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
virtual |
Reimplemented from NLPInterfacePack::NLPSerialPreprocessExplJac.
Definition at line 151 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
virtual |
Implements NLPInterfacePack::NLP.
Definition at line 156 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
virtual |
Reimplemented from NLPInterfacePack::NLPSerialPreprocess.
Definition at line 376 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Reimplemented from NLPInterfacePack::NLPSerialPreprocess.
Definition at line 163 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 168 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 173 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 179 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 185 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 190 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 195 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 200 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 205 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 210 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 215 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 226 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 271 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocess.
Definition at line 295 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Reimplemented from NLPInterfacePack::NLPSerialPreprocess.
Definition at line 306 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Reimplemented from NLPInterfacePack::NLPSerialPreprocess.
Definition at line 365 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocessExplJac.
Definition at line 383 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocessExplJac.
Definition at line 388 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocessExplJac.
Definition at line 393 of file NLPInterfacePack_ExampleNLPBanded.cpp.
|
protectedvirtual |
Implements NLPInterfacePack::NLPSerialPreprocessExplJac.
Definition at line 472 of file NLPInterfacePack_ExampleNLPBanded.cpp.
1.8.6