cobyla_test Subroutine

public subroutine cobyla_test()

Test routine for cobyla.

From: Report DAMTP 1992/NA5.

Arguments

None

Calls

proc~~cobyla_test~~CallsGraph proc~cobyla_test cobyla_test proc~cobyla cobyla proc~cobyla_test->proc~cobyla proc~cobylb cobylb proc~cobyla->proc~cobylb proc~trstlp trstlp proc~cobylb->proc~trstlp

Source Code

    subroutine cobyla_test ()

        implicit none

        real(wp),dimension(10) :: x,xopt
        integer :: nprob,n,m,i,icase,iprint,maxfun
        real(wp) :: rhobeg,rhoend,temp,tempa,tempb,tempc,tempd

        do nprob = 1, 10

            if (nprob == 1) then
    !
    !     minimization of a simple quadratic function of two variables.
    !
                print 10
10              format (/ 7 x, 'Output from test problem 1 (Simple quadratic)')
                n = 2
                m = 0
                xopt (1) = - 1.0_wp
                xopt (2) = 0.0_wp

            else if (nprob == 2) then
    !
    !     Easy two dimensional minimization in unit circle.
    !
                print 20
20              format (/ 7 x, 'Output from test problem 2 (2D unit circle ',&
                'calculation)')
                n = 2
                m = 1
                xopt (1) = sqrt (0.5_wp)
                xopt (2) = - xopt (1)

            else if (nprob == 3) then
    !
    !     Easy three dimensional minimization in ellipsoid.
    !
                print 30
30              format (/ 7 x, 'Output from test problem 3 (3D ellipsoid ',&
                'calculation)')
                n = 3
                m = 1
                xopt (1) = 1.0_wp / sqrt (3.0_wp)
                xopt (2) = 1.0_wp / sqrt (6.0_wp)
                xopt (3) = - 1.0_wp / 3.0_wp

            else if (nprob == 4) then
    !
    !     Weak version of Rosenbrock's problem.
    !
                print 40
40              format (/ 7 x, 'Output from test problem 4 (Weak Rosenbrock)')
                n = 2
                m = 0
                xopt (1) = - 1.0_wp
                xopt (2) = 1.0_wp

            else if (nprob == 5) then
    !
    !     Intermediate version of Rosenbrock's problem.
    !
                print 50
50              format (/ 7 x, 'Output from test problem 5 (Intermediate ', 'Rosenbrock)')
                n = 2
                m = 0
                xopt (1) = - 1.0_wp
                xopt (2) = 1.0_wp

            else if (nprob == 6) then
    !
    !     This problem is taken from Fletcher's book Practical Methods of
    !     Optimization and has the equation number (9.1.15).
    !
                print 60
60              format (/ 7 x, 'Output from test problem 6 (Equation ',&
                '(9.1.15) in Fletcher)')
                n = 2
                m = 2
                xopt (1) = sqrt (0.5_wp)
                xopt (2) = xopt (1)

            else if (nprob == 7) then
    !
    !     This problem is taken from Fletcher's book Practical Methods of
    !     Optimization and has the equation number (14.4.2).
    !
                print 70
70              format (/ 7 x, 'Output from test problem 7 (Equation ',&
                '(14.4.2) in Fletcher)')
                n = 3
                m = 3
                xopt (1) = 0.0_wp
                xopt (2) = - 3.0_wp
                xopt (3) = - 3.0_wp

            else if (nprob == 8) then
    !
    !     This problem is taken from page 66 of Hock and Schittkowski's book Test
    !     Examples for Nonlinear Programming Codes. It is their test problem Number
    !     43, and has the name Rosen-Suzuki.
    !
                print 80
80              format (/ 7 x, 'Output from test problem 8 (Rosen-Suzuki)')
                n = 4
                m = 3
                xopt (1) = 0.0_wp
                xopt (2) = 1.0_wp
                xopt (3) = 2.0_wp
                xopt (4) = - 1.0_wp

            else if (nprob == 9) then
    !
    !     This problem is taken from page 111 of Hock and Schittkowski's
    !     book Test Examples for Nonlinear Programming Codes. It is their
    !     test problem Number 100.
    !
                print 90
90              format (/ 7 x, 'Output from test problem 9 (Hock and ',&
                'Schittkowski 100)')
                n = 7
                m = 4
                xopt (1) = 2.330499_wp
                xopt (2) = 1.951372_wp
                xopt (3) = - 0.4775414_wp
                xopt (4) = 4.365726_wp
                xopt (5) = - 0.624487_wp
                xopt (6) = 1.038131_wp
                xopt (7) = 1.594227_wp

            else if (nprob == 10) then
    !
    !     This problem is taken from page 415 of Luenberger's book Applied
    !     Nonlinear Programming. It is to maximize the area of a hexagon of
    !     unit diameter.
    !
                print 100
100             format (/ 7 x, 'Output from test problem 10 (Hexagon area)')
                n = 9
                m = 14
            end if
            do icase = 1, 2
                do i = 1, n
                    x (i) = 1.0_wp
                end do
                rhobeg = 0.5_wp
                rhoend = 0.001_wp
                if (icase == 2) rhoend = 0.0001_wp
                iprint = 1
                maxfun = 2000
                call cobyla (n, m, x, rhobeg, rhoend, iprint, maxfun, calcfc)
                if (nprob == 10) then
                    tempa = x (1) + x (3) + x (5) + x (7)
                    tempb = x (2) + x (4) + x (6) + x (8)
                    tempc = 0.5_wp / sqrt (tempa*tempa+tempb*tempb)
                    tempd = tempc * sqrt (3.0_wp)
                    xopt (1) = tempd * tempa + tempc * tempb
                    xopt (2) = tempd * tempb - tempc * tempa
                    xopt (3) = tempd * tempa - tempc * tempb
                    xopt (4) = tempd * tempb + tempc * tempa
                    do i = 1, 4
                        xopt (i+4) = xopt (i)
                    end do
                    xopt (9) = 0.0_wp
                end if
                temp = 0.0_wp
                do i = 1, n
                    temp = temp + (x(i)-xopt(i)) ** 2
                end do
                print 150, sqrt (temp)
150             format (/ 5 x, 'Least squares error in variables =', 1 pe16.6)
            end do
            print 170
170         format (2 x, '----------------------------------------------',&
            '--------------------')

        end do

    contains

        subroutine calcfc (n, m, x, f, con)

            implicit none

            integer,intent(in)                :: n
            integer,intent(in)                :: m
            real(wp),dimension(*),intent(in)  :: x
            real(wp),intent(out)              :: f
            real(wp),dimension(*),intent(out) :: con

            if (nprob == 1) then
    !
    !     Test problem 1 (Simple quadratic)
    !
                f = 10.0_wp * (x(1)+1.0_wp) ** 2 + x (2) ** 2

            else if (nprob == 2) then
    !
    !    Test problem 2 (2D unit circle calculation)
    !
                f = x (1) * x (2)
                con (1) = 1.0_wp - x (1) ** 2 - x (2) ** 2

            else if (nprob == 3) then
    !
    !     Test problem 3 (3D ellipsoid calculation)
    !
                f = x (1) * x (2) * x (3)
                con (1) = 1.0_wp - x (1) ** 2 - 2.0_wp * x (2) ** 2 - 3.0_wp * x (3) ** 2

            else if (nprob == 4) then
    !
    !     Test problem 4 (Weak Rosenbrock)
    !
                f = (x(1)**2-x(2)) ** 2 + (1.0_wp+x(1)) ** 2

            else if (nprob == 5) then
    !
    !     Test problem 5 (Intermediate Rosenbrock)
    !
                f = 10.0_wp * (x(1)**2-x(2)) ** 2 + (1.0_wp+x(1)) ** 2

            else if (nprob == 6) then
    !
    !     Test problem 6 (Equation (9.1.15) in Fletcher's book)
    !
                f = - x (1) - x (2)
                con (1) = x (2) - x (1) ** 2
                con (2) = 1.0_wp - x (1) ** 2 - x (2) ** 2

            else if (nprob == 7) then
    !
    !     Test problem 7 (Equation (14.4.2) in Fletcher's book)
    !
                f = x (3)
                con (1) = 5.0_wp * x (1) - x (2) + x (3)
                con (2) = x (3) - x (1) ** 2 - x (2) ** 2 - 4.0_wp * x (2)
                con (3) = x (3) - 5.0_wp * x (1) - x (2)

            else if (nprob == 8) then
    !
    !     Test problem 8 (Rosen-Suzuki)
    !
                f = x (1) ** 2 + x (2) ** 2 + 2.0_wp * x (3) ** 2 + x (4) ** 2 - 5.0_wp * &
               & x (1) - 5.0_wp * x (2) - 21.0_wp * x (3) + 7.0_wp * x (4)
                con (1) = 8.0_wp - x (1) ** 2 - x (2) ** 2 - x (3) ** 2 - x (4) ** 2 - x &
               & (1) + x (2) - x (3) + x (4)
                con (2) = 10.0_wp - x (1) ** 2 - 2.0_wp * x (2) ** 2 - x (3) ** 2 - &
               & 2.0_wp * x (4) ** 2 + x (1) + x (4)
                con (3) = 5.0_wp - 2.0_wp * x (1) ** 2 - x (2) ** 2 - x (3) ** 2 - 2.0_wp &
               & * x (1) + x (2) + x (4)

            else if (nprob == 9) then
    !
    !     Test problem 9 (Hock and Schittkowski 100)
    !
                f = (x(1)-10.0_wp) ** 2 + 5.0_wp * (x(2)-12.0_wp) ** 2 + x (3) ** 4 + &
               & 3.0_wp * (x(4)-11.0_wp) ** 2 + 10.0_wp * x (5) ** 6 + 7.0_wp * x (6) ** &
               & 2 + x (7) ** 4 - 4.0_wp * x (6) * x (7) - 10.0_wp * x (6) - 8.0_wp * x &
               & (7)
                con (1) = 127.0_wp - 2.0_wp * x (1) ** 2 - 3.0_wp * x (2) ** 4 - x (3) - &
               & 4.0_wp * x (4) ** 2 - 5.0_wp * x (5)
                con (2) = 282.0_wp - 7.0_wp * x (1) - 3.0_wp * x (2) - 10.0_wp * x (3) ** &
               & 2 - x (4) + x (5)
                con (3) = 196.0_wp - 23.0_wp * x (1) - x (2) ** 2 - 6.0_wp * x (6) ** 2 + &
               & 8.0_wp * x (7)
                con (4) = - 4.0_wp * x (1) ** 2 - x (2) ** 2 + 3.0_wp * x (1) * x (2) - &
               & 2.0_wp * x (3) ** 2 - 5.0_wp * x (6) + 11.0_wp * x (7)

            else if (nprob == 10) then
    !
    !     Test problem 10 (Hexagon area)
    !
                f = - 0.5_wp * &
                   (x(1)*x(4)-x(2)*x(3)+x(3)*x(9)-x(5)*x(9)+x(5)*x(8)-x(6)*x(7))
                con (1) = 1.0_wp - x (3) ** 2 - x (4) ** 2
                con (2) = 1.0_wp - x (9) ** 2
                con (3) = 1.0_wp - x (5) ** 2 - x (6) ** 2
                con (4) = 1.0_wp - x (1) ** 2 - (x(2)-x(9)) ** 2
                con (5) = 1.0_wp - (x(1)-x(5)) ** 2 - (x(2)-x(6)) ** 2
                con (6) = 1.0_wp - (x(1)-x(7)) ** 2 - (x(2)-x(8)) ** 2
                con (7) = 1.0_wp - (x(3)-x(5)) ** 2 - (x(4)-x(6)) ** 2
                con (8) = 1.0_wp - (x(3)-x(7)) ** 2 - (x(4)-x(8)) ** 2
                con (9) = 1.0_wp - x (7) ** 2 - (x(8)-x(9)) ** 2
                con (10) = x (1) * x (4) - x (2) * x (3)
                con (11) = x (3) * x (9)
                con (12) = - x (5) * x (9)
                con (13) = x (5) * x (8) - x (6) * x (7)
                con (14) = x (9)

            end if

        end subroutine calcfc

    end subroutine cobyla_test