L-BFGS-B  3.0
Large-scale Bound-constrained Optimization
bmv.f
Go to the documentation of this file.
1 c> \file bmv.f
2 
3 c> \brief This subroutine computes the product of the 2m x 2m middle matrix
4 c> in the compact L-BFGS formula of B and a 2m vector v.
5 c>
6 c> This subroutine computes the product of the 2m x 2m middle matrix
7 c> in the compact L-BFGS formula of B and a 2m vector v;
8 c> it returns the product in p.
9 c>
10 c> @param m On entry m is the maximum number of variable metric corrections
11 c> used to define the limited memory matrix.<br/>
12 c> On exit m is unchanged.
13 c>
14 c> @param sy On entry sy specifies the matrix S'Y.<br/>
15 c> On exit sy is unchanged.
16 c>
17 c> @param wt On entry wt specifies the upper triangular matrix J' which is
18 c> the Cholesky factor of (thetaS'S+LD^(-1)L').<br/>
19 c> On exit wt is unchanged.
20 c>
21 c> @param col On entry col specifies the number of s-vectors (or y-vectors)
22 c> stored in the compact L-BFGS formula.<br/>
23 c> On exit col is unchanged.
24 c>
25 c> @param v On entry v specifies vector v.<br/>
26 c> On exit v is unchanged.
27 c>
28 c> @param p On entry p is unspecified.<br/>
29 c> On exit p is the product Mv.
30 c>
31 c> @param info On entry info is unspecified.<br/>
32 c> On exit info = <ul><li>0 for normal return,</li>
33 c> <li>nonzero for abnormal return when the system
34 c> to be solved by dtrsl is singular.</li></ul>
35  subroutine bmv(m, sy, wt, col, v, p, info)
36 
37  integer m, col, info
38  double precision sy(m, m), wt(m, m), v(2*col), p(2*col)
39 
40 c
41 c * * *
42 c
43 c NEOS, November 1994. (Latest revision June 1996.)
44 c Optimization Technology Center.
45 c Argonne National Laboratory and Northwestern University.
46 c Written by
47 c Ciyou Zhu
48 c in collaboration with R.H. Byrd, P. Lu-Chen and J. Nocedal.
49 c
50 c
51 c ************
52 
53  integer i,k,i2
54  double precision sum
55 
56  if (col .eq. 0) return
57 
58 c PART I: solve [ D^(1/2) O ] [ p1 ] = [ v1 ]
59 c [ -L*D^(-1/2) J ] [ p2 ] [ v2 ].
60 
61 c solve Jp2=v2+LD^(-1)v1.
62  p(col + 1) = v(col + 1)
63  do 20 i = 2, col
64  i2 = col + i
65  sum = 0.0d0
66  do 10 k = 1, i - 1
67  sum = sum + sy(i,k)*v(k)/sy(k,k)
68  10 continue
69  p(i2) = v(i2) + sum
70  20 continue
71 c Solve the triangular system
72  !call dtrsl(wt,m,col,p(col+1),11,info)
73  call dtrsm('l','u','t','n',col,1,one,wt,m,p(col+1),col)
74  info = 0
75  !if (info .ne. 0) return
76 
77 c solve D^(1/2)p1=v1.
78  do 30 i = 1, col
79  p(i) = v(i)/sqrt(sy(i,i))
80  30 continue
81 
82 c PART II: solve [ -D^(1/2) D^(-1/2)*L' ] [ p1 ] = [ p1 ]
83 c [ 0 J' ] [ p2 ] [ p2 ].
84 
85 c solve J^Tp2=p2.
86  !call dtrsl(wt,m,col,p(col+1),01,info)
87  !if (info .ne. 0) return
88  call dtrsm('l','u','n','n',col,1,one,wt,m,p(col+1),col)
89  info = 0
90 
91 c compute p1=-D^(-1/2)(p1-D^(-1/2)L'p2)
92 c =-D^(-1/2)p1+D^(-1)L'p2.
93  do 40 i = 1, col
94  p(i) = -p(i)/sqrt(sy(i,i))
95  40 continue
96  do 60 i = 1, col
97  sum = 0.d0
98  do 50 k = i + 1, col
99  sum = sum + sy(k,i)*p(col+k)/sy(i,i)
100  50 continue
101  p(i) = p(i) + sum
102  60 continue
103 
104  return
105 
106  end
subroutine bmv(m, sy, wt, col, v, p, info)
This subroutine computes the product of the 2m x 2m middle matrix in the compact L-BFGS formula of B ...
Definition: bmv.f:36