First steps

Using NLPModels standard modeling interface

NLPModelsAlgencan is naturally built as an NLPModels interface. Hence, one needs an implementation of an AbstractNLPModel to call the Algencan solver.

For instance, the NLPModels package provides an automatic differentiation NLPModel, called ADNLPModel, using ForwardDiff to compute the derivatives. For further details, it is recommended to access this tutorial. As an example, the following optimization problem is defined as:

\[\begin{aligned} \min \quad & f(x_1, x_2) = x_1x_2 + 5 \\ \textrm{s.t.} \quad & x_1+x_2 \leq 5, \\ & x_1^2+x_2^2=10, \\ & 0 \leq x_1, x_2 \leq 5, \end{aligned}\]

And could be modeled as an ADNLModel as follows:

using NLPModels, NLPModelsAlgencan

x0   = [1.0; 1.0]
f(x) = x[1]*x[2] + 5
lvar = [0.0; 0.0]
uvar = [5.0; 5.0]
c(x) = [x[1] + x[2]; x[1]^2 + x[2]^2]
lcon = [-Inf; 10.0]
ucon = [5.0;  10.0]

nlp  = ADNLPModel(f, x0, lvar, uvar, c, lcon, ucon)

Moreover, it is solved with Algencan and displayed the problem statistics, which is provided by a return instance of a GenericExecutionStats:

stats = algencan(nlp)
print(stats)

 ==============================================================================
 This is ALGENCAN 3.1.1.
 ALGENCAN, an Augmented Lagrangian method for nonlinear programming, is part of
 the TANGO Project: Trustable Algorithms for Nonlinear General Optimization.
 See http://www.ime.usp.br/~egbirgin/tango/ for details.
 ==============================================================================

 There are no strings to be processed in the array of parameters.

 The specification file is not being used.

 Available HSL subroutines = NONE

 ALGENCAN PARAMETERS:

 firstde                =                    T
 seconde                =                    T
 truehpr                =                    T
 hptype in TN           =               TRUEHP
 lsslvr in TR           =            NONE/NONE
 lsslvr in NW           =            NONE/NONE
 lsslvr in ACCPROC      =            NONE/NONE
 innslvr                =                   TN
 accproc                =                    F
 rmfixv                 =                    T
 slacks                 =                    F
 scale                  =                    T
 epsfeas                =           1.0000D-08
 epsopt                 =           1.0000D-08
 efstain                =           1.0000D-04
 eostain                =           1.0000D-12
 efacc                  =           1.0000D-04
 eoacc                  =           1.0000D-04
 iprint                 =                   10
 ncomp                  =                    6

 Specification filename =                   ''
 Output filename        =                   ''
 Solution filename      =                   ''

 Number of variables               :       2
 Number of equality constraints    :       1
 Number of inequality constraints  :       1
 Number of bound constraints       :       4
 Number of fixed variables         :       0

 There are no fixed variables to be removed.

 Objective function scale factor   : 1.0D+00
 Smallest constraints scale factor : 5.0D-01

 Entry to ALGENCAN.
 Number of variables  :       2
 Number of constraints:       2

out penalt  objective infeas  scaled    scaled infeas norm   |Grad| inner Newton
ite         function  ibilty  obj-funct infeas +compl graLag infeas totit forKKT
  0         6.000D+00 8.D+00  6.000D+00 4.D+00 4.D+00 1.D+00 4.D+00    0   0   0
  1 8.D+00  9.867D+00 3.D-01  9.867D+00 1.D-01 1.D-01 1.D-08 3.D-01    5C  0   0
  2 1.D+02  1.000D+01 1.D-10  1.000D+01 6.D-11 6.D-11 5.D-11 1.D-10    8C  0   0

 Flag of ALGENCAN: Solution was found.

 User-provided subroutines calls counters:

 Subroutine fsub     (coded=F):        0
 Subroutine gsub     (coded=F):        0
 Subroutine hsub     (coded=F):        0
 Subroutine csub     (coded=F):        0 (       0 calls per constraint in avg)
 Subroutine jacsub   (coded=F):        0 (       0 calls per constraint in avg)
 Subroutine hcsub    (coded=F):        0 (       0 calls per constraint in avg)
 Subroutine fcsub    (coded=T):       18
 Subroutine gjacsub  (coded=T):       16
 Subroutine gjacpsub (coded=F):        0
 Subroutine hlsub    (coded=T):        8
 Subroutine hlpsub   (coded=F):        0


 Total CPU time in seconds =     3.25
Generic Execution stats
  status: first-order stationary
  objective value: 9.999999999942318
  primal feasibility: 1.1536371857800987e-10
  dual feasibility: 5.768185928900493e-11
  solution: [2.2360679774868917  2.2360679774868917]
  multipliers: [0.0  -0.499999999989881]
  iterations: -1
  elapsed time: 3.292833567

Using JuMP modeling interface

Alternatively, one can use the JuMP interface to model the problem and then solve it with Algencan. This solution is provided by the NLPModelsJuMP package. To model the same problem as before, it could be done as follows:

using JuMP, NLPModelsJuMP, NLPModelsAlgencan

model = Model()

@variable(model, 0 ≤ x[1:2] ≤ 5)
set_start_value(x[1], 1.0)
set_start_value(x[2], 1.0)

@NLobjective(model, Min, x[1]*x[2] + 5)

@constraint(model, x[1] + x[2] ≤ 5)
@NLconstraint(model, x[1]^2 + x[2]^2 == 10)

nlp = MathOptNLPModel(model)

stats = algencan(nlp)
print(stats)

 ==============================================================================
 This is ALGENCAN 3.1.1.
 ALGENCAN, an Augmented Lagrangian method for nonlinear programming, is part of
 the TANGO Project: Trustable Algorithms for Nonlinear General Optimization.
 See http://www.ime.usp.br/~egbirgin/tango/ for details.
 ==============================================================================

 There are no strings to be processed in the array of parameters.

 The specification file is not being used.

 Available HSL subroutines = NONE

 ALGENCAN PARAMETERS:

 firstde                =                    T
 seconde                =                    T
 truehpr                =                    T
 hptype in TN           =               TRUEHP
 lsslvr in TR           =            NONE/NONE
 lsslvr in NW           =            NONE/NONE
 lsslvr in ACCPROC      =            NONE/NONE
 innslvr                =                   TN
 accproc                =                    F
 rmfixv                 =                    T
 slacks                 =                    F
 scale                  =                    T
 epsfeas                =           1.0000D-08
 epsopt                 =           1.0000D-08
 efstain                =           1.0000D-04
 eostain                =           1.0000D-12
 efacc                  =           1.0000D-04
 eoacc                  =           1.0000D-04
 iprint                 =                   10
 ncomp                  =                    6

 Specification filename =                   ''
 Output filename        =                   ''
 Solution filename      =                   ''

 Number of variables               :       2
 Number of equality constraints    :       1
 Number of inequality constraints  :       1
 Number of bound constraints       :       4
 Number of fixed variables         :       0

 There are no fixed variables to be removed.

 Objective function scale factor   : 1.0D+00
 Smallest constraints scale factor : 5.0D-01

 Entry to ALGENCAN.
 Number of variables  :       2
 Number of constraints:       2

out penalt  objective infeas  scaled    scaled infeas norm   |Grad| inner Newton
ite         function  ibilty  obj-funct infeas +compl graLag infeas totit forKKT
  0         6.000D+00 8.D+00  6.000D+00 4.D+00 4.D+00 1.D+00 4.D+00    0   0   0
  1 8.D+00  9.867D+00 3.D-01  9.867D+00 1.D-01 1.D-01 1.D-08 3.D-01    5C  0   0
  2 1.D+02  1.000D+01 1.D-10  1.000D+01 6.D-11 6.D-11 5.D-11 1.D-10    8C  0   0

 Flag of ALGENCAN: Solution was found.

 User-provided subroutines calls counters:

 Subroutine fsub     (coded=F):        0
 Subroutine gsub     (coded=F):        0
 Subroutine hsub     (coded=F):        0
 Subroutine csub     (coded=F):        0 (       0 calls per constraint in avg)
 Subroutine jacsub   (coded=F):        0 (       0 calls per constraint in avg)
 Subroutine hcsub    (coded=F):        0 (       0 calls per constraint in avg)
 Subroutine fcsub    (coded=T):       18
 Subroutine gjacsub  (coded=T):       16
 Subroutine gjacpsub (coded=F):        0
 Subroutine hlsub    (coded=T):        8
 Subroutine hlpsub   (coded=F):        0


 Total CPU time in seconds =     2.39
Generic Execution stats
  status: first-order stationary
  objective value: 9.999999999942318
  primal feasibility: 1.1536371857800987e-10
  dual feasibility: 5.768185928900493e-11
  solution: [2.2360679774868917  2.2360679774868917]
  multipliers: [0.0  -0.499999999989881]
  iterations: -1
  elapsed time: 2.435221733