lbfgsb_module Module

This is a modernized version L-BFGS-B.

It is based on the the modified version of L-BFGS-B described in the following paper:

  • Jorge Nocedal and Jose Luis Morales, Remark on "Algorithm 778: L-BFGS-B: Fortran Subroutines for Large-Scale Bound Constrained Optimization" (2011). ACM Transactions on Mathematical Software, Volume 38, Issue 1, Article No.: 7, pp 1-4

The paper describes an improvement and a correction to Algorithm 778. It is shown that the performance of the algorithm can be improved significantly by making a relatively simple modication to the subspace minimization phase. The correction concerns an error caused by the use of routine dpmeps to estimate machine precision.

License

L-BFGS-B is released under the "New BSD License" (aka "Modified BSD License" or "3-clause license") Please read attached file License.txt

Version

  • Based on: L-BFGS-B (version 3.0). April 25, 2011
  • Refactored and modernized by Jacob Williams, 2023

Original Authors

  • J. Nocedal Department of Electrical Engineering and Computer Science. Northwestern University. Evanston, IL. USA
  • J.L Morales Departamento de Matematicas, Instituto Tecnologico Autonomo de Mexico Mexico D.F. Mexico.

Todo

make a high-level wrapper so the user doesn't have to call the reverse communication routine directly.


Uses

  • module~~lbfgsb_module~~UsesGraph module~lbfgsb_module lbfgsb_module iso_fortran_env iso_fortran_env module~lbfgsb_module->iso_fortran_env module~lbfgsb_blas_module lbfgsb_blas_module module~lbfgsb_module->module~lbfgsb_blas_module module~lbfgsb_kinds_module lbfgsb_kinds_module module~lbfgsb_module->module~lbfgsb_kinds_module module~lbfgsb_linpack_module lbfgsb_linpack_module module~lbfgsb_module->module~lbfgsb_linpack_module module~lbfgsb_blas_module->module~lbfgsb_kinds_module module~lbfgsb_kinds_module->iso_fortran_env module~lbfgsb_linpack_module->module~lbfgsb_blas_module module~lbfgsb_linpack_module->module~lbfgsb_kinds_module

Variables

Type Visibility Attributes Name Initial
integer, public, parameter :: lbfgsp_wp = wp
real(kind=wp), private, parameter :: zero = 0.0_wp
real(kind=wp), private, parameter :: one = 1.0_wp
real(kind=wp), private, parameter :: two = 2.0_wp
real(kind=wp), private, parameter :: three = 3.0_wp

Subroutines

public subroutine setulb(n, m, x, l, u, Nbd, f, g, Factr, Pgtol, Wa, Iwa, Task, Iprint, Csave, Lsave, Isave, Dsave, iteration_file)

The main routine.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the problem (the number of variables).

integer, intent(in) :: m

the maximum number of variable metric corrections used to define the limited memory matrix. Values of m < 3 are not recommended, and large values of m can result in excessive computing time. The range 3 <= m <= 20 is recommended.

real(kind=wp), intent(inout) :: x(n)

On initial entry it must be set by the user to the values of the initial estimate of the solution vector. Upon successful exit, it contains the values of the variables at the best point found (usually an approximate solution).

real(kind=wp), intent(in) :: l(n)

the lower bound on x. If the i-th variable has no lower bound, l(i) need not be defined.

real(kind=wp), intent(in) :: u(n)

the upper bound on x. If the i-th variable has no upper bound, u(i) need not be defined.

integer, intent(in) :: Nbd(n)

nbd represents the type of bounds imposed on the variables, and must be specified as follows:

Read more…
real(kind=wp), intent(inout) :: f

On first entry f is unspecified. On final exit f is the value of the function at x. If the setulb returns with task(1:2)= 'FG', then f must be set by the user to contain the value of the function at the point x.

real(kind=wp), intent(inout) :: g(n)

On first entry g is unspecified. On final exit g is the value of the gradient at x. If the setulb returns with task(1:2)= 'FG', then g must be set by the user to contain the components of the gradient at the point x.

real(kind=wp), intent(in) :: Factr

A tolerance in the termination test for the algorithm. The iteration will stop when:

Read more…
real(kind=wp), intent(in) :: Pgtol

The iteration will stop when:

Read more…
real(kind=wp) :: Wa(2*m*n+5*n+11*m*m+8*m)

working array of length (2mmax + 5)nmax + 11mmax^2 + 8mmax. This array must not be altered by the user.

integer :: Iwa(3*n)

an integer working array of length 3nmax. This array must not be altered by the user.

character(len=60), intent(inout) :: Task

Indicates the current job when entering and quitting this subroutine:

Read more…
integer, intent(in) :: Iprint

Controls the frequency and type of output generated:

Read more…
character(len=60) :: Csave

working string

logical :: Lsave(4)

A logical working array of dimension 4. On exit with task = 'NEW_X', the following information is available:

Read more…
integer :: Isave(44)

An integer working array of dimension 44. On exit with 'task' = NEW_X, the following information is available:

Read more…
real(kind=wp) :: Dsave(29)

A real(wp) working array of dimension 29. On exit with 'task' = NEW_X, the following information is available:

Read more…
character(len=*), intent(in), optional :: iteration_file

The name of the iteration file if used (Iprint>=1). If not specified, then 'iterate.dat' is used.

private subroutine mainlb(n, m, x, l, u, Nbd, f, g, Factr, Pgtol, Ws, Wy, Sy, Ss, Wt, Wn, Snd, z, r, d, t, Xp, Wa, Index, Iwhere, Indx2, Task, Iprint, Csave, Lsave, Isave, Dsave, iteration_file)

This subroutine solves bound constrained optimization problems by using the compact formula of the limited memory BFGS updates.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the number of variables.

integer, intent(in) :: m

the maximum number of variable metric corrections allowed in the limited memory matrix.

real(kind=wp), intent(inout) :: x(n)

On entry x is an approximation to the solution. On exit x is the current approximation.

real(kind=wp), intent(in) :: l(n)

the lower bound of x.

real(kind=wp), intent(in) :: u(n)

the upper bound of x.

integer, intent(in) :: Nbd(n)

the type of bounds imposed on the variables, and must be specified as follows:

Read more…
real(kind=wp), intent(inout) :: f

On first entry f is unspecified. On final exit f is the value of the function at x.

real(kind=wp), intent(inout) :: g(n)

On first entry g is unspecified. On final exit g is the value of the gradient at x.

real(kind=wp), intent(in) :: Factr

On entry factr >= 0 is specified by the user. The iteration will stop when

Read more…
real(kind=wp), intent(in) :: Pgtol

On entry pgtol >= 0 is specified by the user. The iteration will stop when

Read more…
real(kind=wp) :: Ws(n,m)

working array used to store information defining the limited memory BFGS matrix: stores S, the matrix of s-vectors;

real(kind=wp) :: Wy(n,m)

working array used to store information defining the limited memory BFGS matrix: stores Y, the matrix of y-vectors;

real(kind=wp) :: Sy(m,m)

working array used to store information defining the limited memory BFGS matrix: stores S'Y;

real(kind=wp) :: Ss(m,m)

working array used to store information defining the limited memory BFGS matrix: stores S'S;

real(kind=wp) :: Wt(m,m)

working array used to store information defining the limited memory BFGS matrix: stores the Cholesky factorization of (theta*S'S+LD^(-1)L'); see eq. (2.26) in [3].

real(kind=wp) :: Wn(2*m,2*m)

working array used to store the LEL^T factorization of the indefinite matrix

Read more…
real(kind=wp) :: Snd(2*m,2*m)

working array used to store the lower triangular part of

Read more…
real(kind=wp) :: z(n)

working array. used at different times to store the Cauchy point and the Newton point.

real(kind=wp) :: r(n)

working array.

real(kind=wp) :: d(n)

working array.

real(kind=wp) :: t(n)

working array.

real(kind=wp) :: Xp(n)

working array. used to safeguard the projected Newton direction

real(kind=wp) :: Wa(8*m)

working array.

integer :: Index(n)

working array. In subroutine freev, index is used to store the free and fixed variables at the Generalized Cauchy Point (GCP).

integer :: Iwhere(n)

working array used to record the status of the vector x for GCP computation:

Read more…
integer :: Indx2(n)

integer working array. Within subroutine cauchy, indx2 corresponds to the array iorder. In subroutine freev, a list of variables entering and leaving the free set is stored in indx2, and it is passed on to subroutine formk with this information.

character(len=60), intent(inout) :: Task

indicates the current job when entering and leaving this subroutine.

integer, intent(in) :: Iprint

Controls the frequency and type of output generated:

Read more…
character(len=60) :: Csave

working string

logical :: Lsave(4)

working array

integer :: Isave(23)

working array

real(kind=wp) :: Dsave(29)

working array

character(len=*), intent(in), optional :: iteration_file

The name of the iteration file if used (Iprint>=1). If not specified, then 'iterate.dat' is used.

private subroutine active(n, l, u, Nbd, x, Iwhere, Iprint, Prjctd, Cnstnd, Boxed)

This subroutine initializes iwhere and projects the initial x to the feasible set if necessary.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the problem (the number of variables).

real(kind=wp), intent(in) :: l(n)

the lower bound on x.

real(kind=wp), intent(in) :: u(n)

the upper bound on x.

integer, intent(in) :: Nbd(n)
real(kind=wp), intent(inout) :: x(n)
integer, intent(out) :: Iwhere(n)

In cauchy, iwhere is given finer gradations.

Read more…
integer, intent(in) :: Iprint

Controls the frequency and type of output generated

logical, intent(out) :: Prjctd
logical, intent(out) :: Cnstnd
logical, intent(out) :: Boxed

private subroutine bmv(m, Sy, Wt, Col, v, p, Info)

This subroutine computes the product of the 2m x 2m middle matrix in the compact L-BFGS formula of B and a 2m vector v; it returns the product in p.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: m

the maximum number of variable metric corrections used to define the limited memory matrix.

real(kind=wp), intent(in) :: Sy(m,m)

specifies the matrix S'Y.

real(kind=wp), intent(in) :: Wt(m,m)

specifies the upper triangular matrix J' which is the Cholesky factor of (thetaS'S+LD^(-1)L').

integer, intent(in) :: Col

specifies the number of s-vectors (or y-vectors) stored in the compact L-BFGS formula.

real(kind=wp), intent(in) :: v(2*Col)

specifies vector v.

real(kind=wp), intent(out) :: p(2*Col)

the product Mv

integer, intent(out) :: Info

On exit:

Read more…

private subroutine cauchy(n, x, l, u, Nbd, g, Iorder, Iwhere, t, d, Xcp, m, Wy, Ws, Sy, Wt, Theta, Col, Head, p, c, Wbp, v, Nseg, Iprint, Sbgnrm, Info, Epsmch)

For given x, l, u, g (with sbgnrm > 0), and a limited memory BFGS matrix B defined in terms of matrices WY, WS, WT, and scalars head, col, and theta, this subroutine computes the generalized Cauchy point (GCP), defined as the first local minimizer of the quadratic

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the problem.

real(kind=wp), intent(in) :: x(n)

the starting point for the GCP computation.

real(kind=wp), intent(in) :: l(n)

the lower bound of x.

real(kind=wp), intent(in) :: u(n)

the upper bound of x.

integer, intent(in) :: Nbd(n)

On entry nbd represents the type of bounds imposed on the variables, and must be specified as follows:

Read more…
real(kind=wp), intent(in) :: g(n)

the gradient of f(x). g must be a nonzero vector.

integer :: Iorder(n)

working array used to store the breakpoints in the piecewise linear path and free variables encountered. On exit:

Read more…
integer, intent(inout) :: Iwhere(n)

On entry iwhere indicates only the permanently fixed (iwhere=3) or free (iwhere= -1) components of x.

Read more…
real(kind=wp) :: t(n)

used to store the break points.

real(kind=wp) :: d(n)

used to store the Cauchy direction P(x-tg)-x.

real(kind=wp), intent(out) :: Xcp(n)

used to return the GCP on exit.

integer, intent(in) :: m

the maximum number of variable metric corrections used to define the limited memory matrix.

real(kind=wp), intent(in) :: Wy(n,Col)

stores information that defines the limited memory BFGS matrix: wy(n,m) stores Y, a set of y-vectors;

real(kind=wp), intent(in) :: Ws(n,Col)

stores information that defines the limited memory BFGS matrix: ws(n,m) stores S, a set of s-vectors;

real(kind=wp), intent(in) :: Sy(m,m)

stores information that defines the limited memory BFGS matrix: sy(m,m) stores S'Y;

real(kind=wp), intent(in) :: Wt(m,m)

stores information that defines the limited memory BFGS matrix: wt(m,m) stores the Cholesky factorization of (theta*S'S+LD^(-1)L').

real(kind=wp), intent(in) :: Theta

the scaling factor specifying B_0 = theta I.

integer, intent(in) :: Col

the actual number of variable metric corrections stored so far.

integer, intent(in) :: Head

the location of the first s-vector (or y-vector) in S (or Y).

real(kind=wp) :: p(2*m)

working array used to store the vector p = W^(T)d.

real(kind=wp) :: c(2*m)

working array used to store the vector c = W^(T)(xcp-x).

real(kind=wp) :: Wbp(2*m)

working array used to store the row of W corresponding to a breakpoint.

real(kind=wp) :: v(2*m)

working array

integer, intent(out) :: Nseg

records the number of quadratic segments explored in searching for the GCP.

integer, intent(in) :: Iprint

controls the frequency and type of output generated:

Read more…
real(kind=wp), intent(in) :: Sbgnrm

the norm of the projected gradient at x.

integer, intent(inout) :: Info

On entry info is 0. On exit:

Read more…
real(kind=wp), intent(in) :: Epsmch

machine precision

private subroutine cmprlb(n, m, x, g, Ws, Wy, Sy, Wt, z, r, Wa, Index, Theta, Col, Head, Nfree, Cnstnd, Info)

This subroutine computes r=-Z'B(xcp-xk)-Z'g by using wa(2m+1)=W'(xcp-x) from subroutine cauchy.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: n
integer :: m
real(kind=wp) :: x(n)
real(kind=wp) :: g(n)
real(kind=wp) :: Ws(n,m)
real(kind=wp) :: Wy(n,m)
real(kind=wp) :: Sy(m,m)
real(kind=wp) :: Wt(m,m)
real(kind=wp) :: z(n)
real(kind=wp) :: r(n)
real(kind=wp) :: Wa(4*m)
integer :: Index(n)
real(kind=wp) :: Theta
integer :: Col
integer :: Head
integer :: Nfree
logical :: Cnstnd
integer :: Info

private subroutine errclb(n, m, Factr, l, u, Nbd, Task, Info, k)

This subroutine checks the validity of the input data.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n
integer, intent(in) :: m
real(kind=wp), intent(in) :: Factr
real(kind=wp), intent(in) :: l(n)

the lower bound on x.

real(kind=wp), intent(in) :: u(n)

the upper bound on x.

integer, intent(in) :: Nbd(n)
character(len=60) :: Task
integer, intent(inout) :: Info
integer, intent(out) :: k

private subroutine formk(n, Nsub, Ind, Nenter, Ileave, Indx2, Iupdat, Updatd, Wn, Wn1, m, Ws, Wy, Sy, Theta, Col, Head, Info)

This subroutine forms the LEL^T factorization of the indefinite matrix:

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the problem.

integer, intent(in) :: Nsub

the number of subspace variables in free set.

integer, intent(in) :: Ind(n)

specifies the indices of subspace variables.

integer, intent(in) :: Nenter

the number of variables entering the free set.

integer, intent(in) :: Ileave

indx2(ileave),...,indx2(n) are the variables leaving the free set.

integer, intent(in) :: Indx2(n)

On entry indx2(1),...,indx2(nenter) are the variables entering the free set, while indx2(ileave),...,indx2(n) are the variables leaving the free set.

integer, intent(in) :: Iupdat

the total number of BFGS updates made so far.

logical, intent(in) :: Updatd

true if the L-BFGS matrix is updatd.

real(kind=wp) :: Wn(2*m,2*m)

On exit the upper triangle of wn stores the LEL^T factorization of the 2*col x 2*col indefinite matrix

Read more…
real(kind=wp) :: Wn1(2*m,2*m)

On entry wn1 stores the lower triangular part of

Read more…
integer, intent(in) :: m

the maximum number of variable metric corrections used to define the limited memory matrix.

real(kind=wp), intent(in) :: Ws(n,m)

information defining the limited memory BFGS matrix: ws(n,m) stores S, a set of s-vectors

real(kind=wp), intent(in) :: Wy(n,m)

information defining the limited memory BFGS matrix: wy(n,m) stores Y, a set of y-vectors

real(kind=wp), intent(in) :: Sy(m,m)

information defining the limited memory BFGS matrix: sy(m,m) stores S'Y

real(kind=wp), intent(in) :: Theta

information defining the limited memory BFGS matrix: the scaling factor specifying B_0 = theta I

integer, intent(in) :: Col

information defining the limited memory BFGS matrix: the number of variable metric corrections stored

integer, intent(in) :: Head

information defining the limited memory BFGS matrix: the location of the 1st s- (or y-) vector in S (or Y)

integer, intent(out) :: Info

On exit:

Read more…

private subroutine formt(m, Wt, Sy, Ss, Col, Theta, Info)

This subroutine forms the upper half of the pos. def. and symm. T = thetaSS + LD^(-1)L', stores T in the upper triangle of the array wt, and performs the Cholesky factorization of T to produce JJ', with J' stored in the upper triangle of wt.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: m
real(kind=wp) :: Wt(m,m)
real(kind=wp) :: Sy(m,m)
real(kind=wp) :: Ss(m,m)
integer :: Col
real(kind=wp), intent(in) :: Theta
integer :: Info

private subroutine freev(n, Nfree, Index, Nenter, Ileave, Indx2, Iwhere, Wrk, Updatd, Cnstnd, Iprint, Iter)

This subroutine counts the entering and leaving variables when iter > 0, and finds the index set of free and active variables at the GCP.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n
integer, intent(inout) :: Nfree
integer, intent(inout) :: Index(n)

for i=nfree+1,...,n, index(i) are the indices of bound variables

Read more…
integer, intent(out) :: Nenter
integer, intent(out) :: Ileave
integer, intent(inout) :: Indx2(n)

On exit with iter>0, indx2 indicates which variables have changed status since the previous iteration.

Read more…
integer, intent(in) :: Iwhere(n)
logical, intent(out) :: Wrk
logical, intent(in) :: Updatd
logical, intent(in) :: Cnstnd

indicating whether bounds are present

integer, intent(in) :: Iprint

Controls the frequency and type of output generated

integer, intent(in) :: Iter

private subroutine hpsolb(n, t, Iorder, Iheap)

This subroutine sorts out the least element of t, and puts the remaining elements of t in a heap.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the arrays t and iorder.

real(kind=wp), intent(inout) :: t(n)

On entry t stores the elements to be sorted, On exit t(n) stores the least elements of t, and t(1) to t(n-1) stores the remaining elements in the form of a heap.

integer, intent(inout) :: Iorder(n)

On entry iorder(i) is the index of t(i). On exit iorder(i) is still the index of t(i), but iorder may be permuted in accordance with t.

integer, intent(in) :: Iheap

iheap == 0 if t(1) to t(n) is not in the form of a heap, iheap /= 0 if otherwise.

private subroutine lnsrlb(n, l, u, Nbd, x, f, Fold, Gd, Gdold, g, d, r, t, z, Stp, Dnorm, Dtd, Xstep, Stpmx, Iter, Ifun, Iback, Nfgv, Info, Task, Boxed, Cnstnd, Csave, Isave, Dsave)

This subroutine calls subroutine dcsrch from the Minpack2 library to perform the line search. Subroutine dscrch is safeguarded so that all trial points lie within the feasible region.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: n
real(kind=wp) :: l(n)
real(kind=wp) :: u(n)
integer :: Nbd(n)
real(kind=wp) :: x(n)
real(kind=wp) :: f
real(kind=wp) :: Fold
real(kind=wp) :: Gd
real(kind=wp) :: Gdold
real(kind=wp) :: g(n)
real(kind=wp) :: d(n)
real(kind=wp) :: r(n)
real(kind=wp) :: t(n)
real(kind=wp) :: z(n)
real(kind=wp) :: Stp
real(kind=wp) :: Dnorm
real(kind=wp) :: Dtd
real(kind=wp) :: Xstep
real(kind=wp) :: Stpmx
integer :: Iter
integer :: Ifun
integer :: Iback
integer :: Nfgv
integer :: Info
character(len=60) :: Task
logical :: Boxed
logical :: Cnstnd
character(len=60) :: Csave
integer :: Isave(2)
real(kind=wp) :: Dsave(13)

private subroutine matupd(n, m, Ws, Wy, Sy, Ss, d, r, Itail, Iupdat, Col, Head, Theta, Rr, Dr, Stp, Dtd)

This subroutine updates matrices WS and WY, and forms the middle matrix in B.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: n
integer :: m
real(kind=wp) :: Ws(n,m)
real(kind=wp) :: Wy(n,m)
real(kind=wp) :: Sy(m,m)
real(kind=wp) :: Ss(m,m)
real(kind=wp) :: d(n)
real(kind=wp) :: r(n)
integer :: Itail
integer :: Iupdat
integer :: Col
integer :: Head
real(kind=wp) :: Theta
real(kind=wp) :: Rr
real(kind=wp) :: Dr
real(kind=wp) :: Stp
real(kind=wp) :: Dtd

private subroutine prn1lb(n, m, l, u, x, Iprint, Itfile, Epsmch)

This subroutine prints the input data, initial point, upper and lower bounds of each variable, machine precision, as well as the headings of the output.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the problem

integer, intent(in) :: m

the maximum number of variable metric corrections allowed in the limited memory matrix.

real(kind=wp), intent(in) :: l(n)

the lower bound on x.

real(kind=wp), intent(in) :: u(n)

the upper bound on x.

real(kind=wp), intent(in) :: x(n)

initial solution point

integer, intent(in) :: Iprint

Controls the frequency and type of output generated

integer, intent(in) :: Itfile

iteration print file unit

real(kind=wp), intent(in) :: Epsmch

machine precision

private subroutine prn2lb(n, x, f, g, Iprint, Itfile, Iter, Nfgv, Nact, Sbgnrm, Nseg, Word, Iword, Iback, Stp, Xstep)

This subroutine prints out new information after a successful line search.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: n
real(kind=wp) :: x(n)
real(kind=wp) :: f
real(kind=wp) :: g(n)
integer :: Iprint
integer :: Itfile
integer :: Iter
integer :: Nfgv
integer :: Nact
real(kind=wp) :: Sbgnrm
integer :: Nseg
character(len=3), intent(out) :: Word

records the status of subspace solutions.

integer :: Iword
integer :: Iback
real(kind=wp) :: Stp
real(kind=wp) :: Xstep

private subroutine prn3lb(n, x, f, Task, Iprint, Info, Itfile, Iter, Nfgv, Nintol, Nskip, Nact, Sbgnrm, Time, Nseg, Word, Iback, Stp, Xstep, k, Cachyt, Sbtime, Lnscht)

This subroutine prints out information when either a built-in convergence test is satisfied or when an error message is generated.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: n
real(kind=wp) :: x(n)
real(kind=wp) :: f
character(len=60) :: Task
integer :: Iprint
integer :: Info
integer :: Itfile
integer :: Iter
integer :: Nfgv
integer :: Nintol
integer :: Nskip
integer :: Nact
real(kind=wp) :: Sbgnrm
real(kind=wp) :: Time
integer :: Nseg
character(len=3) :: Word
integer :: Iback
real(kind=wp) :: Stp
real(kind=wp) :: Xstep
integer :: k
real(kind=wp) :: Cachyt
real(kind=wp) :: Sbtime
real(kind=wp) :: Lnscht

private subroutine projgr(n, l, u, Nbd, x, g, Sbgnrm)

This subroutine computes the infinity norm of the projected gradient.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the problem (the number of variables).

real(kind=wp), intent(in) :: l(n)

the lower bound on x.

real(kind=wp), intent(in) :: u(n)

the upper bound on x.

integer, intent(in) :: Nbd(n)
real(kind=wp), intent(in) :: x(n)
real(kind=wp), intent(in) :: g(n)
real(kind=wp), intent(out) :: Sbgnrm

infinity norm of the projected gradient

private subroutine subsm(n, m, Nsub, Ind, l, u, Nbd, x, d, Xp, Ws, Wy, Theta, Xx, Gg, Col, Head, Iword, Wv, Wn, Iprint, Info)

This routine contains the major changes in the updated version. The changes are described in the accompanying paper

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the dimension of the problem

integer, intent(in) :: m

the maximum number of variable metric corrections used to define the limited memory matrix.

integer, intent(in) :: Nsub

the number of free variables

integer, intent(in) :: Ind(Nsub)

specifies the coordinate indices of free variables.

real(kind=wp), intent(in) :: l(n)

the lower bound of x.

real(kind=wp), intent(in) :: u(n)

the upper bound of x

integer, intent(in) :: Nbd(n)

represents the type of bounds imposed on the variables, and must be specified as follows:

Read more…
real(kind=wp), intent(inout) :: x(n)

On entry x specifies the Cauchy point xcp. On exit x(i) is the minimizer of Q over the subspace of free variables.

real(kind=wp), intent(inout) :: d(n)

On entry d is the reduced gradient of Q at xcp. On exit d is the Newton direction of Q.

real(kind=wp) :: Xp(n)

used to safeguard the projected Newton direction

real(kind=wp), intent(in) :: Ws(n,m)

variable that stores the information defining the limited memory BFGS matrix: ws(n,m) stores S, a set of s-vectors

real(kind=wp), intent(in) :: Wy(n,m)

variable that stores the information defining the limited memory BFGS matrix:

real(kind=wp), intent(in) :: Theta

variable that stores the information defining the limited memory BFGS matrix: theta is the scaling factor specifying B_0 = theta I

real(kind=wp), intent(in) :: Xx(n)

the current iterate

real(kind=wp), intent(in) :: Gg(n)

the gradient at the current iterate

integer, intent(in) :: Col

variable that stores the information defining the limited memory BFGS matrix: col is the number of variable metric corrections stored

integer, intent(in) :: Head

variable that stores the information defining the limited memory BFGS matrix: head is the location of the 1st s- (or y-) vector in S (or Y)

integer, intent(out) :: Iword

On exit iword specifies the status of the subspace solution:

Read more…
real(kind=wp) :: Wv(2*m)

a real(wp) working array of dimension 2m

real(kind=wp), intent(in) :: Wn(2*m,2*m)

The upper triangle of wn stores the LEL^T factorization of the indefinite matrix

Read more…
integer, intent(in) :: Iprint

If iprint >= 99, some minimum debug printing is done.

integer, intent(out) :: Info

On exit:

Read more…

private subroutine dcsrch(f, g, Stp, Ftol, Gtol, Xtol, Stpmin, Stpmax, Task, Isave, Dsave)

This subroutine finds a step that satisfies a sufficient decrease condition and a curvature condition.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(inout) :: f Read more…
real(kind=wp), intent(inout) :: g Read more…
real(kind=wp), intent(inout) :: Stp Read more…
real(kind=wp), intent(in) :: Ftol

ftol specifies a nonnegative tolerance for the sufficient decrease condition.

real(kind=wp), intent(in) :: Gtol

gtol specifies a nonnegative tolerance for the curvature condition.

real(kind=wp), intent(in) :: Xtol

xtol specifies a nonnegative relative tolerance for an acceptable step. The subroutine exits with a warning if the relative difference between sty and stx is less than xtol.

real(kind=wp), intent(in) :: Stpmin

a nonnegative lower bound for the step.

real(kind=wp), intent(in) :: Stpmax

a nonnegative upper bound for the step.

character(len=*), intent(inout) :: Task

task is a character variable of length at least 60:

Read more…
integer :: Isave(2)

integer work array

real(kind=wp) :: Dsave(13)

real work array

private subroutine dcstep(Stx, Fx, Dx, Sty, Fy, Dy, Stp, Fp, Dp, Brackt, Stpmin, Stpmax)

This subroutine computes a safeguarded step for a search procedure and updates an interval that contains a step that satisfies a sufficient decrease and a curvature condition.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(inout) :: Stx

On entry stx is the best step obtained so far and is an endpoint of the interval that contains the minimizer. On exit `stx is the updated best step.

real(kind=wp), intent(inout) :: Fx

On entry fx is the function at stx. On exit fx is the function at stx.

real(kind=wp), intent(inout) :: Dx

On entry dx is the derivative of the function at stx. The derivative must be negative in the direction of the step, that is, dx and stp - stx must have opposite signs. On exit dx is the derivative of the function at stx.

real(kind=wp), intent(inout) :: Sty

On entry sty is the second endpoint of the interval that contains the minimizer. On exit sty is the updated endpoint of the interval that contains the minimizer.

real(kind=wp), intent(inout) :: Fy

On entry fy is the function at sty. On exit fy is the function at sty.

real(kind=wp), intent(inout) :: Dy

On entry dy is the derivative of the function at sty. On exit dy is the derivative of the function at the exit sty.

real(kind=wp), intent(inout) :: Stp

On entry stp is the current step. If brackt is set to .true. then on input stp must be between stx and sty. On exit stp is a new trial step.

real(kind=wp), intent(in) :: Fp

the function at stp.

real(kind=wp), intent(in) :: Dp

the derivative of the function at stp.

logical, intent(inout) :: Brackt

On entry brackt specifies if a minimizer has been bracketed. Initially brackt must be set to .false. On exit brackt specifies if a minimizer has been bracketed. When a minimizer is bracketed brackt is set to .true.

real(kind=wp), intent(in) :: Stpmin

a lower bound for the step.

real(kind=wp), intent(in) :: Stpmax

an upper bound for the step.