psqp_module Module

PSQP: SQP variable metric method for general nonlinear programming problems.

History


Uses

  • module~~psqp_module~~UsesGraph module~psqp_module psqp_module module~psqp_kind_module psqp_kind_module module~psqp_module->module~psqp_kind_module module~psqp_matrix_module psqp_matrix_module module~psqp_module->module~psqp_matrix_module iso_fortran_env iso_fortran_env module~psqp_kind_module->iso_fortran_env module~psqp_matrix_module->module~psqp_kind_module

Variables

Type Visibility Attributes Name Initial
integer, public, parameter :: psqp_wp = wp

export the working precision


Abstract Interfaces

abstract interface

  • private subroutine obj_func(me, nf, x, ff)

    objective function interface

    Arguments

    Type IntentOptional Attributes Name
    class(psqp_class), intent(inout) :: me
    integer :: nf

    the number of variables

    real(kind=wp) :: x(nf)

    a vector of variables

    real(kind=wp) :: ff

    the value of the objective function

abstract interface

  • private subroutine dobj_func(me, nf, x, gf)

    gradient of the objective function interface

    Arguments

    Type IntentOptional Attributes Name
    class(psqp_class), intent(inout) :: me
    integer :: nf

    the number of variables

    real(kind=wp) :: x(nf)

    a vector of variables

    real(kind=wp) :: gf(nf)

    the gradient of the objective function

abstract interface

  • private subroutine con_func(me, nf, kc, x, fc)

    constraint function interface

    Arguments

    Type IntentOptional Attributes Name
    class(psqp_class), intent(inout) :: me
    integer :: nf

    the number of variables

    integer :: kc

    the index of the constraint function

    real(kind=wp) :: x(nf)

    a vector of variables

    real(kind=wp) :: fc

    the value of the constraint function

abstract interface

  • private subroutine dcon_func(me, nf, kc, x, gc)

    gradient of the constraint function interface

    Arguments

    Type IntentOptional Attributes Name
    class(psqp_class), intent(inout) :: me
    integer :: nf

    the number of variables

    integer :: kc

    the index of the constraint function

    real(kind=wp) :: x(nf)

    a vector of variables and

    real(kind=wp) :: gc(nf)

    the gradient of the constraint function

abstract interface

  • private subroutine report_f(me, iter, x, j, f)

    Report function to call once per iteration to report the solution.

    Arguments

    Type IntentOptional Attributes Name
    class(psqp_class), intent(inout) :: me
    integer, intent(in) :: iter

    Iteration number

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

    optimization variables

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

    Objective function value

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

    Constraint functions


Derived Types

type, public ::  psqp_class

The main PSQP class to use.

Components

Type Visibility Attributes Name Initial
integer, public :: nres = 0

number of restarts.

integer, public :: ndec = 0

number of matrix decomposition.

integer, public :: nrem = 0

number of constraint deletions.

integer, public :: nadd = 0

number of constraint additions.

integer, public :: nit = 0

number of iterations.

integer, public :: nfv = 0

number of function evaluations.

integer, public :: nfg = 0

number of gradient evaluations.

integer, public :: nfh = 0

number of hessian evaluations.

integer, private :: mtyp = 0
integer, private :: mode = 0
integer, private :: mes1 = 0
integer, private :: mes2 = 0
real(kind=wp), private :: rl = 0.0_wp
real(kind=wp), private :: fl = 0.0_wp
real(kind=wp), private :: ru = 0.0_wp
real(kind=wp), private :: fu = 0.0_wp
real(kind=wp), private :: ri = 0.0_wp
real(kind=wp), private :: fi = 0.0_wp
procedure(obj_func), private, pointer :: obj => null()

objective function

procedure(dobj_func), private, pointer :: dobj => null()

gradient of the objective function

procedure(con_func), private, pointer :: con => null()

constraint function

procedure(dcon_func), private, pointer :: dcon => null()

gradient of the constraint function

procedure(report_f), private, pointer :: report => null()

iteration report function

Type-Bound Procedures

procedure, public :: psqpn
procedure, public :: psqp
procedure, private :: compute_obj_and_dobj
procedure, private :: dual_range_space_quad_prog
procedure, private :: ops_after_constr_deletion
procedure, private :: compute_con_and_dcon
procedure, private :: extended_line_search

Subroutines

private subroutine psqpn(me, nf, nb, nc, x, bound_type, xl, xu, cf, constraint_type, cl, cu, ipar, rpar, f, gmax, cmax, iprnt, iterm, obj, dobj, con, dcon, report)

Date
97/01/22

easy to use subroutine for general nonlinear programming problems.

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
integer, intent(in) :: nf

number of variables

integer, intent(in) :: nb

choice of simple bounds.

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

number of general nonlinear constraints.

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

x(nf) vector of variables.

integer, intent(in), dimension(nf) :: bound_type

ix(nf) vector containing types of bounds.

Read more…
real(kind=wp), intent(in), dimension(nf) :: xl

xl(nf) vector containing lower bounds for variables.

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

xu(nf) vector containing upper bounds for variables.

real(kind=wp), intent(out), dimension(nc+1) :: cf

cf(nc+1) vector containing values of the constraint functions.

integer, intent(in), dimension(nc) :: constraint_type

ic(nc) vector containing types of constraints:

Read more…
real(kind=wp), intent(in), dimension(nc) :: cl

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp), intent(in), dimension(nc) :: cu

cu(nc) vector containing upper bounds for constraint functions.

integer, intent(in), dimension(6) :: ipar

integer paremeters:

Read more…
real(kind=wp), intent(in), dimension(5) :: rpar

real parameters:

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

value of the objective function.

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

maximum partial derivative of the lagrangian function.

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

maximum constraint violation.

integer, intent(in) :: iprnt

print specification:

Read more…
integer, intent(out) :: iterm

variable that indicates the cause of termination.

Read more…
procedure(obj_func) :: obj

computation of the value of the objective function

procedure(dobj_func) :: dobj

computation of the gradient of the objective function

procedure(con_func) :: con

computation of the value of the constraint function

procedure(dcon_func) :: dcon

computation of the gradient of the constraint function

procedure(report_f), optional :: report

iteration report function. Note: this is independent of iprnt. If this function is associated, each iteration will be reported

private subroutine psqp(me, nf, nb, nc, x, ix, xl, xu, cf, ic, cl, cu, cg, cfo, cfd, gc, ica, cr, cz, cp, gf, g, h, s, xo, go, xmax, tolx, tolc, tolg, rpf, cmax, gmax, f, mit, mfv, met, mec, iprnt, iterm)

Date
97/01/22

recursive quadratic programming method with the bfgs variable metric update for general nonlinear programming problems.

Read more…

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
integer :: nf

number of variables.

integer :: nb

choice of simple bounds.

Read more…
integer :: nc

number of linear constraints.

real(kind=wp) :: x(*)

x(nf) vector of variables.

integer, intent(inout) :: ix(*)

ix(nf) vector containing types of bounds.

Read more…
real(kind=wp) :: xl(*)

xl(nf) vector containing lower bounds for variables.

real(kind=wp) :: xu(*)

xu(nf) vector containing upper bounds for variables.

real(kind=wp) :: cf(*)

cf(nc+1) vector containing values of the constraint functions.

integer, intent(inout) :: ic(*)

ic(nc) vector containing types of constraints.

Read more…
real(kind=wp) :: cl(*)

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp) :: cu(*)

cu(nc) vector containing upper bounds for constraint functions.

real(kind=wp) :: cg(*)

cg(nf*nc) matrix whose columns are normals of the linear constraints.

real(kind=wp) :: cfo(*)

cfo(nc) vector containing saved values of the constraint functions.

real(kind=wp) :: cfd(*)

cfd(nc) vector containing increments of the constraint functions.

real(kind=wp) :: gc(*)

gc(nf) gradient of the selected constraint function.

integer :: ica(*)

ica(nf) vector containing indices of active constraints.

real(kind=wp) :: cr(*)

cr(nf*(nf+1)/2) triangular decomposition of kernel of the orthogonal projection.

real(kind=wp) :: cz(*)

cz(nf) vector of lagrange multipliers.

real(kind=wp) :: cp(*)
real(kind=wp) :: gf(*)

gf(nf) gradient of the model function.

real(kind=wp) :: g(*)

g(nf) gradient of the objective function.

real(kind=wp) :: h(*)

h(nf*(nf+1)/2) triangular decomposition or inversion of the hessian matrix approximation.

real(kind=wp) :: s(*)

s(nf) direction vector.

real(kind=wp) :: xo(*)

xo(nf) vectors of variables difference.

real(kind=wp) :: go(*)

go(nf) gradients difference.

real(kind=wp) :: xmax

maximum stepsize.

real(kind=wp) :: tolx

tolerance for change of variables.

real(kind=wp) :: tolc

tolerance for constraint violations.

real(kind=wp) :: tolg

tolerance for the gradient of the lagrangian function.

real(kind=wp) :: rpf

penalty coefficient.

real(kind=wp) :: cmax

maximum constraint violation.

real(kind=wp) :: gmax

maximum partial derivative of the lagrangian function.

real(kind=wp) :: f

value of the objective function.

integer :: mit

maximum number of iterations.

integer :: mfv

maximum number of function evaluations.

integer :: met

variable metric update used.

Read more…
integer :: mec

correction if the negative curvature occurs.

Read more…
integer :: iprnt

print specification.

Read more…
integer :: iterm

variable that indicates the cause of termination.

Read more…

private subroutine compute_con_and_dcon(me, nf, nc, x, fc, cf, cl, cu, ic, gc, cg, cmax, kd, ld)

Date
97/12/01

computation of the value and the gradient of the constraint function.

Read more…

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
integer :: nf

number of variables.

integer :: nc

number of constraints.

real(kind=wp) :: x(nf)

x(nf) vector of variables.

real(kind=wp) :: fc

value of the selected constraint function.

real(kind=wp) :: cf(*)

cf(nc) vector containing values of constraint functions.

real(kind=wp) :: cl(*)

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp) :: cu(*)

cu(nc) vector containing upper bounds for constraint functions.

integer :: ic(*)

ic(nc) vector containing types of constraints.

real(kind=wp) :: gc(nf)

gc(nf) gradient of the selected constraint function.

real(kind=wp) :: cg(*)

cg(nf*nc) matrix whose columns are gradients of constraint functions.

real(kind=wp) :: cmax

maximum constraint violation.

integer :: kd

degree of required derivatives.

integer :: ld

degree of previously computed derivatives.

private subroutine compute_obj_and_dobj(me, nf, x, gf, g, ff, f, kd, ld, iext)

Date
97/12/01

computation of the value and the gradient of the objective function.

Read more…

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
integer, intent(in) :: nf

number of variables.

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

x(nf) vector of variables.

real(kind=wp), intent(out) :: gf(nf)

gf(nf) gradient of the model function.

real(kind=wp), intent(out) :: g(nf)

g(nf) gradient of the objective function.

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

value of the model function.

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

value of the objective function.

integer, intent(in) :: kd

degree of required derivatives.

integer, intent(inout) :: ld

degree of previously computed derivatives.

integer, intent(in) :: iext

type of extremum.

Read more…

private subroutine dual_range_space_quad_prog(me, nf, nc, x, ix, xl, xu, cf, cfd, ic, ica, cl, cu, cg, cr, cz, g, go, h, s, mfp, kbf, kbc, idecf, eta2, eta9, eps7, eps9, umax, gmax, n, iterq)

Date
97/12/01

dual range space quadratic programming method.

Read more…

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
integer :: nf

number of variables.

integer :: nc

number of linear constraints.

real(kind=wp) :: x(*)

x(nf) vector of variables.

integer :: ix(*)

ix(nf) vector containing types of bounds.

real(kind=wp) :: xl(*)

xl(nf) vector containing lower bounds for variables.

real(kind=wp) :: xu(*)

xu(nf) vector containing upper bounds for variables.

real(kind=wp) :: cf(*)

cf(nf) vector containing values of the constraint functions.

real(kind=wp) :: cfd(*)

cfd(nc) vector containing increments of the constraint functions.

integer :: ic(*)

ic(nc) vector containing types of constraints.

integer :: ica(*)

ica(nf) vector containing indices of active constraints.

real(kind=wp) :: cl(*)

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp) :: cu(*)

cu(nc) vector containing upper bounds for constraint functions.

real(kind=wp) :: cg(*)

cg(nf*nc) matrix whose columns are normals of the linear constraints.

real(kind=wp) :: cr(*)

cr(nf*(nf+1)/2) triangular decomposition of kernel of the orthogonal projection.

real(kind=wp) :: cz(*)

cz(nf) vector of lagrange multipliers.

real(kind=wp) :: g(*)

g(nf) gradient of the lagrangian function.

real(kind=wp) :: go(*)

go(nf) saved gradient of the objective function.

real(kind=wp) :: h(*)

h(nf*(nf+1)/2) triangular decomposition or inversion of the hessian matrix approximation.

real(kind=wp) :: s(*)

s(nf) direction vector.

integer :: mfp

type of feasible point.

Read more…
integer :: kbf

specification of simple bounds.

Read more…
integer :: kbc

specification of linear constraints.

Read more…
integer :: idecf

decomposition indicator.

Read more…
real(kind=wp) :: eta2

tolerance for positive definiteness of the hessian matrix.

real(kind=wp) :: eta9

maximum for real numbers.

real(kind=wp) :: eps7

tolerance for linear independence of constraints.

real(kind=wp) :: eps9

tolerance for activity of constraints.

real(kind=wp) :: umax

maximum absolute value of a negative lagrange multiplier.

real(kind=wp) :: gmax

maximum absolute value of a partial derivative.

integer :: n

dimension of the manifold defined by active constraints.

integer :: iterq

type of feasible point.

Read more…

private subroutine update_tri_decomp_general(nf, n, ica, cg, cr, h, s, g, eps7, gmax, umax, idecf, inew, nadd, ier, job)

Date
97/12/01

triangular decomposition of kernel of the general projection is updated after constraint addition.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: nf

declared number of variables.

integer :: n

actual number of variables.

integer :: ica(*)

ica(nf) vector containing indices of active constraints.

real(kind=wp) :: cg(*)

cg(nf*nc) matrix whose columns are normals of the linear constraints.

real(kind=wp) :: cr(*)

cr(nf*(nf+1)/2) triangular decomposition of kernel of the orthogonal projection.

real(kind=wp) :: h(*)

h(nf*(nf+1)/2) triangular decomposition or inversion of the hessian matrix approximation.

real(kind=wp) :: s(*)

s(nf) auxiliary vector.

real(kind=wp) :: g(*)

g(nf) vector used in the dual range space quadratic programming method.

real(kind=wp) :: eps7

tolerance for linear independence of constraints.

real(kind=wp) :: gmax

maximum absolute value of a partial derivative.

real(kind=wp) :: umax

maximum absolute value of a negative lagrange multiplier.

integer :: idecf

decomposition indicator.

Read more…
integer :: inew

index of the new active constraint.

integer :: nadd

number of constraint additions.

integer :: ier

error indicator.

integer :: job

specification of computation. output vector g is not or is computed in case when n<=0 if job=0 or job=1 respectively.

private subroutine determine_new_active_linear_constr(nf, nc, cf, cfd, ic, cl, cu, cg, s, eps9, par, kbc, inew, knew)

Date
97/12/01

determination of the new active linear constraint.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: nf

number of variables.

integer :: nc

number of constraints.

real(kind=wp) :: cf(*)

cf(nc) vector containing values of the constraint functions.

real(kind=wp) :: cfd(*)

cfd(nc) vector containing increments of the constraint functions.

integer :: ic(*)

ic(nc) vector containing types of constraints.

real(kind=wp) :: cl(*)

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp) :: cu(*)

cu(nc) vector containing upper bounds for constraint functions.

real(kind=wp) :: cg(*)

cg(nf*nc) matrix whose columns are normals of the linear constraints.

real(kind=wp) :: s(*)

s(nf) direction vector.

real(kind=wp) :: eps9

tolerance for active constraints.

real(kind=wp) :: par

auxiliary variable.

integer :: kbc

specification of linear constraints.

Read more…
integer :: inew

index of the new active constraint.

integer :: knew

signum of the new active normal.

private subroutine determine_new_active_simple_bound(nf, ix, xo, xl, xu, s, kbf, inew, knew, eps9, par)

Date
91/12/01

determination of the new active simple bound.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: nf

declared number of variables.

integer :: ix(*)

ix(nf) vector containing types of bounds.

real(kind=wp) :: xo(*)

xo(nf) saved vector of variables.

real(kind=wp) :: xl(*)

xl(nf) vector containing lower bounds for variables.

real(kind=wp) :: xu(*)

xu(nf) vector containing upper bounds for variables.

real(kind=wp) :: s(*)

s(nf) direction vector.

integer :: kbf

specification of simple bounds.

Read more…
integer :: inew

index of the new active constraint.

integer :: knew

signum of the new normal.

real(kind=wp) :: eps9

tolerance for active constraints.

real(kind=wp) :: par

auxiliary variable.

private subroutine test_simple_bound(x, ix, xl, xu, eps9, i)

Date
97/12/01

test on activity of a given simple bound.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x(*)

x(nf) vector of variables.

integer, intent(inout) :: ix(*)

ix(nf) vector containing types of bounds.

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

xl(nf) vector containing lower bounds for variables.

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

xu(nf) vector containing upper bounds for variables.

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

tolerance for active constraints.

integer, intent(in) :: i

index of tested simple bound.

private subroutine transform_incompatible_qp_subproblem(nc, cf, ic, cl, cu, kbc)

Date
98/12/01

transformation of the incompatible quadratic programming subproblem.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: nc

number of current linear constraints.

real(kind=wp) :: cf(*)

cf(nf) vector containing values of the constraint functions.

integer :: ic(nc)

ic(nc) vector containing types of constraints.

real(kind=wp) :: cl(*)

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp) :: cu(*)

cu(nc) vector containing upper bounds for constraint functions.

integer :: kbc

specification of linear constraints.

Read more…

private subroutine ops_after_constr_deletion(me, nf, nc, ix, ia, iaa, ar, ic, s, n, iold, krem, ier)

Date
91/12/01

operations after constraint deletion.

Read more…

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
integer :: nf

declared number of variables.

integer :: nc

number of constraints.

integer :: ix(*)

ix(nf) vector containing types of bounds.

integer :: ia(*)

ia(na) vector containing types of deviations.

integer :: iaa(*)

iaa(nf+1) vector containing indices of active functions.

real(kind=wp) :: ar(*)

ar((nf+1)*(nf+2)/2) triangular decomposition of kernel of the orthogonal projection.

integer :: ic(*)

ic(nc) vector containing types of constraints.

real(kind=wp) :: s(*)

s(nf+1) auxiliary vector.

integer :: n

actual number of variables.

integer :: iold

index of the old active constraint.

integer :: krem

auxiliary variable.

integer :: ier

error indicator.

private subroutine update_tri_decomp_orthogonal(nf, ica, cr, g, n, iold, krem, ier)

Date
91/12/01

triangular decomposition of kernel of the orthogonal projection is updated after constraint deletion.

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: nf

declared number of variables.

integer :: ica(*)

ica(nf) vector containing indices of active constraints.

real(kind=wp) :: cr(*)

cr(nf*(nf+1)/2) triangular decomposition of kernel of the orthogonal projection.

real(kind=wp) :: g(*)

g(nf) auxiliary vector.

integer :: n

actual number of variables.

integer :: iold

index of the old active constraint.

integer :: krem

auxiliary variable.

integer :: ier

error indicator.

private subroutine line_search_interpolation(ro, rl, ru, ri, fo, fl, fu, fi, po, r, mode, mtyp, merr)

Date
91/12/01

extrapolation or interpolation for line search without directional derivatives.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=wp) :: ro

initial value of the stepsize parameter.

real(kind=wp) :: rl

lower value of the stepsize parameter.

real(kind=wp) :: ru

upper value of the stepsize parameter.

real(kind=wp) :: ri

inner value of the stepsize parameter.

real(kind=wp) :: fo

value of the objective function for r=ro.

real(kind=wp) :: fl

value of the objective function for r=rl.

real(kind=wp) :: fu

value of the objective function for r=ru.

real(kind=wp) :: fi

value of the objective function for r=ri.

real(kind=wp) :: po

initial value of the directional derivative.

real(kind=wp) :: r

value of the stepsize parameter obtained.

integer :: mode

mode of line search.

integer :: mtyp

method selection

Read more…
integer :: merr

error indicator. merr=0 for normal return.

private pure subroutine compute_augmented_lagrangian(nf, n, nc, cf, ic, ica, cl, cu, cz, rpf, fc, f)

Date
97/12/01

computation of value of the augmented lagrangian function.

Read more…

Arguments

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

number of variables.

integer, intent(in) :: n

dimension of the constraint null space.

integer, intent(in) :: nc

number of constraints.

real(kind=wp), intent(in) :: cf(*)

cf(nc+1) vector containing values of the constraints.

integer, intent(in) :: ic(*)

ic(nc) vector containing types of constraints.

integer, intent(in) :: ica(*)

ica(nf) vector containing indices of active constraints.

real(kind=wp), intent(in) :: cl(*)

cl(nc) vector containing lower bounds for constraint functions.

real(kind=wp), intent(in) :: cu(*)

cu(nc) vector containing upper bounds for constraint functions.

real(kind=wp), intent(in) :: cz(*)

cz(nc) vector of lagrange multipliers.

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

penalty coefficient.

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

value of the penalty term.

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

value of the penalty function.

private pure subroutine compute_new_penalty_parameters(nf, n, nc, ica, cz, cp)

Date
97/12/01

computation of the new penalty parameters.

Read more…

Arguments

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

declared number of variables.

integer, intent(in) :: n

actual number of variables.

integer, intent(in) :: nc

number of constraints.

integer, intent(in) :: ica(*)

vector containing indices of active constraints.

real(kind=wp), intent(in) :: cz(*)

vector of lagrange multipliers.

real(kind=wp), intent(inout) :: cp(*)

vector containing penalty parameters.

private subroutine extended_line_search(me, r, ro, rp, f, fo, fp, po, pp, fmin, fmax, rmin, rmax, tols, kd, ld, nit, kit, nred, mred, maxst, iest, inits, iters, kters, mes, isys)

Date
97/12/01

extended line search without directional derivatives.

Read more…

Arguments

Type IntentOptional Attributes Name
class(psqp_class), intent(inout) :: me
real(kind=wp) :: r

value of the stepsize parameter.

real(kind=wp) :: ro

initial value of the stepsize parameter.

real(kind=wp) :: rp

previous value of the stepsize parameter.

real(kind=wp) :: f

value of the objective function.

real(kind=wp) :: fo

initial value of the objective function.

real(kind=wp) :: fp

previous value of the objective function.

real(kind=wp) :: po

initial value of the directional derivative.

real(kind=wp) :: pp

previous value of the directional derivative.

real(kind=wp) :: fmin

lower bound for value of the objective function.

real(kind=wp) :: fmax

upper bound for value of the objective function.

real(kind=wp) :: rmin

minimum value of the stepsize parameter

real(kind=wp) :: rmax

maximum value of the stepsize parameter

real(kind=wp) :: tols

termination tolerance for line search (in test on the change of the function value).

integer :: kd

degree of required dervatives.

integer :: ld

degree of previously computed derivatives.

integer :: nit

actual number of iterations.

integer :: kit

number of the iteration after last restart.

integer :: nred

actual number of extrapolations or interpolations.

integer :: mred

maximum number of extrapolations or interpolations.

integer :: maxst

maximum stepsize indicator. maxst=0 or maxst=1 if maximum stepsize was not or was reached.

integer :: iest

lower bound specification. iest=0 or iest=1 if lower bound is not or is given.

integer :: inits

choice of the initial stepsize.

Read more…
integer :: iters

termination indicator.

Read more…
integer :: kters

termination selection.

Read more…
integer :: mes

method selection.

Read more…
integer :: isys

control parameter.

private subroutine bfgs_variable_metric_update(n, h, g, s, xo, go, r, po, nit, kit, iterh, met, met1, mec)

Date
92/12/01

variable metric update of a dense symmetric positive definite matrix using the factorization b=ldtrans(l).

Read more…

Arguments

Type IntentOptional Attributes Name
integer :: n

actual number of variables.

real(kind=wp) :: h(*)

h(m) factorization b=ldtrans(l) of a positive definite approximation of the hessian matrix.

real(kind=wp) :: g(*)

g(nf) gradient of the objective function.

real(kind=wp) :: s(*)

s(nf) auxiliary vector.

real(kind=wp) :: xo(*)

xo(nf) vectors of variables difference.

real(kind=wp) :: go(*)

go(nf) gradients difference.

real(kind=wp) :: r

value of the stepsize parameter.

real(kind=wp) :: po

old value of the directional derivative.

integer :: nit

actual number of iterations.

integer :: kit

number of the iteration after last restart.

integer :: iterh

termination indicator.

Read more…
integer :: met
integer :: met1

selection of self scaling.

Read more…
integer :: mec

correction if the negative curvature occurs.

Read more…

private subroutine dual_range_space_qp(nf, n, x, xo, ica, cg, cz, g, go, r, f, fo, p, po, cmax, cmaxo, dmax, kd, ld, iters)

Date
91/12/01

dual range space quadratic programming method for minimax approximation.

Read more…

Arguments

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

declared number of variables.

integer, intent(inout) :: n

actual number of variables.

real(kind=wp) :: x(*)

x(nf) vector of variables.

real(kind=wp) :: xo(*)

xo(nf) saved vector of variables.

integer, intent(in) :: ica(*)

ica(nf) vector containing indices of active constraints.

real(kind=wp) :: cg(*)

cg(nf*nc) matrix whose columns are normals of the linear constraints.

real(kind=wp) :: cz(*)

cz(nf) vector of lagrange multipliers.

real(kind=wp) :: g(*)

g(nf) gradient of the lagrangian function.

real(kind=wp) :: go(*)

go(nf) saved gradient of the lagrangian function.

real(kind=wp) :: r

value of the stepsize parameter.

real(kind=wp) :: f

new value of the objective function.

real(kind=wp) :: fo

old value of the objective function.

real(kind=wp) :: p

new value of the directional derivative.

real(kind=wp) :: po

old value of the directional derivative.

real(kind=wp) :: cmax

value of the constraint violation.

real(kind=wp) :: cmaxo

saved value of the constraint violation.

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

maximum relative difference of variables.

integer :: kd
integer :: ld
integer :: iters

termination indicator for steplength determination. iters=0 for zero step.