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

DQLT01

Purpose:
 DQLT01 tests DGEQLF, which computes the QL factorization of an m-by-n
 matrix A, and partially tests DORGQL which forms the m-by-m
 orthogonal matrix Q.

 DQLT01 compares L with Q'*A, 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 QL factorization of A, as returned by DGEQLF.
          See DGEQLF for further details.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDA,M)
          The m-by-m orthogonal matrix Q.
[out]L
          L is DOUBLE PRECISION array, dimension (LDA,max(M,N))
[in]LDA
          LDA is INTEGER
          The leading dimension of the arrays A, AF, Q and R.
          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 DGEQLF.
[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 - Q'*A ) / ( M * norm(A) * EPS )
          RESULT(2) = norm( I - Q'*Q ) / ( M * EPS )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 128 of file dqlt01.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: