active Subroutine

private subroutine active(n, l, u, Nbd, x, Iwhere, Iprint, Prjctd, Cnstnd, Boxed)

This subroutine initializes iwhere and projects the initial x to the feasible set if necessary.

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(inout) :: x(n)
integer, intent(out) :: Iwhere(n)
  • iwhere(i)=-1 if x(i) has no bounds
  • iwhere(i)=3 if l(i)=u(i)
  • iwhere(i)=0 otherwise.

In cauchy, iwhere is given finer gradations.

integer, intent(in) :: Iprint

Controls the frequency and type of output generated

logical, intent(out) :: Prjctd
logical, intent(out) :: Cnstnd
logical, intent(out) :: Boxed

Called by

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

Source Code

      subroutine active(n,l,u,Nbd,x,Iwhere,Iprint,Prjctd,Cnstnd,Boxed)
      implicit none

      integer,intent(in) :: n !! the dimension of the problem (the number of variables).
      real(wp),intent(in) :: l(n) !! the lower bound on `x`.
      real(wp),intent(in) :: u(n) !! the upper bound on `x`.
      integer,intent(in) :: Nbd(n)
      real(wp),intent(inout) :: x(n)
      integer,intent(out) :: Iwhere(n) !! * `iwhere(i)=-1` if `x(i)` has no bounds
                                       !! * `iwhere(i)=3`  if `l(i)=u(i)`
                                       !! * `iwhere(i)=0`  otherwise.
                                       !!
                                       !! In [[cauchy]], `iwhere` is given finer gradations.
      integer,intent(in) :: Iprint !! Controls the frequency and type of output generated
      logical,intent(out) :: Prjctd
      logical,intent(out) :: Cnstnd
      logical,intent(out) :: Boxed

      integer nbdd , i

      ! Initialize nbdd, prjctd, cnstnd and boxed.

      nbdd = 0
      Prjctd = .false.
      Cnstnd = .false.
      Boxed = .true.

      ! Project the initial x to the feasible set if necessary.

      do i = 1 , n
         if ( Nbd(i)>0 ) then
            if ( Nbd(i)<=2 .and. x(i)<=l(i) ) then
               if ( x(i)<l(i) ) then
                  Prjctd = .true.
                  x(i) = l(i)
               endif
               nbdd = nbdd + 1
            else if ( Nbd(i)>=2 .and. x(i)>=u(i) ) then
               if ( x(i)>u(i) ) then
                  Prjctd = .true.
                  x(i) = u(i)
               endif
               nbdd = nbdd + 1
            endif
         endif
      enddo

      ! Initialize iwhere and assign values to cnstnd and boxed.

      do i = 1 , n
         if ( Nbd(i)/=2 ) Boxed = .false.
         if ( Nbd(i)==0 ) then
            ! this variable is always free
            Iwhere(i) = -1
            ! otherwise set x(i)=mid(x(i), u(i), l(i)).
         else
            Cnstnd = .true.
            if ( Nbd(i)==2 .and. u(i)-l(i)<=zero ) then
               ! this variable is always fixed
               Iwhere(i) = 3
            else
               Iwhere(i) = 0
            endif
         endif
      enddo

      if ( Iprint>=0 ) then
         if ( Prjctd ) write (output_unit,*) &
            'The initial X is infeasible.  Restart with its projection.'
         if ( .not.Cnstnd ) write (output_unit,*) 'This problem is unconstrained.'
      endif

      if ( Iprint>0 ) write (output_unit,'(/,a,i9,a)') &
            'At X0 ', nbdd, ' variables are exactly at the bounds'

      end subroutine active