newuoa_test Subroutine

public subroutine newuoa_test()

The Chebyquad test problem (Fletcher, 1965) for N = 2,4,6 and 8, with NPT = 2N+1.

Arguments

None

Calls

proc~~newuoa_test~~CallsGraph proc~newuoa_test newuoa_test proc~newuoa newuoa proc~newuoa_test->proc~newuoa proc~newuob newuob proc~newuoa->proc~newuob proc~bigden bigden proc~newuob->proc~bigden proc~biglag biglag proc~newuob->proc~biglag proc~trsapp trsapp proc~newuob->proc~trsapp proc~update update proc~newuob->proc~update den den proc~bigden->den denex denex proc~bigden->denex par par proc~bigden->par

Source Code

    subroutine newuoa_test ()

        implicit none
        
        real(wp) :: x (10)
        integer :: iprint,maxfun,n,npt,i
        real(wp) :: rhoend,rhobeg
        
        iprint = 2
        maxfun = 5000
        rhoend = 1.0e-6_wp
        do n = 2, 8, 2
            npt = 2 * n + 1
            do i = 1, n
                x (i) = real (i, wp) / real (n+1, wp)
            end do
            rhobeg = 0.2_wp * x (1)
            print 20, n, npt
20          format (/ / 4 x, 'Results with N =', i2, ' and NPT =', i3)
            call newuoa (n, npt, x, rhobeg, rhoend, iprint, maxfun, calfun)
        end do
 
    contains
 
        subroutine calfun (n, x, f)
        
            implicit none
            
            integer :: n
            real (wp) :: x (*)
            real (wp) :: f
            
            real(wp) :: y (10, 10)
            integer :: i,j,np,iw
            real(wp) :: sum
            
            do j = 1, n
                y (1, j) = 1.0_wp
                y (2, j) = 2.0_wp * x (j) - 1.0_wp
            end do
            do i = 2, n
                do j = 1, n
                    y (i+1, j) = 2.0_wp * y (2, j) * y (i, j) - y (i-1, j)
                end do
            end do
            f = 0.0_wp
            np = n + 1
            iw = 1
            do i = 1, np
                sum = 0.0_wp
                do j = 1, n
                    sum = sum + y (i, j)
                end do
                sum = sum / real (n, wp)
                if (iw > 0) sum = sum + 1.0_wp / real (i*i-2*i, wp)
                iw = - iw
                f = f + sum * sum
            end do
        end subroutine calfun
 
    end subroutine newuoa_test