LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ cgrqts()

subroutine cgrqts ( integer  m,
integer  p,
integer  n,
complex, dimension( lda, * )  a,
complex, dimension( lda, * )  af,
complex, dimension( lda, * )  q,
complex, dimension( lda, * )  r,
integer  lda,
complex, dimension( * )  taua,
complex, dimension( ldb, * )  b,
complex, dimension( ldb, * )  bf,
complex, dimension( ldb, * )  z,
complex, dimension( ldb, * )  t,
complex, dimension( ldb, * )  bwk,
integer  ldb,
complex, dimension( * )  taub,
complex, dimension( lwork )  work,
integer  lwork,
real, dimension( * )  rwork,
real, dimension( 4 )  result 
)

CGRQTS

Purpose:
 CGRQTS tests CGGRQF, which computes the GRQ factorization of an
 M-by-N matrix A and a P-by-N matrix B: A = R*Q and B = Z*T*Q.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]P
          P is INTEGER
          The number of rows of the matrix B.  P >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrices A and B.  N >= 0.
[in]A
          A is COMPLEX array, dimension (LDA,N)
          The M-by-N matrix A.
[out]AF
          AF is COMPLEX array, dimension (LDA,N)
          Details of the GRQ factorization of A and B, as returned
          by CGGRQF, see CGGRQF for further details.
[out]Q
          Q is COMPLEX array, dimension (LDA,N)
          The N-by-N unitary matrix Q.
[out]R
          R is COMPLEX array, dimension (LDA,MAX(M,N))
[in]LDA
          LDA is INTEGER
          The leading dimension of the arrays A, AF, R and Q.
          LDA >= max(M,N).
[out]TAUA
          TAUA is COMPLEX array, dimension (min(M,N))
          The scalar factors of the elementary reflectors, as returned
          by SGGQRC.
[in]B
          B is COMPLEX array, dimension (LDB,N)
          On entry, the P-by-N matrix A.
[out]BF
          BF is COMPLEX array, dimension (LDB,N)
          Details of the GQR factorization of A and B, as returned
          by CGGRQF, see CGGRQF for further details.
[out]Z
          Z is REAL array, dimension (LDB,P)
          The P-by-P unitary matrix Z.
[out]T
          T is COMPLEX array, dimension (LDB,max(P,N))
[out]BWK
          BWK is COMPLEX array, dimension (LDB,N)
[in]LDB
          LDB is INTEGER
          The leading dimension of the arrays B, BF, Z and T.
          LDB >= max(P,N).
[out]TAUB
          TAUB is COMPLEX array, dimension (min(P,N))
          The scalar factors of the elementary reflectors, as returned
          by SGGRQF.
[out]WORK
          WORK is COMPLEX array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK, LWORK >= max(M,P,N)**2.
[out]RWORK
          RWORK is REAL array, dimension (M)
[out]RESULT
          RESULT is REAL array, dimension (4)
          The test ratios:
            RESULT(1) = norm( R - A*Q' ) / ( MAX(M,N)*norm(A)*ULP)
            RESULT(2) = norm( T*Q - Z'*B ) / (MAX(P,N)*norm(B)*ULP)
            RESULT(3) = norm( I - Q'*Q ) / ( N*ULP )
            RESULT(4) = norm( I - Z'*Z ) / ( P*ULP )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 174 of file cgrqts.f.

176*
177* -- LAPACK test routine --
178* -- LAPACK is a software package provided by Univ. of Tennessee, --
179* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
180*
181* .. Scalar Arguments ..
182 INTEGER LDA, LDB, LWORK, M, P, N
183* ..
184* .. Array Arguments ..
185 REAL RESULT( 4 ), RWORK( * )
186 COMPLEX A( LDA, * ), AF( LDA, * ), R( LDA, * ),
187 $ Q( LDA, * ), B( LDB, * ), BF( LDB, * ),
188 $ T( LDB, * ), Z( LDB, * ), BWK( LDB, * ),
189 $ TAUA( * ), TAUB( * ), WORK( LWORK )
190* ..
191*
192* =====================================================================
193*
194* .. Parameters ..
195 REAL ZERO, ONE
196 parameter( zero = 0.0e+0, one = 1.0e+0 )
197 COMPLEX CZERO, CONE
198 parameter( czero = ( 0.0e+0, 0.0e+0 ),
199 $ cone = ( 1.0e+0, 0.0e+0 ) )
200 COMPLEX CROGUE
201 parameter( crogue = ( -1.0e+10, 0.0e+0 ) )
202* ..
203* .. Local Scalars ..
204 INTEGER INFO
205 REAL ANORM, BNORM, ULP, UNFL, RESID
206* ..
207* .. External Functions ..
208 REAL SLAMCH, CLANGE, CLANHE
209 EXTERNAL slamch, clange, clanhe
210* ..
211* .. External Subroutines ..
212 EXTERNAL cgemm, cggrqf, clacpy, claset, cungqr,
213 $ cungrq, cherk
214* ..
215* .. Intrinsic Functions ..
216 INTRINSIC max, min, real
217* ..
218* .. Executable Statements ..
219*
220 ulp = slamch( 'Precision' )
221 unfl = slamch( 'Safe minimum' )
222*
223* Copy the matrix A to the array AF.
224*
225 CALL clacpy( 'Full', m, n, a, lda, af, lda )
226 CALL clacpy( 'Full', p, n, b, ldb, bf, ldb )
227*
228 anorm = max( clange( '1', m, n, a, lda, rwork ), unfl )
229 bnorm = max( clange( '1', p, n, b, ldb, rwork ), unfl )
230*
231* Factorize the matrices A and B in the arrays AF and BF.
232*
233 CALL cggrqf( m, p, n, af, lda, taua, bf, ldb, taub, work,
234 $ lwork, info )
235*
236* Generate the N-by-N matrix Q
237*
238 CALL claset( 'Full', n, n, crogue, crogue, q, lda )
239 IF( m.LE.n ) THEN
240 IF( m.GT.0 .AND. m.LT.n )
241 $ CALL clacpy( 'Full', m, n-m, af, lda, q( n-m+1, 1 ), lda )
242 IF( m.GT.1 )
243 $ CALL clacpy( '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 clacpy( 'Lower', n-1, n-1, af( m-n+2, 1 ), lda,
248 $ q( 2, 1 ), lda )
249 END IF
250 CALL cungrq( n, n, min( m, n ), q, lda, taua, work, lwork, info )
251*
252* Generate the P-by-P matrix Z
253*
254 CALL claset( 'Full', p, p, crogue, crogue, z, ldb )
255 IF( p.GT.1 )
256 $ CALL clacpy( 'Lower', p-1, n, bf( 2,1 ), ldb, z( 2,1 ), ldb )
257 CALL cungqr( p, p, min( p,n ), z, ldb, taub, work, lwork, info )
258*
259* Copy R
260*
261 CALL claset( 'Full', m, n, czero, czero, r, lda )
262 IF( m.LE.n )THEN
263 CALL clacpy( 'Upper', m, m, af( 1, n-m+1 ), lda, r( 1, n-m+1 ),
264 $ lda )
265 ELSE
266 CALL clacpy( 'Full', m-n, n, af, lda, r, lda )
267 CALL clacpy( '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 claset( 'Full', p, n, czero, czero, t, ldb )
274 CALL clacpy( 'Upper', p, n, bf, ldb, t, ldb )
275*
276* Compute R - A*Q'
277*
278 CALL cgemm( 'No transpose', 'Conjugate transpose', m, n, n, -cone,
279 $ a, lda, q, lda, cone, r, lda )
280*
281* Compute norm( R - A*Q' ) / ( MAX(M,N)*norm(A)*ULP ) .
282*
283 resid = clange( '1', m, n, r, lda, rwork )
284 IF( anorm.GT.zero ) THEN
285 result( 1 ) = ( ( resid / real(max(1,m,n) ) ) / anorm ) / ulp
286 ELSE
287 result( 1 ) = zero
288 END IF
289*
290* Compute T*Q - Z'*B
291*
292 CALL cgemm( 'Conjugate transpose', 'No transpose', p, n, p, cone,
293 $ z, ldb, b, ldb, czero, bwk, ldb )
294 CALL cgemm( 'No transpose', 'No transpose', p, n, n, cone, t, ldb,
295 $ q, lda, -cone, bwk, ldb )
296*
297* Compute norm( T*Q - Z'*B ) / ( MAX(P,N)*norm(A)*ULP ) .
298*
299 resid = clange( '1', p, n, bwk, ldb, rwork )
300 IF( bnorm.GT.zero ) THEN
301 result( 2 ) = ( ( resid / real( max( 1,p,m ) ) )/bnorm ) / ulp
302 ELSE
303 result( 2 ) = zero
304 END IF
305*
306* Compute I - Q*Q'
307*
308 CALL claset( 'Full', n, n, czero, cone, r, lda )
309 CALL cherk( 'Upper', 'No Transpose', n, n, -one, q, lda, one, r,
310 $ lda )
311*
312* Compute norm( I - Q'*Q ) / ( N * ULP ) .
313*
314 resid = clanhe( '1', 'Upper', n, r, lda, rwork )
315 result( 3 ) = ( resid / real( max( 1,n ) ) ) / ulp
316*
317* Compute I - Z'*Z
318*
319 CALL claset( 'Full', p, p, czero, cone, t, ldb )
320 CALL cherk( 'Upper', 'Conjugate transpose', p, p, -one, z, ldb,
321 $ one, t, ldb )
322*
323* Compute norm( I - Z'*Z ) / ( P*ULP ) .
324*
325 resid = clanhe( '1', 'Upper', p, t, ldb, rwork )
326 result( 4 ) = ( resid / real( max( 1,p ) ) ) / ulp
327*
328 RETURN
329*
330* End of CGRQTS
331*
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:188
subroutine cggrqf(m, p, n, a, lda, taua, b, ldb, taub, work, lwork, info)
CGGRQF
Definition cggrqf.f:214
subroutine cherk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
CHERK
Definition cherk.f:173
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clange(norm, m, n, a, lda, work)
CLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition clange.f:115
real function clanhe(norm, uplo, n, a, lda, work)
CLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clanhe.f:124
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:106
subroutine cungqr(m, n, k, a, lda, tau, work, lwork, info)
CUNGQR
Definition cungqr.f:128
subroutine cungrq(m, n, k, a, lda, tau, work, lwork, info)
CUNGRQ
Definition cungrq.f:128
Here is the call graph for this function:
Here is the caller graph for this function: