L-BFGS-B
3.0
Large-scale Bound-constrained Optimization
|
Go to the source code of this file.
Functions/Subroutines | |
subroutine | subsm (n, m, nsub, ind, l, u, nbd, x, d, xp, ws, wy, theta, xx, gg, col, head, iword, wv, wn, iprint, info) |
Performs the subspace minimization. More... | |
subroutine subsm | ( | integer | n, |
integer | m, | ||
integer | nsub, | ||
integer, dimension(nsub) | ind, | ||
double precision, dimension(n) | l, | ||
double precision, dimension(n) | u, | ||
integer, dimension(n) | nbd, | ||
double precision, dimension(n) | x, | ||
double precision, dimension(n) | d, | ||
double precision, dimension(n) | xp, | ||
double precision, dimension(n, m) | ws, | ||
double precision, dimension(n, m) | wy, | ||
double precision | theta, | ||
double precision, dimension(n) | xx, | ||
double precision, dimension(n) | gg, | ||
integer | col, | ||
integer | head, | ||
integer | iword, | ||
double precision, dimension(2*m) | wv, | ||
double precision, dimension(2*m, 2*m) | wn, | ||
integer | iprint, | ||
integer | info | ||
) |
Given xcp, l, u, r, an index set that specifies the active set at xcp, and an l-BFGS matrix B (in terms of WY, WS, SY, WT, head, col, and theta), this subroutine computes an approximate solution of the subspace problem
(P) min Q(x) = r'(x-xcp) + 1/2 (x-xcp)' B (x-xcp)
subject to l<=x<=u x_i=xcp_i for all i in A(xcp)
along the subspace unconstrained Newton direction
d = -(Z'BZ)^(-1) r.
The formula for the Newton direction, given the L-BFGS matrix and the Sherman-Morrison formula, is
d = (1/theta)r + (1/theta*2) Z'WK^(-1)W'Z r.
where K = [-D -Y'ZZ'Y/theta L_a'-R_z' ] [L_a -R_z theta*S'AA'S ]
Note that this procedure for computing d differs from that described in [1]. One can show that the matrix K is equal to the matrix M^[-1]N in that paper.
n | On entry n is the dimension of the problem. On exit n is unchanged. |
m | On entry m is the maximum number of variable metric corrections used to define the limited memory matrix. On exit m is unchanged. |
nsub | On entry nsub is the number of free variables. On exit nsub is unchanged. |
ind | On entry ind specifies the coordinate indices of free variables. On exit ind is unchanged. |
l | On entry l is the lower bound of x. On exit l is unchanged. |
u | On entry u is the upper bound of x. On exit u is unchanged. |
nbd | On entry nbd represents the type of bounds imposed on the variables, and must be specified as follows: nbd(i)=
|
x | On entry x specifies the Cauchy point xcp. On exit x(i) is the minimizer of Q over the subspace of free variables. |
d | On entry d is the reduced gradient of Q at xcp. On exit d is the Newton direction of Q. |
xp | used to safeguard the projected Newton direction |
xx | On entry it holds the current iterate. On output it is unchanged. |
gg | On entry it holds the gradient at the current iterate. On output it is unchanged. |
ws | On entry this stores S, a set of s-vectors, that defines the limited memory BFGS matrix. On exit this array is unchanged. |
wy | On entry this stores Y, a set of y-vectors, that defines the limited memory BFGS matrix. On exit this array is unchanged. |
theta | On entry theta is the scaling factor specifying B_0 = theta I. On exit theta is unchanged. |
col | On entry col is the actual number of variable metric corrections stored so far. On exit col is unchanged. |
head | On entry head is the location of the first s-vector (or y-vector) in S (or Y). On exit col is unchanged. |
iword | On entry iword is unspecified. On exit iword specifies the status of the subspace solution. iword =
|
wv | working array |
wn | On entry the upper triangle of wn stores the LEL^T factorization of the indefinite matrix K = [-D -Y'ZZ'Y/theta L_a'-R_z' ] [L_a -R_z theta*S'AA'S ] where E = [-I 0] [ 0 I] On exit wn is unchanged. |
iprint | must be set by the user; It controls the frequency and type of output generated:
|
info | On entry info is unspecified. On exit info =
|
Definition at line 124 of file subsm.f.
Referenced by mainlb().