chiark / gitweb /
strip
[nlopt.git] / test / t_fortran.f90
1 program main
2   external myfunc, myconstraint
3   double precision lb(2)
4   integer*8 opt
5   double precision d1(2), d2(2)
6   double precision x(2), minf
7   integer ires
8   include 'nlopt.f'
9
10   opt=0
11   call nlo_create(opt, NLOPT_LD_MMA, 2)
12   call nlo_get_lower_bounds(ires, opt, lb)
13   lb(2) = 0.0
14   call nlo_set_lower_bounds(ires, opt, lb)
15   call nlo_set_min_objective(ires, opt, myfunc, 0)
16
17   d1(1) = 2.
18   d1(2) = 0.
19   call nlo_add_inequality_constraint(ires, opt, myconstraint, d1, 1.D-8)
20   d2(1) = -1.
21   d2(2) = 1.
22   call nlo_add_inequality_constraint(ires, opt, myconstraint, d2, 1.D-8)
23
24   call nlo_set_xtol_rel(ires, opt, 1.D-4)
25
26   x(1) = 1.234
27   x(2) = 5.678
28   call nlo_optimize(ires, opt, x, minf)
29   if (ires.lt.0) then
30     write(*,*) 'nlopt failed!'
31     stop 1
32   else
33     write(*,*) 'found min at ', x(1), x(2)
34     write(*,*) 'min val = ', minf
35   endif
36
37   call nlo_destroy(opt)
38
39   end 
40
41   subroutine myfunc(val, n, x, grad, need_gradient, f_data)
42   double precision val, x(n), grad(n)
43   integer n, need_gradient
44   if (need_gradient.ne.0) then
45      grad(1) = 0.0
46      grad(2) = 0.5 / dsqrt(x(2))
47   endif
48   val = dsqrt(x(2))
49   end 
50
51   subroutine myconstraint(val, n, x, grad, need_gradient, d)
52   integer need_gradient
53   double precision val, x(n), grad(n), d(2), a, b
54   a = d(1)
55   b = d(2)
56   if (need_gradient.ne.0) then
57     grad(1) = 3. * a * (a*x(1) + b)**2
58     grad(2) = -1.0
59   endif
60   val = (a*x(1) + b)**3 - x(2)
61   end