LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
dqrt02.f
Go to the documentation of this file.
1 *> \brief \b DQRT02
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE DQRT02( M, N, K, A, AF, Q, R, LDA, TAU, WORK, LWORK,
12 * RWORK, RESULT )
13 *
14 * .. Scalar Arguments ..
15 * INTEGER K, LDA, LWORK, M, N
16 * ..
17 * .. Array Arguments ..
18 * DOUBLE PRECISION A( LDA, * ), AF( LDA, * ), Q( LDA, * ),
19 * $ R( LDA, * ), RESULT( * ), RWORK( * ), TAU( * ),
20 * $ WORK( LWORK )
21 * ..
22 *
23 *
24 *> \par Purpose:
25 * =============
26 *>
27 *> \verbatim
28 *>
29 *> DQRT02 tests DORGQR, which generates an m-by-n matrix Q with
30 *> orthonornmal columns that is defined as the product of k elementary
31 *> reflectors.
32 *>
33 *> Given the QR factorization of an m-by-n matrix A, DQRT02 generates
34 *> the orthogonal matrix Q defined by the factorization of the first k
35 *> columns of A; it compares R(1:n,1:k) with Q(1:m,1:n)'*A(1:m,1:k),
36 *> and checks that the columns of Q are orthonormal.
37 *> \endverbatim
38 *
39 * Arguments:
40 * ==========
41 *
42 *> \param[in] M
43 *> \verbatim
44 *> M is INTEGER
45 *> The number of rows of the matrix Q to be generated. M >= 0.
46 *> \endverbatim
47 *>
48 *> \param[in] N
49 *> \verbatim
50 *> N is INTEGER
51 *> The number of columns of the matrix Q to be generated.
52 *> M >= N >= 0.
53 *> \endverbatim
54 *>
55 *> \param[in] K
56 *> \verbatim
57 *> K is INTEGER
58 *> The number of elementary reflectors whose product defines the
59 *> matrix Q. N >= K >= 0.
60 *> \endverbatim
61 *>
62 *> \param[in] A
63 *> \verbatim
64 *> A is DOUBLE PRECISION array, dimension (LDA,N)
65 *> The m-by-n matrix A which was factorized by DQRT01.
66 *> \endverbatim
67 *>
68 *> \param[in] AF
69 *> \verbatim
70 *> AF is DOUBLE PRECISION array, dimension (LDA,N)
71 *> Details of the QR factorization of A, as returned by DGEQRF.
72 *> See DGEQRF for further details.
73 *> \endverbatim
74 *>
75 *> \param[out] Q
76 *> \verbatim
77 *> Q is DOUBLE PRECISION array, dimension (LDA,N)
78 *> \endverbatim
79 *>
80 *> \param[out] R
81 *> \verbatim
82 *> R is DOUBLE PRECISION array, dimension (LDA,N)
83 *> \endverbatim
84 *>
85 *> \param[in] LDA
86 *> \verbatim
87 *> LDA is INTEGER
88 *> The leading dimension of the arrays A, AF, Q and R. LDA >= M.
89 *> \endverbatim
90 *>
91 *> \param[in] TAU
92 *> \verbatim
93 *> TAU is DOUBLE PRECISION array, dimension (N)
94 *> The scalar factors of the elementary reflectors corresponding
95 *> to the QR factorization in AF.
96 *> \endverbatim
97 *>
98 *> \param[out] WORK
99 *> \verbatim
100 *> WORK is DOUBLE PRECISION array, dimension (LWORK)
101 *> \endverbatim
102 *>
103 *> \param[in] LWORK
104 *> \verbatim
105 *> LWORK is INTEGER
106 *> The dimension of the array WORK.
107 *> \endverbatim
108 *>
109 *> \param[out] RWORK
110 *> \verbatim
111 *> RWORK is DOUBLE PRECISION array, dimension (M)
112 *> \endverbatim
113 *>
114 *> \param[out] RESULT
115 *> \verbatim
116 *> RESULT is DOUBLE PRECISION array, dimension (2)
117 *> The test ratios:
118 *> RESULT(1) = norm( R - Q'*A ) / ( M * norm(A) * EPS )
119 *> RESULT(2) = norm( I - Q'*Q ) / ( M * EPS )
120 *> \endverbatim
121 *
122 * Authors:
123 * ========
124 *
125 *> \author Univ. of Tennessee
126 *> \author Univ. of California Berkeley
127 *> \author Univ. of Colorado Denver
128 *> \author NAG Ltd.
129 *
130 *> \date November 2011
131 *
132 *> \ingroup double_lin
133 *
134 * =====================================================================
135  SUBROUTINE dqrt02( M, N, K, A, AF, Q, R, LDA, TAU, WORK, LWORK,
136  $ rwork, result )
137 *
138 * -- LAPACK test routine (version 3.4.0) --
139 * -- LAPACK is a software package provided by Univ. of Tennessee, --
140 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
141 * November 2011
142 *
143 * .. Scalar Arguments ..
144  INTEGER k, lda, lwork, m, n
145 * ..
146 * .. Array Arguments ..
147  DOUBLE PRECISION a( lda, * ), af( lda, * ), q( lda, * ),
148  $ r( lda, * ), result( * ), rwork( * ), tau( * ),
149  $ work( lwork )
150 * ..
151 *
152 * =====================================================================
153 *
154 * .. Parameters ..
155  DOUBLE PRECISION zero, one
156  parameter( zero = 0.0d+0, one = 1.0d+0 )
157  DOUBLE PRECISION rogue
158  parameter( rogue = -1.0d+10 )
159 * ..
160 * .. Local Scalars ..
161  INTEGER info
162  DOUBLE PRECISION anorm, eps, resid
163 * ..
164 * .. External Functions ..
165  DOUBLE PRECISION dlamch, dlange, dlansy
166  EXTERNAL dlamch, dlange, dlansy
167 * ..
168 * .. External Subroutines ..
169  EXTERNAL dgemm, dlacpy, dlaset, dorgqr, dsyrk
170 * ..
171 * .. Intrinsic Functions ..
172  INTRINSIC dble, max
173 * ..
174 * .. Scalars in Common ..
175  CHARACTER*32 srnamt
176 * ..
177 * .. Common blocks ..
178  common / srnamc / srnamt
179 * ..
180 * .. Executable Statements ..
181 *
182  eps = dlamch( 'Epsilon' )
183 *
184 * Copy the first k columns of the factorization to the array Q
185 *
186  CALL dlaset( 'Full', m, n, rogue, rogue, q, lda )
187  CALL dlacpy( 'Lower', m-1, k, af( 2, 1 ), lda, q( 2, 1 ), lda )
188 *
189 * Generate the first n columns of the matrix Q
190 *
191  srnamt = 'DORGQR'
192  CALL dorgqr( m, n, k, q, lda, tau, work, lwork, info )
193 *
194 * Copy R(1:n,1:k)
195 *
196  CALL dlaset( 'Full', n, k, zero, zero, r, lda )
197  CALL dlacpy( 'Upper', n, k, af, lda, r, lda )
198 *
199 * Compute R(1:n,1:k) - Q(1:m,1:n)' * A(1:m,1:k)
200 *
201  CALL dgemm( 'Transpose', 'No transpose', n, k, m, -one, q, lda, a,
202  $ lda, one, r, lda )
203 *
204 * Compute norm( R - Q'*A ) / ( M * norm(A) * EPS ) .
205 *
206  anorm = dlange( '1', m, k, a, lda, rwork )
207  resid = dlange( '1', n, k, r, lda, rwork )
208  IF( anorm.GT.zero ) THEN
209  result( 1 ) = ( ( resid / dble( max( 1, m ) ) ) / anorm ) / eps
210  ELSE
211  result( 1 ) = zero
212  END IF
213 *
214 * Compute I - Q'*Q
215 *
216  CALL dlaset( 'Full', n, n, zero, one, r, lda )
217  CALL dsyrk( 'Upper', 'Transpose', n, m, -one, q, lda, one, r,
218  $ lda )
219 *
220 * Compute norm( I - Q'*Q ) / ( M * EPS ) .
221 *
222  resid = dlansy( '1', 'Upper', n, r, lda, rwork )
223 *
224  result( 2 ) = ( resid / dble( max( 1, m ) ) ) / eps
225 *
226  return
227 *
228 * End of DQRT02
229 *
230  END