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

◆ dlqt02()

subroutine dlqt02 ( integer  m,
integer  n,
integer  k,
double precision, dimension( lda, * )  a,
double precision, dimension( lda, * )  af,
double precision, dimension( lda, * )  q,
double precision, dimension( lda, * )  l,
integer  lda,
double precision, dimension( * )  tau,
double precision, dimension( lwork )  work,
integer  lwork,
double precision, dimension( * )  rwork,
double precision, dimension( * )  result 
)

DLQT02

Purpose:
 DLQT02 tests DORGLQ, which generates an m-by-n matrix Q with
 orthonormal rows that is defined as the product of k elementary
 reflectors.

 Given the LQ factorization of an m-by-n matrix A, DLQT02 generates
 the orthogonal matrix Q defined by the factorization of the first k
 rows of A; it compares L(1:k,1:m) with A(1:k,1:n)*Q(1:m,1:n)', and
 checks that the rows of Q are orthonormal.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix Q to be generated.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix Q to be generated.
          N >= M >= 0.
[in]K
          K is INTEGER
          The number of elementary reflectors whose product defines the
          matrix Q. M >= K >= 0.
[in]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The m-by-n matrix A which was factorized by DLQT01.
[in]AF
          AF is DOUBLE PRECISION array, dimension (LDA,N)
          Details of the LQ factorization of A, as returned by DGELQF.
          See DGELQF for further details.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDA,N)
[out]L
          L is DOUBLE PRECISION array, dimension (LDA,M)
[in]LDA
          LDA is INTEGER
          The leading dimension of the arrays A, AF, Q and L. LDA >= N.
[in]TAU
          TAU is DOUBLE PRECISION array, dimension (M)
          The scalar factors of the elementary reflectors corresponding
          to the LQ factorization in AF.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (M)
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (2)
          The test ratios:
          RESULT(1) = norm( L - A*Q' ) / ( N * norm(A) * EPS )
          RESULT(2) = norm( I - Q*Q' ) / ( N * EPS )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 133 of file dlqt02.f.

135*
136* -- LAPACK test routine --
137* -- LAPACK is a software package provided by Univ. of Tennessee, --
138* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
139*
140* .. Scalar Arguments ..
141 INTEGER K, LDA, LWORK, M, N
142* ..
143* .. Array Arguments ..
144 DOUBLE PRECISION A( LDA, * ), AF( LDA, * ), L( LDA, * ),
145 $ Q( LDA, * ), RESULT( * ), RWORK( * ), TAU( * ),
146 $ WORK( LWORK )
147* ..
148*
149* =====================================================================
150*
151* .. Parameters ..
152 DOUBLE PRECISION ZERO, ONE
153 parameter( zero = 0.0d+0, one = 1.0d+0 )
154 DOUBLE PRECISION ROGUE
155 parameter( rogue = -1.0d+10 )
156* ..
157* .. Local Scalars ..
158 INTEGER INFO
159 DOUBLE PRECISION ANORM, EPS, RESID
160* ..
161* .. External Functions ..
162 DOUBLE PRECISION DLAMCH, DLANGE, DLANSY
163 EXTERNAL dlamch, dlange, dlansy
164* ..
165* .. External Subroutines ..
166 EXTERNAL dgemm, dlacpy, dlaset, dorglq, dsyrk
167* ..
168* .. Intrinsic Functions ..
169 INTRINSIC dble, max
170* ..
171* .. Scalars in Common ..
172 CHARACTER*32 SRNAMT
173* ..
174* .. Common blocks ..
175 COMMON / srnamc / srnamt
176* ..
177* .. Executable Statements ..
178*
179 eps = dlamch( 'Epsilon' )
180*
181* Copy the first k rows of the factorization to the array Q
182*
183 CALL dlaset( 'Full', m, n, rogue, rogue, q, lda )
184 CALL dlacpy( 'Upper', k, n-1, af( 1, 2 ), lda, q( 1, 2 ), lda )
185*
186* Generate the first n columns of the matrix Q
187*
188 srnamt = 'DORGLQ'
189 CALL dorglq( m, n, k, q, lda, tau, work, lwork, info )
190*
191* Copy L(1:k,1:m)
192*
193 CALL dlaset( 'Full', k, m, zero, zero, l, lda )
194 CALL dlacpy( 'Lower', k, m, af, lda, l, lda )
195*
196* Compute L(1:k,1:m) - A(1:k,1:n) * Q(1:m,1:n)'
197*
198 CALL dgemm( 'No transpose', 'Transpose', k, m, n, -one, a, lda, q,
199 $ lda, one, l, lda )
200*
201* Compute norm( L - A*Q' ) / ( N * norm(A) * EPS ) .
202*
203 anorm = dlange( '1', k, n, a, lda, rwork )
204 resid = dlange( '1', k, m, l, lda, rwork )
205 IF( anorm.GT.zero ) THEN
206 result( 1 ) = ( ( resid / dble( max( 1, n ) ) ) / anorm ) / eps
207 ELSE
208 result( 1 ) = zero
209 END IF
210*
211* Compute I - Q*Q'
212*
213 CALL dlaset( 'Full', m, m, zero, one, l, lda )
214 CALL dsyrk( 'Upper', 'No transpose', m, n, -one, q, lda, one, l,
215 $ lda )
216*
217* Compute norm( I - Q*Q' ) / ( N * EPS ) .
218*
219 resid = dlansy( '1', 'Upper', m, l, lda, rwork )
220*
221 result( 2 ) = ( resid / dble( max( 1, n ) ) ) / eps
222*
223 RETURN
224*
225* End of DLQT02
226*
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
DSYRK
Definition dsyrk.f:169
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function dlange(norm, m, n, a, lda, work)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition dlange.f:114
double precision function dlansy(norm, uplo, n, a, lda, work)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition dlansy.f:122
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:110
subroutine dorglq(m, n, k, a, lda, tau, work, lwork, info)
DORGLQ
Definition dorglq.f:127
Here is the call graph for this function:
Here is the caller graph for this function: