slsqp_core Module

Core subroutines for the SLSQP optimization method. These are refactoried versions of the original routines.


Uses

  • module~~slsqp_core~~UsesGraph module~slsqp_core slsqp_core ieee_arithmetic ieee_arithmetic module~slsqp_core->ieee_arithmetic module~bvls_module bvls_module module~slsqp_core->module~bvls_module module~slsqp_kinds slsqp_kinds module~slsqp_core->module~slsqp_kinds module~slsqp_support slsqp_support module~slsqp_core->module~slsqp_support module~bvls_module->module~slsqp_kinds module~bvls_module->module~slsqp_support iso_fortran_env iso_fortran_env module~slsqp_kinds->iso_fortran_env module~slsqp_support->module~slsqp_kinds

Used by

  • module~~slsqp_core~~UsedByGraph module~slsqp_core slsqp_core module~slsqp_module slsqp_module module~slsqp_module->module~slsqp_core

Variables

Type Visibility Attributes Name Initial
real(kind=wp), private, parameter :: eps = epsilon(1.0_wp)

machine precision


Derived Types

type, public ::  linmin_data

data formerly saved in linmin routine.

Components

Type Visibility Attributes Name Initial
real(kind=wp), public :: a = zero
real(kind=wp), public :: b = zero
real(kind=wp), public :: d = zero
real(kind=wp), public :: e = zero
real(kind=wp), public :: p = zero
real(kind=wp), public :: q = zero
real(kind=wp), public :: r = zero
real(kind=wp), public :: u = zero
real(kind=wp), public :: v = zero
real(kind=wp), public :: w = zero
real(kind=wp), public :: x = zero
real(kind=wp), public :: m = zero
real(kind=wp), public :: fu = zero
real(kind=wp), public :: fv = zero
real(kind=wp), public :: fw = zero
real(kind=wp), public :: fx = zero
real(kind=wp), public :: tol1 = zero
real(kind=wp), public :: tol2 = zero

Type-Bound Procedures

procedure, public :: destroy => destroy_linmin_data

type, public ::  slsqpb_data

data formerly saved in slsqpb.

Components

Type Visibility Attributes Name Initial
real(kind=wp), public :: t = zero
real(kind=wp), public :: f0 = zero
real(kind=wp), public :: h1 = zero
real(kind=wp), public :: h2 = zero
real(kind=wp), public :: h3 = zero
real(kind=wp), public :: h4 = zero
real(kind=wp), public :: t0 = zero
real(kind=wp), public :: gs = zero
real(kind=wp), public :: tol = zero
real(kind=wp), public :: alpha = zero
integer, public :: line = 0
integer, public :: iexact = 0
integer, public :: incons = 0
integer, public :: ireset = 0
integer, public :: itermx = 0
integer, public :: n1 = 0
integer, public :: n2 = 0
integer, public :: n3 = 0

Type-Bound Procedures

procedure, public :: destroy => destroy_slsqpb_data

Functions

private pure function check_convergence(n, f, f0, x, x0, s, h3, acc, tolf, toldf, toldx, converged, not_converged, inconsistent_linearization) result(mode)

Check for convergence.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n
real(kind=wp), intent(in) :: f
real(kind=wp), intent(in) :: f0
real(kind=wp), intent(in), dimension(:) :: x
real(kind=wp), intent(in), dimension(:) :: x0
real(kind=wp), intent(in), dimension(:) :: s
real(kind=wp), intent(in) :: h3
real(kind=wp), intent(in) :: acc
real(kind=wp), intent(in) :: tolf
real(kind=wp), intent(in) :: toldf
real(kind=wp), intent(in) :: toldx
integer, intent(in) :: converged

mode value if converged

integer, intent(in) :: not_converged

mode value if not converged

logical, intent(in) :: inconsistent_linearization

if the SQP problem is inconsistent (will return not_converged)

Return Value integer

private function linmin(mode, ax, bx, f, tol, a, b, d, e, p, q, r, u, v, w, x, m, fu, fv, fw, fx, tol1, tol2)

Linesearch without derivatives (used by slsqp if linesearch_mode=2). Returns the abscissa approximating the point where f attains a minimum.

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(inout) :: mode

controls reverse communication must be set to 0 initially, returns with intermediate values 1 and 2 which must not be changed by the user, ends with convergence with value 3.

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

left endpoint of initial interval

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

right endpoint of initial interval

real(kind=wp) :: f

function value at linmin which is to be brought in by reverse communication controlled by mode

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

desired length of interval of uncertainty of final result

real(kind=wp), intent(inout) :: a
real(kind=wp), intent(inout) :: b
real(kind=wp), intent(inout) :: d
real(kind=wp), intent(inout) :: e
real(kind=wp), intent(inout) :: p
real(kind=wp), intent(inout) :: q
real(kind=wp), intent(inout) :: r
real(kind=wp), intent(inout) :: u
real(kind=wp), intent(inout) :: v
real(kind=wp), intent(inout) :: w
real(kind=wp), intent(inout) :: x
real(kind=wp), intent(inout) :: m
real(kind=wp), intent(inout) :: fu
real(kind=wp), intent(inout) :: fv
real(kind=wp), intent(inout) :: fw
real(kind=wp), intent(inout) :: fx
real(kind=wp), intent(inout) :: tol1
real(kind=wp), intent(inout) :: tol2

Return Value real(kind=wp)


Subroutines

public subroutine slsqp(m, meq, la, n, x, xl, xu, f, c, g, a, acc, iter, mode, w, l_w, sdat, ldat, alphamin, alphamax, tolf, toldf, toldx, max_iter_ls, nnls_mode, infinite_bound)

slsqp: sequential least squares programming to solve general nonlinear optimization problems

Read more…

Arguments

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

is the total number of constraints,

integer, intent(in) :: meq

is the number of equality constraints,

integer, intent(in) :: la

see a,

integer, intent(in) :: n

is the number of variables,

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

x() stores the current iterate of the n vector x on entry x() must be initialized. on exit x() stores the solution vector x if mode = 0.

real(kind=wp), intent(in), dimension(n) :: xl

xl() stores an n vector of lower bounds xl to x.

real(kind=wp), intent(in), dimension(n) :: xu

xu() stores an n vector of upper bounds xu to x.

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

is the value of the objective function.

real(kind=wp), intent(in), dimension(la) :: c

c() stores the m vector c of constraints, equality constraints (if any) first. dimension of c must be greater or equal la, which must be greater or equal max(1,m).

real(kind=wp), intent(in), dimension(n+1) :: g

g() stores the n vector g of partials of the objective function; dimension of g must be greater or equal n+1.

real(kind=wp), intent(in), dimension(la,n+1) :: a

the la by n + 1 array a() stores the m by n matrix a of constraint normals. a() has first dimensioning parameter la, which must be greater or equal max(1,m).

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

abs(acc) controls the final accuracy. if acc < zero an exact linesearch is performed, otherwise an armijo-type linesearch is used.

integer, intent(inout) :: iter

prescribes the maximum number of iterations. on exit iter indicates the number of iterations.

integer, intent(inout) :: mode

mode controls calculation:

Read more…
real(kind=wp), intent(inout), dimension(l_w) :: w

w() is a one dimensional working space. the first m+n+n*n1/2 elements of w must not be changed between subsequent calls of slsqp. on return w(1) ... w(m) contain the multipliers associated with the general constraints, while w(m+1) ... w(m+n(n+1)/2) store the cholesky factor l*d*l(t) of the approximate hessian of the lagrangian columnwise dense as lower triangular unit matrix l with d in its 'diagonal' and w(m+n(n+1)/2+n+2 ... w(m+n(n+1)/2+n+2+m+2n) contain the multipliers associated with all constraints of the quadratic program finding the search direction to the solution x*

integer, intent(in) :: l_w

the length of w, which should be at least:

Read more…
type(slsqpb_data), intent(inout) :: sdat

data for slsqpb.

type(linmin_data), intent(inout) :: ldat

data for linmin.

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

min for line search

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

max for line search

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

stopping criterion if then stop.

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

stopping criterion if then stop.

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

stopping criterion if then stop.

integer, intent(in) :: max_iter_ls

maximum number of iterations in the nnls problem

integer, intent(in) :: nnls_mode

which NNLS method to use

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

"infinity" for the upper and lower bounds. if xl<=-infinite_bound or xu>=infinite_bound then these bounds are considered nonexistant. If infinite_bound=0 then huge() is used for this.

private subroutine slsqpb(m, meq, la, n, x, xl, xu, f, c, g, a, acc, iter, mode, r, l, x0, mu, s, u, v, w, t, f0, h1, h2, h3, h4, n1, n2, n3, t0, gs, tol, line, alpha, iexact, incons, ireset, itermx, ldat, alphamin, alphamax, tolf, toldf, toldx, max_iter_ls, nnls_mode, infBnd)

nonlinear programming by solving sequentially quadratic programs

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: m
integer, intent(in) :: meq
integer, intent(in) :: la
integer, intent(in) :: n
real(kind=wp), dimension(n) :: x
real(kind=wp), dimension(n) :: xl
real(kind=wp), dimension(n) :: xu
real(kind=wp) :: f
real(kind=wp), dimension(la) :: c
real(kind=wp), dimension(n+1) :: g
real(kind=wp), dimension(la,n+1) :: a
real(kind=wp) :: acc
integer, intent(inout) :: iter

in: maximum number of iterations. out: actual number of iterations.

integer, intent(inout) :: mode
real(kind=wp), dimension(m+n+n+2) :: r
real(kind=wp), dimension((n+1)*(n+2)/2) :: l
real(kind=wp), dimension(n) :: x0
real(kind=wp), dimension(la) :: mu
real(kind=wp), dimension(n+1) :: s
real(kind=wp), dimension(n+1) :: u
real(kind=wp), dimension(n+1) :: v
real(kind=wp), intent(inout), dimension(*) :: w

dim(w) =

Read more…
real(kind=wp), intent(inout) :: t
real(kind=wp), intent(inout) :: f0
real(kind=wp), intent(inout) :: h1
real(kind=wp), intent(inout) :: h2
real(kind=wp), intent(inout) :: h3
real(kind=wp), intent(inout) :: h4
integer, intent(inout) :: n1
integer, intent(inout) :: n2
integer, intent(inout) :: n3
real(kind=wp), intent(inout) :: t0
real(kind=wp), intent(inout) :: gs
real(kind=wp), intent(inout) :: tol
integer, intent(inout) :: line
real(kind=wp), intent(inout) :: alpha
integer, intent(inout) :: iexact
integer, intent(inout) :: incons
integer, intent(inout) :: ireset
integer, intent(inout) :: itermx
type(linmin_data), intent(inout) :: ldat

data for linmin.

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

min for line search

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

max for line search

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

stopping criterion if then stop.

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

stopping criterion if then stop

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

stopping criterion if then stop

integer, intent(in) :: max_iter_ls

maximum number of iterations in the nnls problem

integer, intent(in) :: nnls_mode

which NNLS method to use

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

"infinity" for the upper and lower bounds.

private subroutine lsq(m, meq, n, nl, la, l, g, a, b, xl, xu, x, y, w, mode, max_iter_ls, nnls_mode, infbnd)

Minimize with respect to , with upper triangular matrix , and vector , where the unit lower tridiangular matrix is stored columnwise dense in the array with vector stored in its 'diagonal' thus substituting the one-elements of

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: m
integer, intent(in) :: meq
integer, intent(in) :: n
integer, intent(in) :: nl
integer, intent(in) :: la
real(kind=wp), dimension(nl) :: l
real(kind=wp), dimension(n) :: g
real(kind=wp), dimension(la,n) :: a
real(kind=wp), dimension(la) :: b
real(kind=wp), dimension(n) :: xl
real(kind=wp), dimension(n) :: xu
real(kind=wp), dimension(n) :: x

stores the n-dimensional solution vector

real(kind=wp), dimension(m+n+n) :: y

stores the vector of lagrange multipliers of dimension m+n+n (constraints+lower+upper bounds)

real(kind=wp), dimension(*) :: w
integer :: mode

is a success-failure flag with the following meanings:

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

maximum number of iterations in the nnls problem

integer, intent(in) :: nnls_mode

which NNLS method to use

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

"infinity" for the upper and lower bounds.

private subroutine lsei(c, d, e, f, g, h, lc, mc, le, me, lg, mg, n, x, xnrm, w, mode, max_iter_ls, nnls_mode)

for mode=1, the subroutine returns the solution x of equality & inequality constrained least squares problem lsei :

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(inout), dimension(lc,n) :: c
real(kind=wp), intent(inout), dimension(lc) :: d
real(kind=wp), intent(inout), dimension(le,n) :: e
real(kind=wp), intent(inout), dimension(le) :: f
real(kind=wp), intent(inout), dimension(lg,n) :: g
real(kind=wp), intent(inout), dimension(lg) :: h
integer, intent(in) :: lc
integer, intent(in) :: mc
integer, intent(in) :: le
integer, intent(in) :: me
integer, intent(in) :: lg
integer, intent(in) :: mg
integer, intent(in) :: n
real(kind=wp), intent(out), dimension(n) :: x

stores the solution vector

real(kind=wp), intent(out) :: xnrm

stores the residuum of the solution in euclidian norm

real(kind=wp), intent(inout), dimension(*) :: w

on return, stores the vector of lagrange multipliers in its first mc+mg elements

integer, intent(out) :: mode

is a success-failure flag with the following meanings:

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

maximum number of iterations in the nnls problem

integer, intent(in) :: nnls_mode

which NNLS method to use

private subroutine lsi(e, f, g, h, le, me, lg, mg, n, x, xnorm, w, mode, max_iter_ls, nnls_mode)

for mode=1, the subroutine returns the solution x of inequality constrained linear least squares problem:

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(inout), dimension(le,n) :: e
real(kind=wp), intent(inout), dimension(le) :: f
real(kind=wp), intent(inout), dimension(lg,n) :: g
real(kind=wp), intent(inout), dimension(lg) :: h
integer, intent(in) :: le
integer, intent(in) :: me
integer, intent(in) :: lg
integer, intent(in) :: mg
integer, intent(in) :: n
real(kind=wp), intent(out), dimension(n) :: x

stores the solution vector

real(kind=wp), intent(out) :: xnorm

stores the residuum of the solution in euclidian norm

real(kind=wp), intent(inout), dimension(*) :: w

stores the vector of lagrange multipliers in its first mg elements

integer, intent(out) :: mode

is a success-failure flag with the following meanings:

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

maximum number of iterations in the nnls problem

integer, intent(in) :: nnls_mode

which NNLS method to use

private subroutine ldp(g, mg, m, n, h, x, xnorm, w, mode, max_iter_ls, nnls_mode)

Least distance programming routine. Minimize subject to .

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(mg,n) :: g

on entry g stores the m by n matrix of linear inequality constraints. g has first dimensioning parameter mg

integer, intent(in) :: mg
integer, intent(in) :: m
integer, intent(in) :: n
real(kind=wp), intent(in), dimension(m) :: h

the right side of the inequality system.

real(kind=wp), intent(out), dimension(n) :: x

solution vector x if mode=1.

real(kind=wp), intent(out) :: xnorm

euclidian norm of the solution vector if computation is successful

real(kind=wp), intent(inout), dimension(*) :: w

w is a one dimensional working space, the length of which should be at least (m+2)*(n+1) + 2*m. on exit w stores the lagrange multipliers associated with the constraints. at the solution of problem ldp.

integer, intent(out) :: mode

success-failure flag with the following meanings:

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

maximum number of iterations in the nnls problem

integer, intent(in) :: nnls_mode

which NNLS method to use

private subroutine nnls(a, mda, m, n, b, x, rnorm, w, zz, mode, max_iter)

Nonnegative least squares algorithm.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(inout), dimension(mda,n) :: a

on entry, contains the m by n matrix, a. on exit, contains the product matrix, q*a, where q is an m by m orthogonal matrix generated implicitly by this subroutine.

integer, intent(in) :: mda

first dimensioning parameter for the array a.

integer, intent(in) :: m
integer, intent(in) :: n
real(kind=wp), intent(inout), dimension(m) :: b

on entry, contains the m-vector b. on exit, contains q*b.

real(kind=wp), intent(out), dimension(n) :: x

the solution vector.

real(kind=wp), intent(out) :: rnorm

euclidean norm of the residual vector.

real(kind=wp), intent(inout), dimension(n) :: w

array of working space. on exit w will contain the dual solution vector. w will satisfy w(i) = 0 for all i in set p and w(i) <= 0 for all i in set z.

real(kind=wp), intent(inout), dimension(m) :: zz

an m-array of working space.

integer, intent(out) :: mode

this is a success-failure flag with the following meanings:

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

maximum number of iterations (if <=0, then 3*n is used)

private subroutine hfti(a, mda, m, n, b, mdb, nb, tau, krank, rnorm, h, g)

Rank-deficient least squares algorithm using householder forward triangulation with column interchanges.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(inout), dimension(mda,n) :: a

the array a initially contains the matrix of the least squares problem . either m >= n or m < n is permitted. there is no restriction on the rank of a. the matrix a will be modified by the subroutine.

integer, intent(in) :: mda

the first dimensioning parameter of matrix a (mda >= m).

integer, intent(in) :: m
integer, intent(in) :: n
real(kind=wp), intent(inout), dimension(mdb,nb) :: b

if nb = 0 the subroutine will make no reference to the array b. if nb > 0 the array b must initially contain the m x nb matrix b of the the least squares problem ax = b and on return the array b will contain the n x nb solution x.

integer, intent(in) :: mdb

first dimensioning parameter of matrix b (mdb>=max(m,n))

integer, intent(in) :: nb
real(kind=wp), intent(in) :: tau

absolute tolerance parameter for pseudorank determination, provided by the user.

integer, intent(out) :: krank

pseudorank of a, set by the subroutine.

real(kind=wp), intent(out), dimension(nb) :: rnorm

on exit, rnorm(j) will contain the euclidian norm of the residual vector for the problem defined by the j-th column vector of the array b.

real(kind=wp), intent(inout), dimension(n) :: h

array of working space

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

array of working space

private subroutine h12(mode, lpivot, l1, m, u, iue, up, c, ice, icv, ncv)

Construction and/or application of a single householder transformation .

Read more…

Arguments

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

1 or 2 -- selects algorithm h1 to construct and apply a householder transformation, or algorithm h2 to apply a previously constructed transformation.

integer, intent(in) :: lpivot

the index of the pivot element

integer, intent(in) :: l1

if l1 <= m the transformation will be constructed to zero elements indexed from l1 through m. if l1 > m the subroutine does an identity transformation.

integer, intent(in) :: m

see li.

real(kind=wp), intent(inout), dimension(iue,*) :: u

on entry with mode = 1, u contains the pivot vector. iue is the storage increment between elements. on exit when mode = 1, u and up contain quantities defining the vector u of the householder transformation. on entry with mode = 2, u and up should contain quantities previously computed with mode = 1. these will not be modified during the entry with mode = 2. dimension[u(iue,m)]

integer, intent(in) :: iue

see u.

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

see u.

real(kind=wp), intent(inout), dimension(*) :: c

on entry with mode = 1 or 2, c contains a matrix which will be regarded as a set of vectors to which the householder transformation is to be applied. on exit c contains the set of transformed vectors.

integer, intent(in) :: ice

storage increment between elements of vectors in c.

integer, intent(in) :: icv

storage increment between vectors in c.

integer, intent(in) :: ncv

number of vectors in c to be transformed. if ncv <= 0 no operations will be done on c.

private subroutine g1(a, b, c, s, sig)

Compute orthogonal rotation matrix.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp) :: a
real(kind=wp) :: b
real(kind=wp), intent(out) :: c
real(kind=wp), intent(out) :: s
real(kind=wp) :: sig

private subroutine ldl(n, a, z, sigma, w)

- rank-one - update

Read more…

Arguments

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

order of the coefficient matrix a

real(kind=wp), intent(inout), dimension(*) :: a

In: positive definite matrix of dimension n; only the lower triangle is used and is stored column by column as one dimensional array of dimension n*(n+1)/2.

Read more…
real(kind=wp), intent(inout), dimension(*) :: z

vector of dimension n of updating elements.

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

scalar factor by which the modifying dyade is multiplied.

real(kind=wp), intent(inout), dimension(*) :: w

working array of dimension n (used only if ).

private subroutine enforce_bounds(x, xl, xu, infbnd)

enforce the bound constraints on x.

Arguments

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

optimization variable vector

real(kind=wp), intent(in), dimension(:) :: xl

lower bounds (must be same dimension as x)

real(kind=wp), intent(in), dimension(:) :: xu

upper bounds (must be same dimension as x)

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

"infinity" for the upper and lower bounds. Note that NaN may also be used to indicate no bound.

private subroutine destroy_slsqpb_data(me)

Destructor for slsqpb_data type.

Arguments

Type IntentOptional Attributes Name
class(slsqpb_data), intent(out) :: me

private subroutine destroy_linmin_data(me)

Destructor for linmin_data type.

Arguments

Type IntentOptional Attributes Name
class(linmin_data), intent(out) :: me