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

◆ drqt01()

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

DRQT01

Purpose:
 DRQT01 tests DGERQF, which computes the RQ factorization of an m-by-n
 matrix A, and partially tests DORGRQ which forms the n-by-n
 orthogonal matrix Q.

 DRQT01 compares R with A*Q', and checks that Q is orthogonal.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The m-by-n matrix A.
[out]AF
          AF is DOUBLE PRECISION array, dimension (LDA,N)
          Details of the RQ factorization of A, as returned by DGERQF.
          See DGERQF for further details.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDA,N)
          The n-by-n orthogonal matrix Q.
[out]R
          R is DOUBLE PRECISION array, dimension (LDA,max(M,N))
[in]LDA
          LDA is INTEGER
          The leading dimension of the arrays A, AF, Q and L.
          LDA >= max(M,N).
[out]TAU
          TAU is DOUBLE PRECISION array, dimension (min(M,N))
          The scalar factors of the elementary reflectors, as returned
          by DGERQF.
[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 (max(M,N))
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (2)
          The test ratios:
          RESULT(1) = norm( R - 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 124 of file drqt01.f.

126*
127* -- LAPACK test routine --
128* -- LAPACK is a software package provided by Univ. of Tennessee, --
129* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
130*
131* .. Scalar Arguments ..
132 INTEGER LDA, LWORK, M, N
133* ..
134* .. Array Arguments ..
135 DOUBLE PRECISION A( LDA, * ), AF( LDA, * ), Q( LDA, * ),
136 $ R( LDA, * ), RESULT( * ), RWORK( * ), TAU( * ),
137 $ WORK( LWORK )
138* ..
139*
140* =====================================================================
141*
142* .. Parameters ..
143 DOUBLE PRECISION ZERO, ONE
144 parameter( zero = 0.0d+0, one = 1.0d+0 )
145 DOUBLE PRECISION ROGUE
146 parameter( rogue = -1.0d+10 )
147* ..
148* .. Local Scalars ..
149 INTEGER INFO, MINMN
150 DOUBLE PRECISION ANORM, EPS, RESID
151* ..
152* .. External Functions ..
153 DOUBLE PRECISION DLAMCH, DLANGE, DLANSY
154 EXTERNAL dlamch, dlange, dlansy
155* ..
156* .. External Subroutines ..
157 EXTERNAL dgemm, dgerqf, dlacpy, dlaset, dorgrq, dsyrk
158* ..
159* .. Intrinsic Functions ..
160 INTRINSIC dble, max, min
161* ..
162* .. Scalars in Common ..
163 CHARACTER*32 SRNAMT
164* ..
165* .. Common blocks ..
166 COMMON / srnamc / srnamt
167* ..
168* .. Executable Statements ..
169*
170 minmn = min( m, n )
171 eps = dlamch( 'Epsilon' )
172*
173* Copy the matrix A to the array AF.
174*
175 CALL dlacpy( 'Full', m, n, a, lda, af, lda )
176*
177* Factorize the matrix A in the array AF.
178*
179 srnamt = 'DGERQF'
180 CALL dgerqf( m, n, af, lda, tau, work, lwork, info )
181*
182* Copy details of Q
183*
184 CALL dlaset( 'Full', n, n, rogue, rogue, q, lda )
185 IF( m.LE.n ) THEN
186 IF( m.GT.0 .AND. m.LT.n )
187 $ CALL dlacpy( 'Full', m, n-m, af, lda, q( n-m+1, 1 ), lda )
188 IF( m.GT.1 )
189 $ CALL dlacpy( 'Lower', m-1, m-1, af( 2, n-m+1 ), lda,
190 $ q( n-m+2, n-m+1 ), lda )
191 ELSE
192 IF( n.GT.1 )
193 $ CALL dlacpy( 'Lower', n-1, n-1, af( m-n+2, 1 ), lda,
194 $ q( 2, 1 ), lda )
195 END IF
196*
197* Generate the n-by-n matrix Q
198*
199 srnamt = 'DORGRQ'
200 CALL dorgrq( n, n, minmn, q, lda, tau, work, lwork, info )
201*
202* Copy R
203*
204 CALL dlaset( 'Full', m, n, zero, zero, r, lda )
205 IF( m.LE.n ) THEN
206 IF( m.GT.0 )
207 $ CALL dlacpy( 'Upper', m, m, af( 1, n-m+1 ), lda,
208 $ r( 1, n-m+1 ), lda )
209 ELSE
210 IF( m.GT.n .AND. n.GT.0 )
211 $ CALL dlacpy( 'Full', m-n, n, af, lda, r, lda )
212 IF( n.GT.0 )
213 $ CALL dlacpy( 'Upper', n, n, af( m-n+1, 1 ), lda,
214 $ r( m-n+1, 1 ), lda )
215 END IF
216*
217* Compute R - A*Q'
218*
219 CALL dgemm( 'No transpose', 'Transpose', m, n, n, -one, a, lda, q,
220 $ lda, one, r, lda )
221*
222* Compute norm( R - Q'*A ) / ( N * norm(A) * EPS ) .
223*
224 anorm = dlange( '1', m, n, a, lda, rwork )
225 resid = dlange( '1', m, n, r, lda, rwork )
226 IF( anorm.GT.zero ) THEN
227 result( 1 ) = ( ( resid / dble( max( 1, n ) ) ) / anorm ) / eps
228 ELSE
229 result( 1 ) = zero
230 END IF
231*
232* Compute I - Q*Q'
233*
234 CALL dlaset( 'Full', n, n, zero, one, r, lda )
235 CALL dsyrk( 'Upper', 'No transpose', n, n, -one, q, lda, one, r,
236 $ lda )
237*
238* Compute norm( I - Q*Q' ) / ( N * EPS ) .
239*
240 resid = dlansy( '1', 'Upper', n, r, lda, rwork )
241*
242 result( 2 ) = ( resid / dble( max( 1, n ) ) ) / eps
243*
244 RETURN
245*
246* End of DRQT01
247*
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dgerqf(m, n, a, lda, tau, work, lwork, info)
DGERQF
Definition dgerqf.f:139
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 dorgrq(m, n, k, a, lda, tau, work, lwork, info)
DORGRQ
Definition dorgrq.f:128
Here is the call graph for this function:
Here is the caller graph for this function: