LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dgqrts ( integer  N,
integer  M,
integer  P,
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 
)

DGQRTS

Purpose:
 DGQRTS tests DGGQRF, which computes the GQR factorization of an
 N-by-M matrix A and a N-by-P matrix B: A = Q*R and B = Q*T*Z.
Parameters
[in]N
          N is INTEGER
          The number of rows of the matrices A and B.  N >= 0.
[in]M
          M is INTEGER
          The number of columns of the matrix A.  M >= 0.
[in]P
          P is INTEGER
          The number of columns of the matrix B.  P >= 0.
[in]A
          A is DOUBLE PRECISION array, dimension (LDA,M)
          The N-by-M matrix A.
[out]AF
          AF is DOUBLE PRECISION array, dimension (LDA,N)
          Details of the GQR factorization of A and B, as returned
          by DGGQRF, see SGGQRF for further details.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDA,N)
          The M-by-M 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 DGGQRF.
[in]B
          B is DOUBLE PRECISION array, dimension (LDB,P)
          On entry, the N-by-P matrix A.
[out]BF
          BF is DOUBLE PRECISION array, dimension (LDB,N)
          Details of the GQR factorization of A and B, as returned
          by DGGQRF, see SGGQRF 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(N,M,P)**2.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(N,M,P))
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (4)
          The test ratios:
            RESULT(1) = norm( R - Q'*A ) / ( MAX(M,N)*norm(A)*ULP)
            RESULT(2) = norm( T*Z - Q'*B ) / (MAX(P,N)*norm(B)*ULP)
            RESULT(3) = norm( I - Q'*Q ) / ( M*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 dgqrts.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, n, p
186 * ..
187 * .. Array Arguments ..
188  DOUBLE PRECISION a( lda, * ), af( lda, * ), b( ldb, * ),
189  $ bf( ldb, * ), bwk( ldb, * ), q( lda, * ),
190  $ r( lda, * ), result( 4 ), rwork( * ),
191  $ t( ldb, * ), taua( * ), taub( * ),
192  $ work( lwork ), z( ldb, * )
193 * ..
194 *
195 * =====================================================================
196 *
197 * .. Parameters ..
198  DOUBLE PRECISION zero, one
199  parameter ( zero = 0.0d+0, one = 1.0d+0 )
200  DOUBLE PRECISION rogue
201  parameter ( rogue = -1.0d+10 )
202 * ..
203 * .. Local Scalars ..
204  INTEGER info
205  DOUBLE PRECISION anorm, bnorm, resid, ulp, unfl
206 * ..
207 * .. External Functions ..
208  DOUBLE PRECISION dlamch, dlange, dlansy
209  EXTERNAL dlamch, dlange, dlansy
210 * ..
211 * .. External Subroutines ..
212  EXTERNAL dgemm, dggqrf, dlacpy, dlaset, dorgqr, dorgrq,
213  $ dsyrk
214 * ..
215 * .. Intrinsic Functions ..
216  INTRINSIC dble, max, min
217 * ..
218 * .. Executable Statements ..
219 *
220  ulp = dlamch( 'Precision' )
221  unfl = dlamch( 'Safe minimum' )
222 *
223 * Copy the matrix A to the array AF.
224 *
225  CALL dlacpy( 'Full', n, m, a, lda, af, lda )
226  CALL dlacpy( 'Full', n, p, b, ldb, bf, ldb )
227 *
228  anorm = max( dlange( '1', n, m, a, lda, rwork ), unfl )
229  bnorm = max( dlange( '1', n, p, b, ldb, rwork ), unfl )
230 *
231 * Factorize the matrices A and B in the arrays AF and BF.
232 *
233  CALL dggqrf( n, m, p, af, lda, taua, bf, ldb, taub, work, lwork,
234  $ info )
235 *
236 * Generate the N-by-N matrix Q
237 *
238  CALL dlaset( 'Full', n, n, rogue, rogue, q, lda )
239  CALL dlacpy( 'Lower', n-1, m, af( 2, 1 ), lda, q( 2, 1 ), lda )
240  CALL dorgqr( n, n, min( n, m ), q, lda, taua, work, lwork, info )
241 *
242 * Generate the P-by-P matrix Z
243 *
244  CALL dlaset( 'Full', p, p, rogue, rogue, z, ldb )
245  IF( n.LE.p ) THEN
246  IF( n.GT.0 .AND. n.LT.p )
247  $ CALL dlacpy( 'Full', n, p-n, bf, ldb, z( p-n+1, 1 ), ldb )
248  IF( n.GT.1 )
249  $ CALL dlacpy( 'Lower', n-1, n-1, bf( 2, p-n+1 ), ldb,
250  $ z( p-n+2, p-n+1 ), ldb )
251  ELSE
252  IF( p.GT.1 )
253  $ CALL dlacpy( 'Lower', p-1, p-1, bf( n-p+2, 1 ), ldb,
254  $ z( 2, 1 ), ldb )
255  END IF
256  CALL dorgrq( p, p, min( n, p ), z, ldb, taub, work, lwork, info )
257 *
258 * Copy R
259 *
260  CALL dlaset( 'Full', n, m, zero, zero, r, lda )
261  CALL dlacpy( 'Upper', n, m, af, lda, r, lda )
262 *
263 * Copy T
264 *
265  CALL dlaset( 'Full', n, p, zero, zero, t, ldb )
266  IF( n.LE.p ) THEN
267  CALL dlacpy( 'Upper', n, n, bf( 1, p-n+1 ), ldb, t( 1, p-n+1 ),
268  $ ldb )
269  ELSE
270  CALL dlacpy( 'Full', n-p, p, bf, ldb, t, ldb )
271  CALL dlacpy( 'Upper', p, p, bf( n-p+1, 1 ), ldb, t( n-p+1, 1 ),
272  $ ldb )
273  END IF
274 *
275 * Compute R - Q'*A
276 *
277  CALL dgemm( 'Transpose', 'No transpose', n, m, n, -one, q, lda, a,
278  $ lda, one, r, lda )
279 *
280 * Compute norm( R - Q'*A ) / ( MAX(M,N)*norm(A)*ULP ) .
281 *
282  resid = dlange( '1', n, m, r, lda, rwork )
283  IF( anorm.GT.zero ) THEN
284  result( 1 ) = ( ( resid / dble( max( 1, m, n ) ) ) / anorm ) /
285  $ ulp
286  ELSE
287  result( 1 ) = zero
288  END IF
289 *
290 * Compute T*Z - Q'*B
291 *
292  CALL dgemm( 'No Transpose', 'No transpose', n, p, p, one, t, ldb,
293  $ z, ldb, zero, bwk, ldb )
294  CALL dgemm( 'Transpose', 'No transpose', n, p, n, -one, q, lda, b,
295  $ ldb, one, bwk, ldb )
296 *
297 * Compute norm( T*Z - Q'*B ) / ( MAX(P,N)*norm(A)*ULP ) .
298 *
299  resid = dlange( '1', n, p, bwk, ldb, rwork )
300  IF( bnorm.GT.zero ) THEN
301  result( 2 ) = ( ( resid / dble( max( 1, p, n ) ) ) / bnorm ) /
302  $ ulp
303  ELSE
304  result( 2 ) = zero
305  END IF
306 *
307 * Compute I - Q'*Q
308 *
309  CALL dlaset( 'Full', n, n, zero, one, r, lda )
310  CALL dsyrk( 'Upper', 'Transpose', n, n, -one, q, lda, one, r,
311  $ lda )
312 *
313 * Compute norm( I - Q'*Q ) / ( N * ULP ) .
314 *
315  resid = dlansy( '1', 'Upper', n, r, lda, rwork )
316  result( 3 ) = ( resid / dble( max( 1, n ) ) ) / ulp
317 *
318 * Compute I - Z'*Z
319 *
320  CALL dlaset( 'Full', p, p, zero, one, t, ldb )
321  CALL dsyrk( 'Upper', 'Transpose', p, p, -one, z, ldb, one, t,
322  $ ldb )
323 *
324 * Compute norm( I - Z'*Z ) / ( P*ULP ) .
325 *
326  resid = dlansy( '1', 'Upper', p, t, ldb, rwork )
327  result( 4 ) = ( resid / dble( max( 1, p ) ) ) / ulp
328 *
329  RETURN
330 *
331 * End of DGQRTS
332 *
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:112
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, or the element of largest absolute value of a real symmetric matrix.
Definition: dlansy.f:124
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine dsyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
DSYRK
Definition: dsyrk.f:171
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:116
subroutine dorgrq(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
DORGRQ
Definition: dorgrq.f:130
subroutine dorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
DORGQR
Definition: dorgqr.f:130
subroutine dggqrf(N, M, P, A, LDA, TAUA, B, LDB, TAUB, WORK, LWORK, INFO)
DGGQRF
Definition: dggqrf.f:217

Here is the call graph for this function:

Here is the caller graph for this function: