naginterfaces.library.opt.bounds_​quasi_​deriv_​easy

naginterfaces.library.opt.bounds_quasi_deriv_easy(ibound, funct2, bl, bu, x, liw=None, lw=None, data=None)[source]

bounds_quasi_deriv_easy is an easy-to-use quasi-Newton algorithm for finding a minimum of a function , subject to fixed upper and lower bounds on the independent variables , when first derivatives of are available.

It is intended for functions which are continuous and which have continuous first and second derivatives (although it will usually work even if the derivatives have occasional discontinuities).

For full information please refer to the NAG Library document for e04ky

https://www.nag.com/numeric/nl/nagdoc_29.3/flhtml/e04/e04kyf.html

Parameters
iboundint

Indicates whether the facility for dealing with bounds of special forms is to be used. It must be set to one of the following values:

If you are supplying all the and individually.

If there are no bounds on any .

If all the bounds are of the form .

If and .

funct2callable (fc, gc) = funct2(xc, data=None)

You must supply to calculate the values of the function and its first derivative at any point .

It should be tested separately before being used in conjunction with bounds_quasi_deriv_easy (see the E04 Introduction).

Parameters
xcfloat, ndarray, shape

The point at which the function and derivatives are required.

dataarbitrary, optional, modifiable in place

User-communication data for callback functions.

Returns
fcfloat

The value of the function at the current point .

gcfloat, array-like, shape

must be set to the value of the first derivative at the point , for .

blfloat, array-like, shape

The lower bounds .

If is set to , you must set to , for . (If a lower bound is not specified for a particular , the corresponding should be set to .)

If is set to , you must set to ; bounds_quasi_deriv_easy will then set the remaining elements of equal to .

bufloat, array-like, shape

The upper bounds .

If is set to , you must set to , for . (If an upper bound is not specified for a particular , the corresponding should be set to .)

If is set to , you must set to ; bounds_quasi_deriv_easy will then set the remaining elements of equal to .

xfloat, array-like, shape

must be set to a guess at the th component of the position of the minimum, for . The function checks the gradient at the starting point, and is more likely to detect any error in your programming if the initial are nonzero and mutually distinct.

liwNone or int, optional

Note: if this argument is None then a default value will be used, determined as follows: .

The dimension of the array .

lwNone or int, optional

Note: if this argument is None then a default value will be used, determined as follows: .

The dimension of the array .

dataarbitrary, optional

User-communication data for callback functions.

Returns
blfloat, ndarray, shape

The lower bounds actually used by bounds_quasi_deriv_easy.

bufloat, ndarray, shape

The upper bounds actually used by bounds_quasi_deriv_easy.

xfloat, ndarray, shape

The lowest point found during the calculations. Thus, if no exception or warning is raised on exit, is the th component of the position of the minimum.

ffloat

The value of corresponding to the final point stored in .

gfloat, ndarray, shape

The value of corresponding to the final point stored in , for ; the value of for variables not on a bound should normally be close to zero.

iwint, ndarray, shape

If the function exits successfully or = 3 or 5, the first elements of contain information about which variables are currently on their bounds and which are free. Specifically, if is:

  • fixed on its upper bound, is ;

  • fixed on its lower bound, is ;

  • effectively a constant (i.e., ), is ;

  • free, gives its position in the sequence of free variables.

In addition, contains the number of free variables (i.e., ).

The rest of the array is used as workspace.

wfloat, ndarray, shape

If the function exits successfully or = 3 or 5, contains the th element of the projected gradient vector , for . In addition, contains an estimate of the condition number of the projected Hessian matrix (i.e., ). The rest of the array is used as workspace.

Raises
NagValueError
(errno )

On entry, and for some .

(errno )

On entry, and .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

There have been function evaluations.

(errno )

An overflow occurred during computation.

(errno )

The modulus of a variable has become very large. There may be a mistake in , your problem has no finite solution, or the problem needs rescaling.

(errno )

It is very likely that you have made an error forming the gradient.

Warns
NagAlgorithmicWarning
(errno )

The conditions for a minimum have not all been satisfied, but a lower point could not be found.

(errno )

It is probable that a local minimum has been found, but it cannot be guaranteed.

(errno )

It is possible that a local minimum has been found, but it cannot be guaranteed.

(errno )

It is unlikely that a local minimum has been found.

(errno )

It is very unlikely that a local minimum has been found.

Notes

No equivalent traditional C interface for this routine exists in the NAG Library.

bounds_quasi_deriv_easy is applicable to problems of the form:

when first derivatives are available.

Special provision is made for problems which actually have no bounds on the , problems which have only non-negativity bounds, and problems in which and . You must supply a function to calculate the values of and its first derivatives at any point .

From a starting point you supplied there is generated, on the basis of estimates of the curvature of , a sequence of feasible points which is intended to converge to a local minimum of the constrained function. An attempt is made to verify that the final point is a minimum.

A typical iteration starts at the current point where (say) variables are free from both their bounds. The projected gradient vector , whose elements are the derivatives of with respect to the free variables, is known. A unit lower triangular matrix and a diagonal matrix (both of dimension ), such that is a positive definite approximation of the matrix of second derivatives with respect to the free variables (i.e., the projected Hessian) are also held. The equations

are solved to give a search direction , which is expanded to an -vector by an insertion of appropriate zero elements. Then is found such that is approximately a minimum (subject to the fixed bounds) with respect to ; is replaced by , and the matrices and are updated so as to be consistent with the change produced in the gradient by the step . If any variable actually reaches a bound during the search along , it is fixed and is reduced for the next iteration.

There are two sets of convergence criteria – a weaker and a stronger. Whenever the weaker criteria are satisfied, the Lagrange multipliers are estimated for all the active constraints. If any Lagrange multiplier estimate is significantly negative, then one of the variables associated with a negative Lagrange multiplier estimate is released from its bound and the next search direction is computed in the extended subspace (i.e., is increased). Otherwise minimization continues in the current subspace provided that this is practicable. When it is not, or when the stronger convergence criteria are already satisfied, then, if one or more Lagrange multiplier estimates are close to zero, a slight perturbation is made in the values of the corresponding variables in turn until a lower function value is obtained. The normal algorithm is then resumed from the perturbed point.

If a saddle point is suspected, a local search is carried out with a view to moving away from the saddle point. A local search is also performed when a point is found which is thought to be a constrained minimum.

References

Gill, P E and Murray, W, 1976, Minimization subject to bounds on the variables, NPL Report NAC 72, National Physical Laboratory