NLSL
qrfac.f90
Go to the documentation of this file.
1  subroutine qrfac(m,n,a,lda,pivot,ipvt,lipvt,rdiag,acnorm,wa)
2  use rnddbl
3  integer m,n,lda,lipvt
4  integer ipvt(n)
5  logical pivot
6  double precision a(lda,n),rdiag(n),acnorm(n),wa(n)
7 c **********
8 c
9 c subroutine qrfac
10 c
11 c This subroutine uses Householder transformations with column
12 c pivoting (optional) to compute a QR factorization of the
13 c m by n matrix A. That is, qrfac determines an orthogonal
14 c matrix Q, a permutation matrix P, and an upper trapezoidal
15 c matrix R with diagonal elements of nonincreasing magnitude,
16 c such that A*P = Q*R. The Householder transformation for
17 c column k, k = 1,2,...,min(m,n), is of the form
18 c
19 c T
20 c i - (1/u(k))*u*u
21 c
22 c where u has zeros in the first k-1 positions. The form of
23 c this transformation and the method of pivoting first
24 c appeared in the corresponding LINPACK subroutine.
25 c
26 c The subroutine statement is
27 c
28 c subroutine qrfac(m,n,a,lda,pivot,ipvt,lipvt,rdiag,acnorm,wa)
29 c
30 c where
31 c
32 c m is a positive integer input variable set to the number
33 c of rows of a.
34 c
35 c n is a positive integer input variable set to the number
36 c of columns of a.
37 c
38 c a is an m by n array. On input a contains the matrix for
39 c which the qr factorization is to be computed. On output
40 c the strict upper trapezoidal part of A contains the strict
41 c upper trapezoidal part of R, and the lower trapezoidal
42 c part of a contains a factored form of Q (the non-trivial
43 c elements of the u vectors described above).
44 c
45 c lda is a positive integer input variable not less than m
46 c which specifies the leading dimension of the array a.
47 c
48 c pivot is a logical input variable. If pivot is set true,
49 c then column pivoting is enforced. If pivot is set false,
50 c then no column pivoting is done.
51 c
52 c ipvt is an integer output array of length n. ipvt
53 c defines the permutation matrix P such that A*P = Q*R.
54 c Column j of P is column ipvt(j) of the identity matrix.
55 c If pivot is false, ipvt is not referenced.
56 c
57 c lipvt is a positive integer input variable. If pivot is true
58 c then lipvt specifies the number of leading columns in the
59 c Jacobian that may be pivoted. Columns lipvt+1 through n are
60 c fixed at the right hand side of the matrix.
61 c
62 c rdiag is an output array of length n which contains the
63 c diagonal elements of r.
64 c
65 c acnorm is an output array of length n which contains the
66 c norms of the corresponding columns of the input matrix a.
67 c if this information is not needed, then acnorm can coincide
68 c with rdiag.
69 c
70 c wa is a work array of length n. if pivot is false, then wa
71 c can coincide with rdiag.
72 c
73 c Subprograms called
74 c
75 c MINPACK-supplied ... enorm
76 c
77 c Fortran-supplied ... dmax1,dsqrt,min0
78 c
79 c Argonne National Laboratory. MINPACK project. March 1980.
80 c Burton S. Garbow, Kenneth E. Hillstrom, Jorge J. More
81 c
82 c **********
83  integer i,j,jp1,k,kmax,minmn
84  double precision ajnorm,epsmch,one,p05,sum,temp,zero
85  parameter(epsmch=rndoff)
86  double precision enorm
87  data one,p05,zero /1.0d0,5.0d-2,0.0d0/
88 c
89 c epsmch is the machine precision.
90 c
91  if (lipvt .gt. n) lipvt = n
92 c
93 c----------------------------------------------------------------------
94 c Compute the initial column norms and initialize several arrays.
95 c----------------------------------------------------------------------
96  do j = 1, n
97  acnorm(j) = enorm(m,a(1,j))
98  rdiag(j) = acnorm(j)
99  wa(j) = rdiag(j)
100  if (pivot) ipvt(j) = j
101  end do
102 c
103 c----------------------------------------------------------------------
104 c Reduce A to R with Householder transformations.
105 c----------------------------------------------------------------------
106  minmn = min0(m,n)
107  do j = 1, minmn
108  if (pivot) then
109 c
110 c----------------------------------------------------------------------
111 c Bring the column of largest norm into the pivot position.
112 c----------------------------------------------------------------------
113  kmax = j
114  do k = j, lipvt
115  if (rdiag(k) .gt. rdiag(kmax)) kmax = k
116  end do
117 c
118  if (j.ne.kmax) then
119  do i = 1, m
120  temp = a(i,j)
121  a(i,j) = a(i,kmax)
122  a(i,kmax) = temp
123  end do
124 c
125  rdiag(kmax) = rdiag(j)
126  wa(kmax) = wa(j)
127  k = ipvt(j)
128  ipvt(j) = ipvt(kmax)
129  ipvt(kmax) = k
130  end if
131  end if
132 c
133 c----------------------------------------------------------------------
134 c Compute the Householder transformation to reduce the
135 c j-th column of a to a multiple of the j-th unit vector.
136 c----------------------------------------------------------------------
137  ajnorm = enorm(m-j+1,a(j,j))
138 c
139  if (ajnorm .ne. zero) then
140  if (a(j,j) .lt. zero) ajnorm = -ajnorm
141  do i = j, m
142  a(i,j) = a(i,j)/ajnorm
143  end do
144 c
145  a(j,j) = a(j,j) + one
146 c
147 c----------------------------------------------------------------------
148 c Apply the transformation to the remaining columns
149 c and update the norms.
150 c----------------------------------------------------------------------
151  jp1 = j + 1
152  if (jp1.le.n) then
153  do k = jp1, n
154  sum = zero
155  do i = j, m
156  sum = sum + a(i,j)*a(i,k)
157  end do
158 c
159  temp = sum/a(j,j)
160  do i = j, m
161  a(i,k) = a(i,k) - temp*a(i,j)
162  end do
163 c
164  if (pivot .and. rdiag(k) .ne. zero) then
165  temp = a(j,k)/rdiag(k)
166  rdiag(k) = rdiag(k)*dsqrt(dmax1(zero,one-temp**2))
167  if (p05*(rdiag(k)/wa(k))**2 .le. epsmch) then
168  rdiag(k) = enorm(m-j,a(jp1,k))
169  wa(k) = rdiag(k)
170  end if
171  end if
172  end do
173 c
174  end if
175  end if
176 c
177  rdiag(j) = -ajnorm
178  end do
179 c
180  return
181  end
subroutine qrfac(m, n, a, lda, pivot, ipvt, lipvt, rdiag, acnorm, wa)
Definition: qrfac.f90:2
double precision, parameter rndoff
Definition: rnddbl.f90:86