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

Go to the source code of this file.

Functions/Subroutines

subroutine formk (n, nsub, ind, nenter, ileave, indx2, iupdat, updatd, wn, wn1, m, ws, wy, sy, theta, col, head, info)
 Forms the LEL^T factorization of the indefinite matrix K. More...
 

Function/Subroutine Documentation

◆ formk()

subroutine formk ( integer  n,
integer  nsub,
integer, dimension(n)  ind,
integer  nenter,
integer  ileave,
integer, dimension(n)  indx2,
integer  iupdat,
logical  updatd,
double precision, dimension(2*m, 2*m)  wn,
double precision, dimension(2*m, 2*m)  wn1,
integer  m,
double precision, dimension(n, m)  ws,
double precision, dimension(n, m)  wy,
double precision, dimension(m, m)  sy,
double precision  theta,
integer  col,
integer  head,
integer  info 
)

This subroutine forms 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]

The matrix K can be shown to be equal to the matrix M^[-1]N occurring in section 5.1 of [1], as well as to the matrix Mbar^[-1] Nbar in section 5.3.

Parameters
nOn entry n is the dimension of the problem.
On exit n is unchanged.
nsubOn entry nsub is the number of subspace variables in free set.
On exit nsub is not changed.
indOn entry ind specifies the indices of subspace variables.
On exit ind is unchanged.
nenterOn entry nenter is the number of variables entering the free set.
On exit nenter is unchanged.
ileaveOn entry indx2(ileave),...,indx2(n) are the variables leaving the free set.
On exit ileave is unchanged.
indx2On entry indx2(1),...,indx2(nenter) are the variables entering the free set, while indx2(ileave),...,indx2(n) are the variables leaving the free set.
On exit indx2 is unchanged.
iupdatOn entry iupdat is the total number of BFGS updates made so far.
On exit iupdat is unchanged.
updatdOn entry 'updatd' is true if the L-BFGS matrix is updated.
On exit 'updatd' is unchanged.
wnOn entry wn is unspecified.
On exit the upper triangle of wn stores the LEL^T factorization of the 2*col x 2*col indefinite matrix [-D -Y'ZZ'Y/theta L_a'-R_z' ] [L_a -R_z theta*S'AA'S ]
wn1On entry wn1 stores the lower triangular part of [Y' ZZ'Y L_a'+R_z'] [L_a+R_z S'AA'S ] in the previous iteration.
On exit wn1 stores the corresponding updated matrices.
The purpose of wn1 is just to store these inner products so they can be easily updated and inserted into wn.
mOn entry m is the maximum number of variable metric corrections used to define the limited memory matrix.
On exit m 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.
syOn entry this stores S'Y, 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.
infoOn entry info is unspecified.
On exit info
  • = 0 for normal return;
  • = -1 when the 1st Cholesky factorization failed;
  • = -2 when the 2st Cholesky factorization failed.

Definition at line 89 of file formk.f.

Referenced by mainlb().

Here is the caller graph for this function: