newuoa Subroutine

public subroutine newuoa(n, npt, x, rhobeg, rhoend, iprint, maxfun, calfun)

This subroutine seeks the least value of a function of many variables, by a trust region method that forms quadratic models by interpolation. There can be some freedom in the interpolation conditions, which is taken up by minimizing the Frobenius norm of the change to the second derivative of the quadratic model, beginning with a zero matrix.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

the number of variables. must be at least 2.

integer, intent(in) :: npt

The number of interpolation conditions. Its value must be in the interval [N+2,(N+1)(N+2)/2].

real(kind=wp), intent(inout), dimension(*) :: x

Initial values of the variables must be set in X(1),X(2),...,X(N). They will be changed to the values that give the least calculated F.

real(kind=wp), intent(in) :: rhobeg

RHOBEG and RHOEND must be set to the initial and final values of a trust region radius, so both must be positive with RHOEND<=RHOBEG. Typically RHOBEG should be about one tenth of the greatest expected change to a variable, and RHOEND should indicate the accuracy that is required in the final values of the variables.

real(kind=wp), intent(in) :: rhoend

RHOBEG and RHOEND must be set to the initial and final values of a trust region radius, so both must be positive with RHOEND<=RHOBEG. Typically RHOBEG should be about one tenth of the greatest expected change to a variable, and RHOEND should indicate the accuracy that is required in the final values of the variables.

integer, intent(in) :: iprint

The value of IPRINT should be set to 0, 1, 2 or 3, which controls the amount of printing. Specifically, there is no output if IPRINT=0 and there is output only at the return if IPRINT=1. Otherwise, each new value of RHO is printed, with the best vector of variables so far and the corresponding value of the objective function. Further, each new value of F with its variables are output if IPRINT=3.

integer, intent(in) :: maxfun

an upper bound on the number of calls of CALFUN.

procedure(func) :: calfun

It must set F to the value of the objective function for the variables X(1),X(2),...,X(N).


Calls

proc~~newuoa~~CallsGraph proc~newuoa 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

Called by

proc~~newuoa~~CalledByGraph proc~newuoa newuoa proc~newuoa_test newuoa_test proc~newuoa_test->proc~newuoa

Source Code

    subroutine newuoa (n, npt, x, rhobeg, rhoend, iprint, maxfun, calfun)
    
        implicit none
        
        integer,intent(in)                  :: n       !! the number of variables. must be at least 2.
        integer,intent(in)                  :: npt     !! The number of interpolation conditions.
                                                       !! Its value must be in the interval `[N+2,(N+1)(N+2)/2]`.
        real(wp),dimension(*),intent(inout) :: x       !! Initial values of the variables must be set in X(1),X(2),...,X(N). They
                                                       !! will be changed to the values that give the least calculated F.
        real(wp),intent(in)                 :: rhobeg  !! RHOBEG and RHOEND must be set to the initial and final values of a trust
                                                       !! region radius, so both must be positive with RHOEND<=RHOBEG. Typically
                                                       !! RHOBEG should be about one tenth of the greatest expected change to a
                                                       !! variable, and RHOEND should indicate the accuracy that is required in
                                                       !! the final values of the variables.
        real(wp),intent(in)                 :: rhoend  !! RHOBEG and RHOEND must be set to the initial and final values of a trust
                                                       !! region radius, so both must be positive with RHOEND<=RHOBEG. Typically
                                                       !! RHOBEG should be about one tenth of the greatest expected change to a
                                                       !! variable, and RHOEND should indicate the accuracy that is required in
                                                       !! the final values of the variables. 
        integer,intent(in)                  :: iprint  !! The value of IPRINT should be set to 0, 1, 2 or 3, which controls the
                                                       !! amount of printing. Specifically, there is no output if IPRINT=0 and
                                                       !! there is output only at the return if IPRINT=1. Otherwise, each new
                                                       !! value of RHO is printed, with the best vector of variables so far and
                                                       !! the corresponding value of the objective function. Further, each new
                                                       !! value of F with its variables are output if IPRINT=3.
        integer,intent(in)                  :: maxfun  !! an upper bound on the number of calls of CALFUN.
        procedure(func)                     :: calfun  !! It must set F to the value of the objective function 
                                                       !! for the variables `X(1),X(2),...,X(N)`.

        real(wp),dimension(:),allocatable :: w
        integer :: np,nptm,ndim,ixb,ixo,ixn,ixp,ifv,igq,ihq,ipq,ibmat,izmat,id,ivl,iw
        
        ! Partition the working space array, so that different parts of it can be
        ! treated separately by the subroutine that performs the main calculation.

        np = n + 1
        nptm = npt - np
        if (npt < n+2 .or. npt > ((n+2)*np)/2) then
            write(*,*) 'Return from NEWUOA because NPT is not in the required interval'
            return
        end if
        
        ! The array W will be used for working space
        allocate(w((NPT+13)*(NPT+N)+3*N*(N+3)/2))
        
        ndim = npt + n
        ixb = 1
        ixo = ixb + n
        ixn = ixo + n
        ixp = ixn + n
        ifv = ixp + n * npt
        igq = ifv + npt
        ihq = igq + n
        ipq = ihq + (n*np) / 2
        ibmat = ipq + npt
        izmat = ibmat + ndim * n
        id = izmat + npt * nptm
        ivl = id + n
        iw = ivl + ndim

        ! The above settings provide a partition of W for subroutine NEWUOB.
        ! The partition requires the first NPT*(NPT+N)+5*N*(N+3)/2 elements of
        ! W plus the space that is needed by the last array of NEWUOB.

        call newuob (n, npt, x, rhobeg, rhoend, iprint, maxfun, w(ixb), w(ixo), w(ixn), &
                     w(ixp), w(ifv), w(igq), w(ihq), w(ipq), w(ibmat), w(izmat), ndim, &
                     w(id), w(ivl), w(iw), calfun)

        deallocate(w)

    end subroutine newuoa