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

◆ dgrqts()

subroutine dgrqts ( integer  m,
integer  p,
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( * )  taua,
double precision, dimension( ldb, * )  b,
double precision, dimension( ldb, * )  bf,
double precision, dimension( ldb, * )  z,
double precision, dimension( ldb, * )  t,
double precision, dimension( ldb, * )  bwk,
integer  ldb,
double precision, dimension( * )  taub,
double precision, dimension( lwork )  work,
integer  lwork,
double precision, dimension( * )  rwork,
double precision, dimension( 4 )  result 
)

DGRQTS

Purpose:
 DGRQTS tests DGGRQF, 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 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 GRQ factorization of A and B, as returned
          by DGGRQF, see SGGRQF 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, R and Q.
          LDA >= max(M,N).
[out]TAUA
          TAUA is DOUBLE PRECISION array, dimension (min(M,N))
          The scalar factors of the elementary reflectors, as returned
          by DGGQRC.
[in]B
          B is DOUBLE PRECISION array, dimension (LDB,N)
          On entry, the P-by-N matrix A.
[out]BF
          BF is DOUBLE PRECISION array, dimension (LDB,N)
          Details of the GQR factorization of A and B, as returned
          by DGGRQF, see SGGRQF for further details.
[out]Z
          Z is DOUBLE PRECISION array, dimension (LDB,P)
          The P-by-P orthogonal matrix Z.
[out]T
          T is DOUBLE PRECISION array, dimension (LDB,max(P,N))
[out]BWK
          BWK is DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (min(P,N))
          The scalar factors of the elementary reflectors, as returned
          by DGGRQF.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (LWORK)
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK, LWORK >= max(M,P,N)**2.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (M)
[out]RESULT
          RESULT is DOUBLE PRECISION 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 dgrqts.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, N, P
183* ..
184* .. Array Arguments ..
185 DOUBLE PRECISION A( LDA, * ), AF( LDA, * ), B( LDB, * ),
186 $ BF( LDB, * ), BWK( LDB, * ), Q( LDA, * ),
187 $ R( LDA, * ), RESULT( 4 ), RWORK( * ),
188 $ T( LDB, * ), TAUA( * ), TAUB( * ),
189 $ WORK( LWORK ), Z( LDB, * )
190* ..
191*
192* =====================================================================
193*
194* .. Parameters ..
195 DOUBLE PRECISION ZERO, ONE
196 parameter( zero = 0.0d+0, one = 1.0d+0 )
197 DOUBLE PRECISION ROGUE
198 parameter( rogue = -1.0d+10 )
199* ..
200* .. Local Scalars ..
201 INTEGER INFO
202 DOUBLE PRECISION ANORM, BNORM, RESID, ULP, UNFL
203* ..
204* .. External Functions ..
205 DOUBLE PRECISION DLAMCH, DLANGE, DLANSY
206 EXTERNAL dlamch, dlange, dlansy
207* ..
208* .. External Subroutines ..
209 EXTERNAL dgemm, dggrqf, dlacpy, dlaset, dorgqr, dorgrq,
210 $ dsyrk
211* ..
212* .. Intrinsic Functions ..
213 INTRINSIC dble, max, min
214* ..
215* .. Executable Statements ..
216*
217 ulp = dlamch( 'Precision' )
218 unfl = dlamch( 'Safe minimum' )
219*
220* Copy the matrix A to the array AF.
221*
222 CALL dlacpy( 'Full', m, n, a, lda, af, lda )
223 CALL dlacpy( 'Full', p, n, b, ldb, bf, ldb )
224*
225 anorm = max( dlange( '1', m, n, a, lda, rwork ), unfl )
226 bnorm = max( dlange( '1', p, n, b, ldb, rwork ), unfl )
227*
228* Factorize the matrices A and B in the arrays AF and BF.
229*
230 CALL dggrqf( m, p, n, af, lda, taua, bf, ldb, taub, work, lwork,
231 $ info )
232*
233* Generate the N-by-N matrix Q
234*
235 CALL dlaset( 'Full', n, n, rogue, rogue, q, lda )
236 IF( m.LE.n ) THEN
237 IF( m.GT.0 .AND. m.LT.n )
238 $ CALL dlacpy( 'Full', m, n-m, af, lda, q( n-m+1, 1 ), lda )
239 IF( m.GT.1 )
240 $ CALL dlacpy( 'Lower', m-1, m-1, af( 2, n-m+1 ), lda,
241 $ q( n-m+2, n-m+1 ), lda )
242 ELSE
243 IF( n.GT.1 )
244 $ CALL dlacpy( 'Lower', n-1, n-1, af( m-n+2, 1 ), lda,
245 $ q( 2, 1 ), lda )
246 END IF
247 CALL dorgrq( n, n, min( m, n ), q, lda, taua, work, lwork, info )
248*
249* Generate the P-by-P matrix Z
250*
251 CALL dlaset( 'Full', p, p, rogue, rogue, z, ldb )
252 IF( p.GT.1 )
253 $ CALL dlacpy( 'Lower', p-1, n, bf( 2, 1 ), ldb, z( 2, 1 ), ldb )
254 CALL dorgqr( p, p, min( p, n ), z, ldb, taub, work, lwork, info )
255*
256* Copy R
257*
258 CALL dlaset( 'Full', m, n, zero, zero, r, lda )
259 IF( m.LE.n ) THEN
260 CALL dlacpy( 'Upper', m, m, af( 1, n-m+1 ), lda, r( 1, n-m+1 ),
261 $ lda )
262 ELSE
263 CALL dlacpy( 'Full', m-n, n, af, lda, r, lda )
264 CALL dlacpy( 'Upper', n, n, af( m-n+1, 1 ), lda, r( m-n+1, 1 ),
265 $ lda )
266 END IF
267*
268* Copy T
269*
270 CALL dlaset( 'Full', p, n, zero, zero, t, ldb )
271 CALL dlacpy( 'Upper', p, n, bf, ldb, t, ldb )
272*
273* Compute R - A*Q'
274*
275 CALL dgemm( 'No transpose', 'Transpose', m, n, n, -one, a, lda, q,
276 $ lda, one, r, lda )
277*
278* Compute norm( R - A*Q' ) / ( MAX(M,N)*norm(A)*ULP ) .
279*
280 resid = dlange( '1', m, n, r, lda, rwork )
281 IF( anorm.GT.zero ) THEN
282 result( 1 ) = ( ( resid / dble( max( 1, m, n ) ) ) / anorm ) /
283 $ ulp
284 ELSE
285 result( 1 ) = zero
286 END IF
287*
288* Compute T*Q - Z'*B
289*
290 CALL dgemm( 'Transpose', 'No transpose', p, n, p, one, z, ldb, b,
291 $ ldb, zero, bwk, ldb )
292 CALL dgemm( 'No transpose', 'No transpose', p, n, n, one, t, ldb,
293 $ q, lda, -one, bwk, ldb )
294*
295* Compute norm( T*Q - Z'*B ) / ( MAX(P,N)*norm(A)*ULP ) .
296*
297 resid = dlange( '1', p, n, bwk, ldb, rwork )
298 IF( bnorm.GT.zero ) THEN
299 result( 2 ) = ( ( resid / dble( max( 1, p, m ) ) ) / bnorm ) /
300 $ ulp
301 ELSE
302 result( 2 ) = zero
303 END IF
304*
305* Compute I - Q*Q'
306*
307 CALL dlaset( 'Full', n, n, zero, one, r, lda )
308 CALL dsyrk( 'Upper', 'No Transpose', n, n, -one, q, lda, one, r,
309 $ lda )
310*
311* Compute norm( I - Q'*Q ) / ( N * ULP ) .
312*
313 resid = dlansy( '1', 'Upper', n, r, lda, rwork )
314 result( 3 ) = ( resid / dble( max( 1, n ) ) ) / ulp
315*
316* Compute I - Z'*Z
317*
318 CALL dlaset( 'Full', p, p, zero, one, t, ldb )
319 CALL dsyrk( 'Upper', 'Transpose', p, p, -one, z, ldb, one, t,
320 $ ldb )
321*
322* Compute norm( I - Z'*Z ) / ( P*ULP ) .
323*
324 resid = dlansy( '1', 'Upper', p, t, ldb, rwork )
325 result( 4 ) = ( resid / dble( max( 1, p ) ) ) / ulp
326*
327 RETURN
328*
329* End of DGRQTS
330*
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dggrqf(m, p, n, a, lda, taua, b, ldb, taub, work, lwork, info)
DGGRQF
Definition dggrqf.f:214
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 dorgqr(m, n, k, a, lda, tau, work, lwork, info)
DORGQR
Definition dorgqr.f:128
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: