projgr Subroutine

private subroutine projgr(n, l, u, Nbd, x, g, Sbgnrm)

This subroutine computes the infinity norm of the projected gradient.

Credits

  • NEOS, November 1994. (Latest revision June 1996.) Optimization Technology Center. Argonne National Laboratory and Northwestern University. Written by Ciyou Zhu in collaboration with R.H. Byrd, P. Lu-Chen and J. Nocedal.

Arguments

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

the dimension of the problem (the number of variables).

real(kind=wp), intent(in) :: l(n)

the lower bound on x.

real(kind=wp), intent(in) :: u(n)

the upper bound on x.

integer, intent(in) :: Nbd(n)
real(kind=wp), intent(in) :: x(n)
real(kind=wp), intent(in) :: g(n)
real(kind=wp), intent(out) :: Sbgnrm

infinity norm of the projected gradient


Called by

proc~~projgr~~CalledByGraph proc~projgr lbfgsb_module::projgr proc~mainlb lbfgsb_module::mainlb proc~mainlb->proc~projgr proc~setulb lbfgsb_module::setulb proc~setulb->proc~mainlb

Source Code

      subroutine projgr(n,l,u,Nbd,x,g,Sbgnrm)

      implicit none

      integer,intent(in) :: n !! the dimension of the problem (the number of variables).
      integer,intent(in) :: Nbd(n)
      real(wp),intent(out) :: Sbgnrm !! infinity norm of the projected gradient
      real(wp),intent(in) :: x(n)
      real(wp),intent(in) :: l(n) !! the lower bound on `x`.
      real(wp),intent(in) :: u(n) !! the upper bound on `x`.
      real(wp),intent(in) :: g(n)

      integer :: i
      real(wp) :: gi

      Sbgnrm = zero
      do i = 1 , n
         gi = g(i)
         if ( Nbd(i)/=0 ) then
            if ( gi<zero ) then
               if ( Nbd(i)>=2 ) gi = max((x(i)-u(i)),gi)
            else
               if ( Nbd(i)<=2 ) gi = min((x(i)-l(i)),gi)
            endif
         endif
         Sbgnrm = max(Sbgnrm,abs(gi))
      enddo

      end subroutine projgr