LAPACK  3.10.1
LAPACK: Linear Algebra PACKage

◆ dqrt01p()

subroutine dqrt01p ( 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 
)

DQRT01P

Purpose:
 DQRT01P tests DGEQRFP, which computes the QR factorization of an m-by-n
 matrix A, and partially tests DORGQR which forms the m-by-m
 orthogonal matrix Q.

 DQRT01P compares R 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 QR factorization of A, as returned by DGEQRFP.
          See DGEQRFP for further details.
[out]Q
          Q is DOUBLE PRECISION array, dimension (LDA,M)
          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, 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 DGEQRFP.
[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( R - 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.

Definition at line 124 of file dqrt01p.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, dgeqrfp, dlacpy, dlaset, dorgqr, 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 = 'DGEQRFP'
180  CALL dgeqrfp( m, n, af, lda, tau, work, lwork, info )
181 *
182 * Copy details of Q
183 *
184  CALL dlaset( 'Full', m, m, rogue, rogue, q, lda )
185  CALL dlacpy( 'Lower', m-1, n, af( 2, 1 ), lda, q( 2, 1 ), lda )
186 *
187 * Generate the m-by-m matrix Q
188 *
189  srnamt = 'DORGQR'
190  CALL dorgqr( m, m, minmn, q, lda, tau, work, lwork, info )
191 *
192 * Copy R
193 *
194  CALL dlaset( 'Full', m, n, zero, zero, r, lda )
195  CALL dlacpy( 'Upper', m, n, af, lda, r, lda )
196 *
197 * Compute R - Q'*A
198 *
199  CALL dgemm( 'Transpose', 'No transpose', m, n, m, -one, q, lda, a,
200  $ lda, one, r, lda )
201 *
202 * Compute norm( R - Q'*A ) / ( M * norm(A) * EPS ) .
203 *
204  anorm = dlange( '1', m, n, a, lda, rwork )
205  resid = dlange( '1', m, n, r, lda, rwork )
206  IF( anorm.GT.zero ) THEN
207  result( 1 ) = ( ( resid / dble( max( 1, m ) ) ) / anorm ) / eps
208  ELSE
209  result( 1 ) = zero
210  END IF
211 *
212 * Compute I - Q'*Q
213 *
214  CALL dlaset( 'Full', m, m, zero, one, r, lda )
215  CALL dsyrk( 'Upper', 'Transpose', m, m, -one, q, lda, one, r,
216  $ lda )
217 *
218 * Compute norm( I - Q'*Q ) / ( M * EPS ) .
219 *
220  resid = dlansy( '1', 'Upper', m, r, lda, rwork )
221 *
222  result( 2 ) = ( resid / dble( max( 1, m ) ) ) / eps
223 *
224  RETURN
225 *
226 * End of DQRT01P
227 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
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
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 dsyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
DSYRK
Definition: dsyrk.f:169
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:187
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
subroutine dgeqrfp(M, N, A, LDA, TAU, WORK, LWORK, INFO)
DGEQRFP
Definition: dgeqrfp.f:149
subroutine dorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
DORGQR
Definition: dorgqr.f:128
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
Here is the call graph for this function:
Here is the caller graph for this function: