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