LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dgrqts.f
Go to the documentation of this file.
1 *> \brief \b DGRQTS
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 DGRQTS( M, P, N, A, AF, Q, R, LDA, TAUA, B, BF, Z, T,
12 * BWK, LDB, TAUB, WORK, LWORK, RWORK, RESULT )
13 *
14 * .. Scalar Arguments ..
15 * INTEGER LDA, LDB, LWORK, M, N, P
16 * ..
17 * .. Array Arguments ..
18 * DOUBLE PRECISION A( LDA, * ), AF( LDA, * ), B( LDB, * ),
19 * $ BF( LDB, * ), BWK( LDB, * ), Q( LDA, * ),
20 * $ R( LDA, * ), RESULT( 4 ), RWORK( * ),
21 * $ T( LDB, * ), TAUA( * ), TAUB( * ),
22 * $ WORK( LWORK ), Z( LDB, * )
23 * ..
24 *
25 *
26 *> \par Purpose:
27 * =============
28 *>
29 *> \verbatim
30 *>
31 *> DGRQTS tests DGGRQF, which computes the GRQ factorization of an
32 *> M-by-N matrix A and a P-by-N matrix B: A = R*Q and B = Z*T*Q.
33 *> \endverbatim
34 *
35 * Arguments:
36 * ==========
37 *
38 *> \param[in] M
39 *> \verbatim
40 *> M is INTEGER
41 *> The number of rows of the matrix A. M >= 0.
42 *> \endverbatim
43 *>
44 *> \param[in] P
45 *> \verbatim
46 *> P is INTEGER
47 *> The number of rows of the matrix B. P >= 0.
48 *> \endverbatim
49 *>
50 *> \param[in] N
51 *> \verbatim
52 *> N is INTEGER
53 *> The number of columns of the matrices A and B. N >= 0.
54 *> \endverbatim
55 *>
56 *> \param[in] A
57 *> \verbatim
58 *> A is DOUBLE PRECISION array, dimension (LDA,N)
59 *> The M-by-N matrix A.
60 *> \endverbatim
61 *>
62 *> \param[out] AF
63 *> \verbatim
64 *> AF is DOUBLE PRECISION array, dimension (LDA,N)
65 *> Details of the GRQ factorization of A and B, as returned
66 *> by DGGRQF, see SGGRQF for further details.
67 *> \endverbatim
68 *>
69 *> \param[out] Q
70 *> \verbatim
71 *> Q is DOUBLE PRECISION array, dimension (LDA,N)
72 *> The N-by-N orthogonal matrix Q.
73 *> \endverbatim
74 *>
75 *> \param[out] R
76 *> \verbatim
77 *> R is DOUBLE PRECISION array, dimension (LDA,MAX(M,N))
78 *> \endverbatim
79 *>
80 *> \param[in] LDA
81 *> \verbatim
82 *> LDA is INTEGER
83 *> The leading dimension of the arrays A, AF, R and Q.
84 *> LDA >= max(M,N).
85 *> \endverbatim
86 *>
87 *> \param[out] TAUA
88 *> \verbatim
89 *> TAUA is DOUBLE PRECISION array, dimension (min(M,N))
90 *> The scalar factors of the elementary reflectors, as returned
91 *> by DGGQRC.
92 *> \endverbatim
93 *>
94 *> \param[in] B
95 *> \verbatim
96 *> B is DOUBLE PRECISION array, dimension (LDB,N)
97 *> On entry, the P-by-N matrix A.
98 *> \endverbatim
99 *>
100 *> \param[out] BF
101 *> \verbatim
102 *> BF is DOUBLE PRECISION array, dimension (LDB,N)
103 *> Details of the GQR factorization of A and B, as returned
104 *> by DGGRQF, see SGGRQF for further details.
105 *> \endverbatim
106 *>
107 *> \param[out] Z
108 *> \verbatim
109 *> Z is DOUBLE PRECISION array, dimension (LDB,P)
110 *> The P-by-P orthogonal matrix Z.
111 *> \endverbatim
112 *>
113 *> \param[out] T
114 *> \verbatim
115 *> T is DOUBLE PRECISION array, dimension (LDB,max(P,N))
116 *> \endverbatim
117 *>
118 *> \param[out] BWK
119 *> \verbatim
120 *> BWK is DOUBLE PRECISION array, dimension (LDB,N)
121 *> \endverbatim
122 *>
123 *> \param[in] LDB
124 *> \verbatim
125 *> LDB is INTEGER
126 *> The leading dimension of the arrays B, BF, Z and T.
127 *> LDB >= max(P,N).
128 *> \endverbatim
129 *>
130 *> \param[out] TAUB
131 *> \verbatim
132 *> TAUB is DOUBLE PRECISION array, dimension (min(P,N))
133 *> The scalar factors of the elementary reflectors, as returned
134 *> by DGGRQF.
135 *> \endverbatim
136 *>
137 *> \param[out] WORK
138 *> \verbatim
139 *> WORK is DOUBLE PRECISION array, dimension (LWORK)
140 *> \endverbatim
141 *>
142 *> \param[in] LWORK
143 *> \verbatim
144 *> LWORK is INTEGER
145 *> The dimension of the array WORK, LWORK >= max(M,P,N)**2.
146 *> \endverbatim
147 *>
148 *> \param[out] RWORK
149 *> \verbatim
150 *> RWORK is DOUBLE PRECISION array, dimension (M)
151 *> \endverbatim
152 *>
153 *> \param[out] RESULT
154 *> \verbatim
155 *> RESULT is DOUBLE PRECISION array, dimension (4)
156 *> The test ratios:
157 *> RESULT(1) = norm( R - A*Q' ) / ( MAX(M,N)*norm(A)*ULP)
158 *> RESULT(2) = norm( T*Q - Z'*B ) / (MAX(P,N)*norm(B)*ULP)
159 *> RESULT(3) = norm( I - Q'*Q ) / ( N*ULP )
160 *> RESULT(4) = norm( I - Z'*Z ) / ( P*ULP )
161 *> \endverbatim
162 *
163 * Authors:
164 * ========
165 *
166 *> \author Univ. of Tennessee
167 *> \author Univ. of California Berkeley
168 *> \author Univ. of Colorado Denver
169 *> \author NAG Ltd.
170 *
171 *> \date November 2011
172 *
173 *> \ingroup double_eig
174 *
175 * =====================================================================
176  SUBROUTINE dgrqts( M, P, N, A, AF, Q, R, LDA, TAUA, B, BF, Z, T,
177  $ bwk, ldb, taub, work, lwork, rwork, result )
178 *
179 * -- LAPACK test routine (version 3.4.0) --
180 * -- LAPACK is a software package provided by Univ. of Tennessee, --
181 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
182 * November 2011
183 *
184 * .. Scalar Arguments ..
185  INTEGER LDA, LDB, LWORK, M, N, P
186 * ..
187 * .. Array Arguments ..
188  DOUBLE PRECISION A( lda, * ), AF( lda, * ), B( ldb, * ),
189  $ bf( ldb, * ), bwk( ldb, * ), q( lda, * ),
190  $ r( lda, * ), result( 4 ), rwork( * ),
191  $ t( ldb, * ), taua( * ), taub( * ),
192  $ work( lwork ), z( ldb, * )
193 * ..
194 *
195 * =====================================================================
196 *
197 * .. Parameters ..
198  DOUBLE PRECISION ZERO, ONE
199  parameter ( zero = 0.0d+0, one = 1.0d+0 )
200  DOUBLE PRECISION ROGUE
201  parameter ( rogue = -1.0d+10 )
202 * ..
203 * .. Local Scalars ..
204  INTEGER INFO
205  DOUBLE PRECISION ANORM, BNORM, RESID, ULP, UNFL
206 * ..
207 * .. External Functions ..
208  DOUBLE PRECISION DLAMCH, DLANGE, DLANSY
209  EXTERNAL dlamch, dlange, dlansy
210 * ..
211 * .. External Subroutines ..
212  EXTERNAL dgemm, dggrqf, dlacpy, dlaset, dorgqr, dorgrq,
213  $ dsyrk
214 * ..
215 * .. Intrinsic Functions ..
216  INTRINSIC dble, max, min
217 * ..
218 * .. Executable Statements ..
219 *
220  ulp = dlamch( 'Precision' )
221  unfl = dlamch( 'Safe minimum' )
222 *
223 * Copy the matrix A to the array AF.
224 *
225  CALL dlacpy( 'Full', m, n, a, lda, af, lda )
226  CALL dlacpy( 'Full', p, n, b, ldb, bf, ldb )
227 *
228  anorm = max( dlange( '1', m, n, a, lda, rwork ), unfl )
229  bnorm = max( dlange( '1', p, n, b, ldb, rwork ), unfl )
230 *
231 * Factorize the matrices A and B in the arrays AF and BF.
232 *
233  CALL dggrqf( m, p, n, af, lda, taua, bf, ldb, taub, work, lwork,
234  $ info )
235 *
236 * Generate the N-by-N matrix Q
237 *
238  CALL dlaset( 'Full', n, n, rogue, rogue, q, lda )
239  IF( m.LE.n ) THEN
240  IF( m.GT.0 .AND. m.LT.n )
241  $ CALL dlacpy( 'Full', m, n-m, af, lda, q( n-m+1, 1 ), lda )
242  IF( m.GT.1 )
243  $ CALL dlacpy( 'Lower', m-1, m-1, af( 2, n-m+1 ), lda,
244  $ q( n-m+2, n-m+1 ), lda )
245  ELSE
246  IF( n.GT.1 )
247  $ CALL dlacpy( 'Lower', n-1, n-1, af( m-n+2, 1 ), lda,
248  $ q( 2, 1 ), lda )
249  END IF
250  CALL dorgrq( n, n, min( m, n ), q, lda, taua, work, lwork, info )
251 *
252 * Generate the P-by-P matrix Z
253 *
254  CALL dlaset( 'Full', p, p, rogue, rogue, z, ldb )
255  IF( p.GT.1 )
256  $ CALL dlacpy( 'Lower', p-1, n, bf( 2, 1 ), ldb, z( 2, 1 ), ldb )
257  CALL dorgqr( p, p, min( p, n ), z, ldb, taub, work, lwork, info )
258 *
259 * Copy R
260 *
261  CALL dlaset( 'Full', m, n, zero, zero, r, lda )
262  IF( m.LE.n ) THEN
263  CALL dlacpy( 'Upper', m, m, af( 1, n-m+1 ), lda, r( 1, n-m+1 ),
264  $ lda )
265  ELSE
266  CALL dlacpy( 'Full', m-n, n, af, lda, r, lda )
267  CALL dlacpy( 'Upper', n, n, af( m-n+1, 1 ), lda, r( m-n+1, 1 ),
268  $ lda )
269  END IF
270 *
271 * Copy T
272 *
273  CALL dlaset( 'Full', p, n, zero, zero, t, ldb )
274  CALL dlacpy( 'Upper', p, n, bf, ldb, t, ldb )
275 *
276 * Compute R - A*Q'
277 *
278  CALL dgemm( 'No transpose', 'Transpose', m, n, n, -one, a, lda, q,
279  $ lda, one, r, lda )
280 *
281 * Compute norm( R - A*Q' ) / ( MAX(M,N)*norm(A)*ULP ) .
282 *
283  resid = dlange( '1', m, n, r, lda, rwork )
284  IF( anorm.GT.zero ) THEN
285  result( 1 ) = ( ( resid / dble( max( 1, m, n ) ) ) / anorm ) /
286  $ ulp
287  ELSE
288  result( 1 ) = zero
289  END IF
290 *
291 * Compute T*Q - Z'*B
292 *
293  CALL dgemm( 'Transpose', 'No transpose', p, n, p, one, z, ldb, b,
294  $ ldb, zero, bwk, ldb )
295  CALL dgemm( 'No transpose', 'No transpose', p, n, n, one, t, ldb,
296  $ q, lda, -one, bwk, ldb )
297 *
298 * Compute norm( T*Q - Z'*B ) / ( MAX(P,N)*norm(A)*ULP ) .
299 *
300  resid = dlange( '1', p, n, bwk, ldb, rwork )
301  IF( bnorm.GT.zero ) THEN
302  result( 2 ) = ( ( resid / dble( max( 1, p, m ) ) ) / bnorm ) /
303  $ ulp
304  ELSE
305  result( 2 ) = zero
306  END IF
307 *
308 * Compute I - Q*Q'
309 *
310  CALL dlaset( 'Full', n, n, zero, one, r, lda )
311  CALL dsyrk( 'Upper', 'No Transpose', n, n, -one, q, lda, one, r,
312  $ lda )
313 *
314 * Compute norm( I - Q'*Q ) / ( N * ULP ) .
315 *
316  resid = dlansy( '1', 'Upper', n, r, lda, rwork )
317  result( 3 ) = ( resid / dble( max( 1, n ) ) ) / ulp
318 *
319 * Compute I - Z'*Z
320 *
321  CALL dlaset( 'Full', p, p, zero, one, t, ldb )
322  CALL dsyrk( 'Upper', 'Transpose', p, p, -one, z, ldb, one, t,
323  $ ldb )
324 *
325 * Compute norm( I - Z'*Z ) / ( P*ULP ) .
326 *
327  resid = dlansy( '1', 'Upper', p, t, ldb, rwork )
328  result( 4 ) = ( resid / dble( max( 1, p ) ) ) / ulp
329 *
330  RETURN
331 *
332 * End of DGRQTS
333 *
334  END
subroutine dlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: dlaset.f:112
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine dsyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
DSYRK
Definition: dsyrk.f:171
subroutine dorgrq(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
DORGRQ
Definition: dorgrq.f:130
subroutine dgrqts(M, P, N, A, AF, Q, R, LDA, TAUA, B, BF, Z, T, BWK, LDB, TAUB, WORK, LWORK, RWORK, RESULT)
DGRQTS
Definition: dgrqts.f:178
subroutine dggrqf(M, P, N, A, LDA, TAUA, B, LDB, TAUB, WORK, LWORK, INFO)
DGGRQF
Definition: dggrqf.f:216
subroutine dorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
DORGQR
Definition: dorgqr.f:130