L-BFGS-B  3.0
Large-scale Bound-constrained Optimization
active.f
Go to the documentation of this file.
1 c> \file active.f
2 
3 c> \brief This subroutine initializes iwhere and projects the initial x to
4 c> the feasible set if necessary.
5 c>
6 c> This subroutine initializes iwhere and projects the initial x to
7 c> the feasible set if necessary.
8 c>
9 c> @param n number of parameters
10 c> @param l lower bounds on parameters
11 c> @param u upper bounds on parameters
12 c> @param nbd indicates which bounds are present
13 c> @param x position
14 c> @param iwhere On entry iwhere is unspecified.<br/>
15 c> On exit: iwhere(i)=<ul><li>-1 if x(i) has no bounds</li>
16 c> <li> 3 if l(i)=u(i),</li>
17 c> <li> 0 otherwise.</li></ul>
18 c> In cauchy, iwhere is given finer gradations.
19 c> @param iprint console output flag
20 c> @param prjctd TODO
21 c> @param cnstnd TODO
22 c> @param boxed TODO
23 c
24 c * * *
25 c
26 c NEOS, November 1994. (Latest revision June 1996.)
27 c Optimization Technology Center.
28 c Argonne National Laboratory and Northwestern University.
29 c Written by
30 c Ciyou Zhu
31 c in collaboration with R.H. Byrd, P. Lu-Chen and J. Nocedal.
32  subroutine active(n, l, u, nbd, x, iwhere, iprint,
33  + prjctd, cnstnd, boxed)
34 
35  logical prjctd, cnstnd, boxed
36  integer n, iprint, nbd(n), iwhere(n)
37  double precision x(n), l(n), u(n)
38 
39 c ************
40 
41  integer nbdd,i
42  double precision zero
43  parameter(zero=0.0d0)
44 
45 c Initialize nbdd, prjctd, cnstnd and boxed.
46 
47  nbdd = 0
48  prjctd = .false.
49  cnstnd = .false.
50  boxed = .true.
51 
52 c Project the initial x to the feasible set if necessary.
53 
54  do 10 i = 1, n
55  if (nbd(i) .gt. 0) then
56  if (nbd(i) .le. 2 .and. x(i) .le. l(i)) then
57  if (x(i) .lt. l(i)) then
58  prjctd = .true.
59  x(i) = l(i)
60  endif
61  nbdd = nbdd + 1
62  else if (nbd(i) .ge. 2 .and. x(i) .ge. u(i)) then
63  if (x(i) .gt. u(i)) then
64  prjctd = .true.
65  x(i) = u(i)
66  endif
67  nbdd = nbdd + 1
68  endif
69  endif
70  10 continue
71 
72 c Initialize iwhere and assign values to cnstnd and boxed.
73 
74  do 20 i = 1, n
75  if (nbd(i) .ne. 2) boxed = .false.
76  if (nbd(i) .eq. 0) then
77 c this variable is always free
78  iwhere(i) = -1
79 
80 c otherwise set x(i)=mid(x(i), u(i), l(i)).
81  else
82  cnstnd = .true.
83  if (nbd(i) .eq. 2 .and. u(i) - l(i) .le. zero) then
84 c this variable is always fixed
85  iwhere(i) = 3
86  else
87  iwhere(i) = 0
88  endif
89  endif
90  20 continue
91 
92  if (iprint .ge. 0) then
93  if (prjctd) write (6,*)
94  + 'The initial X is infeasible. Restart with its projection.'
95  if (.not. cnstnd)
96  + write (6,*) 'This problem is unconstrained.'
97  endif
98 
99  if (iprint .gt. 0) write (6,1001) nbdd
100 
101  1001 format (/,'At X0 ',i9,' variables are exactly at the bounds')
102 
103  return
104 
105  end
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.
Definition: active.f:34