L-BFGS-B  3.0
Large-scale Bound-constrained Optimization
subsm.f File Reference

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...
 

Function/Subroutine Documentation

◆ subsm()

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.

Parameters
nOn entry n is the dimension of the problem.
On exit n is unchanged.
mOn entry m is the maximum number of variable metric corrections used to define the limited memory matrix.
On exit m is unchanged.
nsubOn entry nsub is the number of free variables.
On exit nsub is unchanged.
indOn entry ind specifies the coordinate indices of free variables.
On exit ind is unchanged.
lOn entry l is the lower bound of x.
On exit l is unchanged.
uOn entry u is the upper bound of x.
On exit u is unchanged.
nbdOn entry nbd represents the type of bounds imposed on the variables, and must be specified as follows: nbd(i)=
  • 0 if x(i) is unbounded,
  • 1 if x(i) has only a lower bound,
  • 2 if x(i) has both lower and upper bounds, and
  • 3 if x(i) has only an upper bound.
On exit nbd is unchanged.
xOn entry x specifies the Cauchy point xcp.
On exit x(i) is the minimizer of Q over the subspace of free variables.
dOn entry d is the reduced gradient of Q at xcp.
On exit d is the Newton direction of Q.
xpused to safeguard the projected Newton direction
xxOn entry it holds the current iterate.
On output it is unchanged.
ggOn entry it holds the gradient at the current iterate.
On output it is unchanged.
wsOn entry this stores S, a set of s-vectors, that defines the limited memory BFGS matrix.
On exit this array is unchanged.
wyOn entry this stores Y, a set of y-vectors, that defines the limited memory BFGS matrix.
On exit this array is unchanged.
thetaOn entry theta is the scaling factor specifying B_0 = theta I.
On exit theta is unchanged.
colOn entry col is the actual number of variable metric corrections stored so far.
On exit col is unchanged.
headOn entry head is the location of the first s-vector (or y-vector) in S (or Y).
On exit col is unchanged.
iwordOn entry iword is unspecified.
On exit iword specifies the status of the subspace solution. iword =
  • 0 if the solution is in the box,
  • 1 if some bound is encountered.
wvworking array
wnOn 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.
iprintmust be set by the user; It controls the frequency and type of output generated:
  • iprint<0 no output is generated;
  • iprint=0 print only one line at the last iteration;
  • 0<iprint<99 print also f and |proj g| every iprint iterations;
  • iprint=99 print details of every iteration except n-vectors;
  • iprint=100 print also the changes of active set and final x;
  • iprint>100 print details of every iteration including x and g;
When iprint > 0, the file iterate.dat will be created to summarize the iteration.
infoOn entry info is unspecified.
On exit info =
  • 0 for normal return,
  • nonzero for abnormal return when the matrix K is ill-conditioned.

Definition at line 124 of file subsm.f.

Referenced by mainlb().

Here is the caller graph for this function: