126 SUBROUTINE clqt01( M, N, A, AF, Q, L, LDA, TAU, WORK, LWORK,
135 INTEGER LDA, LWORK, M, N
138 REAL RESULT( * ), RWORK( * )
139 COMPLEX A( lda, * ), AF( lda, * ), L( lda, * ),
140 $ q( lda, * ), tau( * ), work( lwork )
147 parameter ( zero = 0.0e+0, one = 1.0e+0 )
149 parameter ( rogue = ( -1.0e+10, -1.0e+10 ) )
153 REAL ANORM, EPS, RESID
156 REAL CLANGE, CLANSY, SLAMCH
157 EXTERNAL clange, clansy, slamch
163 INTRINSIC cmplx, max, min, real
169 COMMON / srnamc / srnamt
174 eps = slamch(
'Epsilon' )
178 CALL clacpy(
'Full', m, n, a, lda, af, lda )
183 CALL cgelqf( m, n, af, lda, tau, work, lwork, info )
187 CALL claset(
'Full', n, n, rogue, rogue, q, lda )
189 $
CALL clacpy(
'Upper', m, n-1, af( 1, 2 ), lda, q( 1, 2 ), lda )
194 CALL cunglq( n, n, minmn, q, lda, tau, work, lwork, info )
198 CALL claset(
'Full', m, n, cmplx( zero ), cmplx( zero ), l, lda )
199 CALL clacpy(
'Lower', m, n, af, lda, l, lda )
203 CALL cgemm(
'No transpose',
'Conjugate transpose', m, n, n,
204 $ cmplx( -one ), a, lda, q, lda, cmplx( one ), l, lda )
208 anorm = clange(
'1', m, n, a, lda, rwork )
209 resid = clange(
'1', m, n, l, lda, rwork )
210 IF( anorm.GT.zero )
THEN
211 result( 1 ) = ( ( resid /
REAL( MAX( 1, N ) ) ) / anorm ) / eps
218 CALL claset(
'Full', n, n, cmplx( zero ), cmplx( one ), l, lda )
219 CALL cherk(
'Upper',
'No transpose', n, n, -one, q, lda, one, l,
224 resid = clansy(
'1',
'Upper', n, l, lda, rwork )
226 result( 2 ) = ( resid /
REAL( MAX( 1, N ) ) ) / eps
subroutine clqt01(M, N, A, AF, Q, L, LDA, TAU, WORK, LWORK, RWORK, RESULT)
CLQT01
subroutine cherk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
CHERK
subroutine cgelqf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
CGELQF
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...
subroutine cunglq(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
CUNGLQ
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM