API

As stated in the home page, we consider the nonlinear optimization problem in the following format:

\[\begin{align*} \min \quad & f(x) \\ & c_L \leq c(x) \leq c_U \\ & \ell \leq x \leq u. \end{align*}\]

To develop an optimization algorithm, we are usually worried not only with $f(x)$ and $c(x)$, but also with their derivatives. Namely,

There are many ways to access some of these values, so here is a little reference guide.

Reference guide

The following naming should be easy enough to follow. If not, click on the link and go to the description.

Feel free to open an issue to suggest other methods that should apply to all NLPModels instances.

FunctionNLPModels function
$f(x)$obj, objgrad, objgrad!, objcons, objcons!
$\nabla f(x)$grad, grad!, objgrad, objgrad!
$\nabla^2 f(x)$hess, hess_op, hess_op!, hess_coord, hprod, hprod!
$c(x)$cons, cons!, objcons, objcons!
$J(x)$jac, jac_op, jac_op!, jac_coord, jprod, jprod!, jtprod, jtprod!
$\nabla^2 L(x,y)$hess, hess_op, hess_coord, hprod, hprod!

AbstractNLPModel functions

NLPModels.objFunction.

f = obj(nlp, x)

Evaluate $f(x)$, the objective function of nlp at x.

source
NLPModels.gradFunction.

g = grad(nlp, x)

Evaluate $\nabla f(x)$, the gradient of the objective function at x.

source
NLPModels.grad!Function.

g = grad!(nlp, x, g)

Evaluate $\nabla f(x)$, the gradient of the objective function at x in place.

source
NLPModels.objgradFunction.

f, g = objgrad(nlp, x)

Evaluate $f(x)$ and $\nabla f(x)$ at x.

source
NLPModels.objgrad!Function.

f, g = objgrad!(nlp, x, g)

Evaluate $f(x)$ and $\nabla f(x)$ at x. g is overwritten with the value of $\nabla f(x)$.

source
NLPModels.consFunction.

c = cons(nlp, x)

Evaluate $c(x)$, the constraints at x.

source
NLPModels.cons!Function.

c = cons!(nlp, x, c)

Evaluate $c(x)$, the constraints at x in place.

source
NLPModels.objconsFunction.

f, c = objcons(nlp, x)

Evaluate $f(x)$ and $c(x)$ at x.

source
NLPModels.objcons!Function.

f = objcons!(nlp, x, c)

Evaluate $f(x)$ and $c(x)$ at x. c is overwritten with the value of $c(x)$.

source
NLPModels.jac_coordFunction.

(rows,cols,vals) = jac_coord(nlp, x)

Evaluate $\nabla c(x)$, the constraint's Jacobian at x in sparse coordinate format.

source
NLPModels.jacFunction.

Jx = jac(nlp, x)

Evaluate $\nabla c(x)$, the constraint's Jacobian at x as a sparse matrix.

source
NLPModels.jac_opFunction.

J = jac_op(nlp, x)

Return the Jacobian at x as a linear operator. The resulting object may be used as if it were a matrix, e.g., J * v or J' * v.

source
NLPModels.jac_op!Function.

J = jac_op!(nlp, x, Jv, Jtv)

Return the Jacobian at x as a linear operator. The resulting object may be used as if it were a matrix, e.g., J * v or J' * v. The values Jv and Jtv are used as preallocated storage for the operations.

source
NLPModels.jprodFunction.

Jv = jprod(nlp, x, v)

Evaluate $\nabla c(x)v$, the Jacobian-vector product at x.

source
NLPModels.jprod!Function.

Jv = jprod!(nlp, x, v, Jv)

Evaluate $\nabla c(x)v$, the Jacobian-vector product at x in place.

source
NLPModels.jtprodFunction.

Jtv = jtprod(nlp, x, v, Jtv)

Evaluate $\nabla c(x)^Tv$, the transposed-Jacobian-vector product at x.

source
NLPModels.jtprod!Function.

Jtv = jtprod!(nlp, x, v, Jtv)

Evaluate $\nabla c(x)^Tv$, the transposed-Jacobian-vector product at x in place.

source
NLPModels.hess_coordFunction.

(rows,cols,vals) = hess_coord(nlp, x; obj_weight=1.0, y=zeros)

Evaluate the Lagrangian Hessian at (x,y) in sparse coordinate format, with objective function scaled by obj_weight, i.e.,

\[ \nabla^2L(x,y) = \sigma * \nabla^2 f(x) + \sum_{i=1}^m y_i\nabla^2 c_i(x), \]

with σ = obj_weight. Only the lower triangle is returned.

source
NLPModels.hessFunction.

Hx = hess(nlp, x; obj_weight=1.0, y=zeros)

Evaluate the Lagrangian Hessian at (x,y) as a sparse matrix, with objective function scaled by obj_weight, i.e.,

\[ \nabla^2L(x,y) = \sigma * \nabla^2 f(x) + \sum_{i=1}^m y_i\nabla^2 c_i(x), \]

with σ = obj_weight. Only the lower triangle is returned.

source
NLPModels.hess_opFunction.

H = hess_op(nlp, x; obj_weight=1.0, y=zeros)

Return the Lagrangian Hessian at (x,y) with objective function scaled by obj_weight as a linear operator. The resulting object may be used as if it were a matrix, e.g., H * v. The linear operator H represents

\[ \nabla^2L(x,y) = \sigma * \nabla^2 f(x) + \sum_{i=1}^m y_i\nabla^2 c_i(x), \]

with σ = obj_weight.

source
NLPModels.hess_op!Function.

H = hess_op!(nlp, x, Hv; obj_weight=1.0, y=zeros)

Return the Lagrangian Hessian at (x,y) with objective function scaled by obj_weight as a linear operator, and storing the result on Hv. The resulting object may be used as if it were a matrix, e.g., w = H * v. The vector Hv is used as preallocated storage for the operation. The linear operator H represents

\[ \nabla^2L(x,y) = \sigma * \nabla^2 f(x) + \sum_{i=1}^m y_i\nabla^2 c_i(x), \]

with σ = obj_weight.

source
NLPModels.hprodFunction.

Hv = hprod(nlp, x, v; obj_weight=1.0, y=zeros)

Evaluate the product of the Lagrangian Hessian at (x,y) with the vector v, with objective function scaled by obj_weight, i.e.,

\[ \nabla^2L(x,y) = \sigma * \nabla^2 f(x) + \sum_{i=1}^m y_i\nabla^2 c_i(x), \]

with σ = obj_weight.

source
NLPModels.hprod!Function.

Hv = hprod!(nlp, x, v, Hv; obj_weight=1.0, y=zeros)

Evaluate the product of the Lagrangian Hessian at (x,y) with the vector v in place, with objective function scaled by obj_weight, i.e.,

\[ \nabla^2L(x,y) = \sigma * \nabla^2 f(x) + \sum_{i=1}^m y_i\nabla^2 c_i(x), \]

with σ = obj_weight.

source
NLPModels.NLPtoMPBFunction.
mp = NLPtoMPB(nlp, solver)

Return a MathProgBase model corresponding to an AbstractNLPModel.

Arguments

  • nlp::AbstractNLPModel

  • solver::AbstractMathProgSolver a solver instance, e.g., IpoptSolver()

Currently, all models are treated as nonlinear models.

Return values

The function returns a MathProgBase model mpbmodel such that it should be possible to call

MathProgBase.optimize!(mpbmodel)
source

reset!(counters)

Reset evaluation counters

source

`reset!(nlp)

Reset evaluation count in nlp

source

Derivative check

Check the first derivatives of the objective at x against centered finite differences.

This function returns a dictionary indexed by components of the gradient for which the relative error exceeds rtol.

source

Check the first derivatives of the constraints at x against centered finite differences.

This function returns a dictionary indexed by (j, i) tuples such that the relative error in the i-th partial derivative of the j-th constraint exceeds rtol.

source

Check the second derivatives of the objective and each constraints at x against centered finite differences. This check does not rely on exactness of the first derivatives, only on objective and constraint values.

The sgn arguments refers to the formulation of the Lagrangian in the problem. It should have a positive value if the Lagrangian is formulated as

L(x,y) = f(x) + ∑ yⱼ cⱼ(x)

e.g., as in JuMPNLPModels, and a negative value if the Lagrangian is formulated as

L(x,y) = f(x) - ∑ yⱼ cⱼ(x)

e.g., as in AmplModels. Only the sign of sgn is important.

This function returns a dictionary indexed by functions. The 0-th function is the objective while the k-th function (for k > 0) is the k-th constraint. The values of the dictionary are dictionaries indexed by tuples (i, j) such that the relative error in the second derivative ∂²fₖ/∂xᵢ∂xⱼ exceeds rtol.

source

Check the second derivatives of the objective and each constraints at x against centered finite differences. This check assumes exactness of the first derivatives.

The sgn arguments refers to the formulation of the Lagrangian in the problem. It should have a positive value if the Lagrangian is formulated as

L(x,y) = f(x) + ∑ yⱼ cⱼ(x)

e.g., as in JuMPNLPModels, and a negative value if the Lagrangian is formulated as

L(x,y) = f(x) - ∑ yⱼ cⱼ(x)

e.g., as in AmplModels. Only the sign of sgn is important.

This function returns a dictionary indexed by functions. The 0-th function is the objective while the k-th function (for k > 0) is the k-th constraint. The values of the dictionary are dictionaries indexed by tuples (i, j) such that the relative error in the second derivative ∂²fₖ/∂xᵢ∂xⱼ exceeds rtol.

source