LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
November 2011

Definition at line 178 of file cgrqts.f.

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, p, n
186 * ..
187 * .. Array Arguments ..
188  REAL result( 4 ), rwork( * )
189  COMPLEX a( lda, * ), af( lda, * ), r( lda, * ),
190  $ q( lda, * ), b( ldb, * ), bf( ldb, * ),
191  $ t( ldb, * ), z( ldb, * ), bwk( ldb, * ),
192  $ taua( * ), taub( * ), work( lwork )
193 * ..
194 *
195 * =====================================================================
196 *
197 * .. Parameters ..
198  REAL zero, one
199  parameter ( zero = 0.0e+0, one = 1.0e+0 )
200  COMPLEX czero, cone
201  parameter ( czero = ( 0.0e+0, 0.0e+0 ),
202  $ cone = ( 1.0e+0, 0.0e+0 ) )
203  COMPLEX crogue
204  parameter ( crogue = ( -1.0e+10, 0.0e+0 ) )
205 * ..
206 * .. Local Scalars ..
207  INTEGER info
208  REAL anorm, bnorm, ulp, unfl, resid
209 * ..
210 * .. External Functions ..
211  REAL slamch, clange, clanhe
212  EXTERNAL slamch, clange, clanhe
213 * ..
214 * .. External Subroutines ..
215  EXTERNAL cgemm, cggrqf, clacpy, claset, cungqr,
216  $ cungrq, cherk
217 * ..
218 * .. Intrinsic Functions ..
219  INTRINSIC max, min, real
220 * ..
221 * .. Executable Statements ..
222 *
223  ulp = slamch( 'Precision' )
224  unfl = slamch( 'Safe minimum' )
225 *
226 * Copy the matrix A to the array AF.
227 *
228  CALL clacpy( 'Full', m, n, a, lda, af, lda )
229  CALL clacpy( 'Full', p, n, b, ldb, bf, ldb )
230 *
231  anorm = max( clange( '1', m, n, a, lda, rwork ), unfl )
232  bnorm = max( clange( '1', p, n, b, ldb, rwork ), unfl )
233 *
234 * Factorize the matrices A and B in the arrays AF and BF.
235 *
236  CALL cggrqf( m, p, n, af, lda, taua, bf, ldb, taub, work,
237  $ lwork, info )
238 *
239 * Generate the N-by-N matrix Q
240 *
241  CALL claset( 'Full', n, n, crogue, crogue, q, lda )
242  IF( m.LE.n ) THEN
243  IF( m.GT.0 .AND. m.LT.n )
244  $ CALL clacpy( 'Full', m, n-m, af, lda, q( n-m+1, 1 ), lda )
245  IF( m.GT.1 )
246  $ CALL clacpy( 'Lower', m-1, m-1, af( 2, n-m+1 ), lda,
247  $ q( n-m+2, n-m+1 ), lda )
248  ELSE
249  IF( n.GT.1 )
250  $ CALL clacpy( 'Lower', n-1, n-1, af( m-n+2, 1 ), lda,
251  $ q( 2, 1 ), lda )
252  END IF
253  CALL cungrq( n, n, min( m, n ), q, lda, taua, work, lwork, info )
254 *
255 * Generate the P-by-P matrix Z
256 *
257  CALL claset( 'Full', p, p, crogue, crogue, z, ldb )
258  IF( p.GT.1 )
259  $ CALL clacpy( 'Lower', p-1, n, bf( 2,1 ), ldb, z( 2,1 ), ldb )
260  CALL cungqr( p, p, min( p,n ), z, ldb, taub, work, lwork, info )
261 *
262 * Copy R
263 *
264  CALL claset( 'Full', m, n, czero, czero, r, lda )
265  IF( m.LE.n )THEN
266  CALL clacpy( 'Upper', m, m, af( 1, n-m+1 ), lda, r( 1, n-m+1 ),
267  $ lda )
268  ELSE
269  CALL clacpy( 'Full', m-n, n, af, lda, r, lda )
270  CALL clacpy( 'Upper', n, n, af( m-n+1, 1 ), lda, r( m-n+1, 1 ),
271  $ lda )
272  END IF
273 *
274 * Copy T
275 *
276  CALL claset( 'Full', p, n, czero, czero, t, ldb )
277  CALL clacpy( 'Upper', p, n, bf, ldb, t, ldb )
278 *
279 * Compute R - A*Q'
280 *
281  CALL cgemm( 'No transpose', 'Conjugate transpose', m, n, n, -cone,
282  $ a, lda, q, lda, cone, r, lda )
283 *
284 * Compute norm( R - A*Q' ) / ( MAX(M,N)*norm(A)*ULP ) .
285 *
286  resid = clange( '1', m, n, r, lda, rwork )
287  IF( anorm.GT.zero ) THEN
288  result( 1 ) = ( ( resid / REAL(MAX(1,M,N) ) ) / anorm ) / ulp
289  ELSE
290  result( 1 ) = zero
291  END IF
292 *
293 * Compute T*Q - Z'*B
294 *
295  CALL cgemm( 'Conjugate transpose', 'No transpose', p, n, p, cone,
296  $ z, ldb, b, ldb, czero, bwk, ldb )
297  CALL cgemm( 'No transpose', 'No transpose', p, n, n, cone, t, ldb,
298  $ q, lda, -cone, bwk, ldb )
299 *
300 * Compute norm( T*Q - Z'*B ) / ( MAX(P,N)*norm(A)*ULP ) .
301 *
302  resid = clange( '1', p, n, bwk, ldb, rwork )
303  IF( bnorm.GT.zero ) THEN
304  result( 2 ) = ( ( resid / REAL( MAX( 1,P,M ) ) )/bnorm ) / ulp
305  ELSE
306  result( 2 ) = zero
307  END IF
308 *
309 * Compute I - Q*Q'
310 *
311  CALL claset( 'Full', n, n, czero, cone, r, lda )
312  CALL cherk( 'Upper', 'No Transpose', n, n, -one, q, lda, one, r,
313  $ lda )
314 *
315 * Compute norm( I - Q'*Q ) / ( N * ULP ) .
316 *
317  resid = clanhe( '1', 'Upper', n, r, lda, rwork )
318  result( 3 ) = ( resid / REAL( MAX( 1,N ) ) ) / ulp
319 *
320 * Compute I - Z'*Z
321 *
322  CALL claset( 'Full', p, p, czero, cone, t, ldb )
323  CALL cherk( 'Upper', 'Conjugate transpose', p, p, -one, z, ldb,
324  $ one, t, ldb )
325 *
326 * Compute norm( I - Z'*Z ) / ( P*ULP ) .
327 *
328  resid = clanhe( '1', 'Upper', p, t, ldb, rwork )
329  result( 4 ) = ( resid / REAL( MAX( 1,P ) ) ) / ulp
330 *
331  RETURN
332 *
333 * End of CGRQTS
334 *
subroutine cggrqf(M, P, N, A, LDA, TAUA, B, LDB, TAUB, WORK, LWORK, INFO)
CGGRQF
Definition: cggrqf.f:216
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, or the element of largest absolute value of a complex Hermitian matrix.
Definition: clanhe.f:126
subroutine cherk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
CHERK
Definition: cherk.f:175
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:117
subroutine cungrq(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
CUNGRQ
Definition: cungrq.f:130
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:108
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine cungqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
CUNGQR
Definition: cungqr.f:130
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:189

Here is the call graph for this function:

Here is the caller graph for this function: