LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
dtzrqf.f
Go to the documentation of this file.
1 *> \brief \b DTZRQF
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DTZRQF + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtzrqf.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtzrqf.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtzrqf.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DTZRQF( M, N, A, LDA, TAU, INFO )
22 *
23 * .. Scalar Arguments ..
24 * INTEGER INFO, LDA, M, N
25 * ..
26 * .. Array Arguments ..
27 * DOUBLE PRECISION A( LDA, * ), TAU( * )
28 * ..
29 *
30 *
31 *> \par Purpose:
32 * =============
33 *>
34 *> \verbatim
35 *>
36 *> This routine is deprecated and has been replaced by routine DTZRZF.
37 *>
38 *> DTZRQF reduces the M-by-N ( M<=N ) real upper trapezoidal matrix A
39 *> to upper triangular form by means of orthogonal transformations.
40 *>
41 *> The upper trapezoidal matrix A is factored as
42 *>
43 *> A = ( R 0 ) * Z,
44 *>
45 *> where Z is an N-by-N orthogonal matrix and R is an M-by-M upper
46 *> triangular matrix.
47 *> \endverbatim
48 *
49 * Arguments:
50 * ==========
51 *
52 *> \param[in] M
53 *> \verbatim
54 *> M is INTEGER
55 *> The number of rows of the matrix A. M >= 0.
56 *> \endverbatim
57 *>
58 *> \param[in] N
59 *> \verbatim
60 *> N is INTEGER
61 *> The number of columns of the matrix A. N >= M.
62 *> \endverbatim
63 *>
64 *> \param[in,out] A
65 *> \verbatim
66 *> A is DOUBLE PRECISION array, dimension (LDA,N)
67 *> On entry, the leading M-by-N upper trapezoidal part of the
68 *> array A must contain the matrix to be factorized.
69 *> On exit, the leading M-by-M upper triangular part of A
70 *> contains the upper triangular matrix R, and elements M+1 to
71 *> N of the first M rows of A, with the array TAU, represent the
72 *> orthogonal matrix Z as a product of M elementary reflectors.
73 *> \endverbatim
74 *>
75 *> \param[in] LDA
76 *> \verbatim
77 *> LDA is INTEGER
78 *> The leading dimension of the array A. LDA >= max(1,M).
79 *> \endverbatim
80 *>
81 *> \param[out] TAU
82 *> \verbatim
83 *> TAU is DOUBLE PRECISION array, dimension (M)
84 *> The scalar factors of the elementary reflectors.
85 *> \endverbatim
86 *>
87 *> \param[out] INFO
88 *> \verbatim
89 *> INFO is INTEGER
90 *> = 0: successful exit
91 *> < 0: if INFO = -i, the i-th argument had an illegal value
92 *> \endverbatim
93 *
94 * Authors:
95 * ========
96 *
97 *> \author Univ. of Tennessee
98 *> \author Univ. of California Berkeley
99 *> \author Univ. of Colorado Denver
100 *> \author NAG Ltd.
101 *
102 *> \date November 2011
103 *
104 *> \ingroup doubleOTHERcomputational
105 *
106 *> \par Further Details:
107 * =====================
108 *>
109 *> \verbatim
110 *>
111 *> The factorization is obtained by Householder's method. The kth
112 *> transformation matrix, Z( k ), which is used to introduce zeros into
113 *> the ( m - k + 1 )th row of A, is given in the form
114 *>
115 *> Z( k ) = ( I 0 ),
116 *> ( 0 T( k ) )
117 *>
118 *> where
119 *>
120 *> T( k ) = I - tau*u( k )*u( k )**T, u( k ) = ( 1 ),
121 *> ( 0 )
122 *> ( z( k ) )
123 *>
124 *> tau is a scalar and z( k ) is an ( n - m ) element vector.
125 *> tau and z( k ) are chosen to annihilate the elements of the kth row
126 *> of X.
127 *>
128 *> The scalar tau is returned in the kth element of TAU and the vector
129 *> u( k ) in the kth row of A, such that the elements of z( k ) are
130 *> in a( k, m + 1 ), ..., a( k, n ). The elements of R are returned in
131 *> the upper triangular part of A.
132 *>
133 *> Z is given by
134 *>
135 *> Z = Z( 1 ) * Z( 2 ) * ... * Z( m ).
136 *> \endverbatim
137 *>
138 * =====================================================================
139  SUBROUTINE dtzrqf( M, N, A, LDA, TAU, INFO )
140 *
141 * -- LAPACK computational routine (version 3.4.0) --
142 * -- LAPACK is a software package provided by Univ. of Tennessee, --
143 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
144 * November 2011
145 *
146 * .. Scalar Arguments ..
147  INTEGER info, lda, m, n
148 * ..
149 * .. Array Arguments ..
150  DOUBLE PRECISION a( lda, * ), tau( * )
151 * ..
152 *
153 * =====================================================================
154 *
155 * .. Parameters ..
156  DOUBLE PRECISION one, zero
157  parameter( one = 1.0d+0, zero = 0.0d+0 )
158 * ..
159 * .. Local Scalars ..
160  INTEGER i, k, m1
161 * ..
162 * .. Intrinsic Functions ..
163  INTRINSIC max, min
164 * ..
165 * .. External Subroutines ..
166  EXTERNAL daxpy, dcopy, dgemv, dger, dlarfg, xerbla
167 * ..
168 * .. Executable Statements ..
169 *
170 * Test the input parameters.
171 *
172  info = 0
173  IF( m.LT.0 ) THEN
174  info = -1
175  ELSE IF( n.LT.m ) THEN
176  info = -2
177  ELSE IF( lda.LT.max( 1, m ) ) THEN
178  info = -4
179  END IF
180  IF( info.NE.0 ) THEN
181  CALL xerbla( 'DTZRQF', -info )
182  return
183  END IF
184 *
185 * Perform the factorization.
186 *
187  IF( m.EQ.0 )
188  $ return
189  IF( m.EQ.n ) THEN
190  DO 10 i = 1, n
191  tau( i ) = zero
192  10 continue
193  ELSE
194  m1 = min( m+1, n )
195  DO 20 k = m, 1, -1
196 *
197 * Use a Householder reflection to zero the kth row of A.
198 * First set up the reflection.
199 *
200  CALL dlarfg( n-m+1, a( k, k ), a( k, m1 ), lda, tau( k ) )
201 *
202  IF( ( tau( k ).NE.zero ) .AND. ( k.GT.1 ) ) THEN
203 *
204 * We now perform the operation A := A*P( k ).
205 *
206 * Use the first ( k - 1 ) elements of TAU to store a( k ),
207 * where a( k ) consists of the first ( k - 1 ) elements of
208 * the kth column of A. Also let B denote the first
209 * ( k - 1 ) rows of the last ( n - m ) columns of A.
210 *
211  CALL dcopy( k-1, a( 1, k ), 1, tau, 1 )
212 *
213 * Form w = a( k ) + B*z( k ) in TAU.
214 *
215  CALL dgemv( 'No transpose', k-1, n-m, one, a( 1, m1 ),
216  $ lda, a( k, m1 ), lda, one, tau, 1 )
217 *
218 * Now form a( k ) := a( k ) - tau*w
219 * and B := B - tau*w*z( k )**T.
220 *
221  CALL daxpy( k-1, -tau( k ), tau, 1, a( 1, k ), 1 )
222  CALL dger( k-1, n-m, -tau( k ), tau, 1, a( k, m1 ), lda,
223  $ a( 1, m1 ), lda )
224  END IF
225  20 continue
226  END IF
227 *
228  return
229 *
230 * End of DTZRQF
231 *
232  END