NLSL
qrsolv.f90
Go to the documentation of this file.
1  subroutine qrsolv(n,r,ldr,ipvt,diag,qtb,x,sdiag,wa)
2  integer n,ldr
3  integer ipvt(n)
4  double precision r(ldr,n),diag(n),qtb(n),x(n),sdiag(n),wa(n)
5 c **********
6 c
7 c subroutine qrsolv
8 c
9 c Given an m by n matrix A, an n by n diagonal matrix D,
10 c and an m-vector B, the problem is to determine an X which
11 c solves the system
12 c
13 c A*X = B , D*X = 0 ,
14 c
15 c in the least squares sense.
16 c
17 c This subroutine completes the solution of the problem
18 c if it is provided with the necessary information from the
19 c QR factorization, with column pivoting, of A. That is, if
20 c A*P = Q*R, where P is a permutationf matrix, Q has orthogonal
21 c columns, and R is an upper triangular matrix with diagonal
22 c elements of nonincreasing magnitude, then qrsolv expects
23 c the full upper triangle of R, the permutation matrix P,
24 c and the first n components of (Q transpose)*B. The system
25 c A*X = B, D*X = 0, is then equivalent to
26 c
27 c t t
28 c R*Z = Q *B , P *D*P*Z = 0 ,
29 c
30 c where X = P*Z. If this system does not have full rank,
31 c then a least squares solution is obtained. On output qrsolv
32 c also provides an upper triangular matrix S such that
33 c
34 c t t t
35 C P *(A *A + D*D)*P = S *S .
36 c
37 c S is computed within qrsolv and may be of separate interest.
38 c
39 c The subroutine statement is
40 c
41 c subroutine qrsolv(n,r,ldr,ipvt,diag,qtb,x,sdiag,wa)
42 c
43 c where
44 c
45 c n is a positive integer input variable set to the order of r.
46 c
47 c r is an n by n array. On input the full upper triangle
48 c must contain the full upper triangle of the matrix r.
49 c On output the full upper triangle is unaltered, and the
50 c strict lower triangle contains the strict upper triangle
51 c (transposed) of the upper triangular matrix S
52 c
53 c ldr is a positive integer input variable not less than n
54 c which specifies the leading dimension of the array r.
55 c
56 c ipvt is an integer input array of length n which defines the
57 c permutation matrix p such that A*P = Q*R. Column j of P
58 c is column ipvt(j) of the identity matrix.
59 c
60 c diag is an input array of length n which must contain the
61 c diagonal elements of the matrix D.
62 c
63 c qtb is an input array of length n which must contain the first
64 c n elements of the vector (Q transpose)*B.
65 c
66 c x is an output array of length n which contains the least
67 c squares solution of the system A*X = B, D*X = 0.
68 c
69 c sdiag is an output array of length n which contains the
70 c diagonal elements of the upper triangular matrix S.
71 c
72 c wa is a work array of length n.
73 c
74 c Subprograms called
75 c
76 c Fortran-supplied ... dabs,dsqrt
77 c
78 c Argonne National Laboratory. MINPACK project. March 1980.
79 c Burton S. Garbow, Kenneth E. Hillstrom, Jorge J. More
80 c
81 c **********
82  integer i,j,jp1,k,kp1,l,nsing
83  double precision cos,cotan,p5,p25,qtbpj,sin,sum,tan,temp,zero
84  data p5,p25,zero /5.0d-1,2.5d-1,0.0d0/
85 c
86 c----------------------------------------------------------------------
87 c Copy R and (Q transpose)*B to preserve input and initialize S.
88 c In particular, save the diagonal elements of R in X.
89 c----------------------------------------------------------------------
90  do 20 j = 1, n
91  do 10 i = j, n
92  r(i,j) = r(j,i)
93  10 continue
94  x(j) = r(j,j)
95  wa(j) = qtb(j)
96  20 continue
97 c
98 c----------------------------------------------------------------------
99 c Eliminate the diagonal matrix D using a Givens rotation.
100 c----------------------------------------------------------------------
101  do 100 j = 1, n
102 c
103 c Prepare the row of D to be eliminated, locating the
104 c diagonal element using P from the QR factorization.
105 c
106  l = ipvt(j)
107  if (diag(l) .ne. zero) then
108  do 30 k = j, n
109  sdiag(k) = zero
110  30 continue
111  sdiag(j) = diag(l)
112 c
113 c The transformations to eliminate the row of D
114 c modify only a single element of (Q transpose)*B
115 c beyond the first n, which is initially zero.
116 c
117  qtbpj = zero
118  do 80 k = j, n
119 c
120 c Determine a Givens rotation which eliminates the
121 c appropriate element in the current row of D
122 c
123  if (sdiag(k) .ne. zero) then
124  if (dabs(r(k,k)) .ge. dabs(sdiag(k))) then
125  tan = sdiag(k)/r(k,k)
126  cos = p5/dsqrt(p25+p25*tan**2)
127  sin = cos*tan
128  else
129  cotan = r(k,k)/sdiag(k)
130  sin = p5/dsqrt(p25+p25*cotan**2)
131  cos = sin*cotan
132  end if
133 c
134 c Compute the modified diagonal element of R and
135 c the modified element of ((Q transpose)*B,0).
136 c
137  r(k,k) = cos*r(k,k) + sin*sdiag(k)
138  temp = cos*wa(k) + sin*qtbpj
139  qtbpj = -sin*wa(k) + cos*qtbpj
140  wa(k) = temp
141 c
142 c Accumulate the tranformation in the row of S.
143 c
144  kp1 = k + 1
145  if (kp1 .le. n) then
146  do 60 i = kp1, n
147  temp = cos*r(i,k) + sin*sdiag(i)
148  sdiag(i) = -sin*r(i,k) + cos*sdiag(i)
149  r(i,k) = temp
150  60 continue
151  end if
152  end if
153  80 continue
154  end if
155 c
156 c Store the diagonal element of S and restore
157 c the corresponding diagonal element of r.
158 c
159  sdiag(j) = r(j,j)
160  r(j,j) = x(j)
161  100 continue
162 c
163 c Solve the triangular system for z. If the system is
164 c singular, then obtain a least squares solution.
165 c
166  nsing = n
167  do 110 j = 1, n
168  if (sdiag(j) .eq. zero .and. nsing .eq. n) nsing = j - 1
169  if (nsing .lt. n) wa(j) = zero
170  110 continue
171 c
172  if (nsing .ge. 1) then
173  do 140 k = 1, nsing
174  j = nsing - k + 1
175  sum = zero
176  jp1 = j + 1
177  if (jp1 .le. nsing) then
178  do 120 i = jp1, nsing
179  sum = sum + r(i,j)*wa(i)
180  120 continue
181  end if
182  wa(j) = (wa(j) - sum)/sdiag(j)
183  140 continue
184  end if
185 c
186 c Permute the components of Z back to components of X.
187 c
188  do 160 j = 1, n
189  l = ipvt(j)
190  x(l) = wa(j)
191  160 continue
192 c
193  return
194  end
subroutine qrsolv(n, r, ldr, ipvt, diag, qtb, x, sdiag, wa)
Definition: qrsolv.f90:2