LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dqrt01p.f
Go to the documentation of this file.
1*> \brief \b DQRT01P
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE DQRT01P( M, N, A, AF, Q, R, LDA, TAU, WORK, LWORK,
12* RWORK, RESULT )
13*
14* .. Scalar Arguments ..
15* INTEGER LDA, LWORK, M, N
16* ..
17* .. Array Arguments ..
18* DOUBLE PRECISION A( LDA, * ), AF( LDA, * ), Q( LDA, * ),
19* $ R( LDA, * ), RESULT( * ), RWORK( * ), TAU( * ),
20* $ WORK( LWORK )
21* ..
22*
23*
24*> \par Purpose:
25* =============
26*>
27*> \verbatim
28*>
29*> DQRT01P tests DGEQRFP, which computes the QR factorization of an m-by-n
30*> matrix A, and partially tests DORGQR which forms the m-by-m
31*> orthogonal matrix Q.
32*>
33*> DQRT01P compares R with Q'*A, and checks that Q is orthogonal.
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] M
40*> \verbatim
41*> M is INTEGER
42*> The number of rows of the matrix A. M >= 0.
43*> \endverbatim
44*>
45*> \param[in] N
46*> \verbatim
47*> N is INTEGER
48*> The number of columns of the matrix A. N >= 0.
49*> \endverbatim
50*>
51*> \param[in] A
52*> \verbatim
53*> A is DOUBLE PRECISION array, dimension (LDA,N)
54*> The m-by-n matrix A.
55*> \endverbatim
56*>
57*> \param[out] AF
58*> \verbatim
59*> AF is DOUBLE PRECISION array, dimension (LDA,N)
60*> Details of the QR factorization of A, as returned by DGEQRFP.
61*> See DGEQRFP for further details.
62*> \endverbatim
63*>
64*> \param[out] Q
65*> \verbatim
66*> Q is DOUBLE PRECISION array, dimension (LDA,M)
67*> The m-by-m orthogonal matrix Q.
68*> \endverbatim
69*>
70*> \param[out] R
71*> \verbatim
72*> R is DOUBLE PRECISION array, dimension (LDA,max(M,N))
73*> \endverbatim
74*>
75*> \param[in] LDA
76*> \verbatim
77*> LDA is INTEGER
78*> The leading dimension of the arrays A, AF, Q and R.
79*> LDA >= max(M,N).
80*> \endverbatim
81*>
82*> \param[out] TAU
83*> \verbatim
84*> TAU is DOUBLE PRECISION array, dimension (min(M,N))
85*> The scalar factors of the elementary reflectors, as returned
86*> by DGEQRFP.
87*> \endverbatim
88*>
89*> \param[out] WORK
90*> \verbatim
91*> WORK is DOUBLE PRECISION array, dimension (LWORK)
92*> \endverbatim
93*>
94*> \param[in] LWORK
95*> \verbatim
96*> LWORK is INTEGER
97*> The dimension of the array WORK.
98*> \endverbatim
99*>
100*> \param[out] RWORK
101*> \verbatim
102*> RWORK is DOUBLE PRECISION array, dimension (M)
103*> \endverbatim
104*>
105*> \param[out] RESULT
106*> \verbatim
107*> RESULT is DOUBLE PRECISION array, dimension (2)
108*> The test ratios:
109*> RESULT(1) = norm( R - Q'*A ) / ( M * norm(A) * EPS )
110*> RESULT(2) = norm( I - Q'*Q ) / ( M * EPS )
111*> \endverbatim
112*
113* Authors:
114* ========
115*
116*> \author Univ. of Tennessee
117*> \author Univ. of California Berkeley
118*> \author Univ. of Colorado Denver
119*> \author NAG Ltd.
120*
121*> \ingroup double_lin
122*
123* =====================================================================
124 SUBROUTINE dqrt01p( M, N, A, AF, Q, R, LDA, TAU, WORK, LWORK,
125 $ RWORK, RESULT )
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*
228 END
subroutine dqrt01p(m, n, a, af, q, r, lda, tau, work, lwork, rwork, result)
DQRT01P
Definition dqrt01p.f:126
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dgeqrfp(m, n, a, lda, tau, work, lwork, info)
DGEQRFP
Definition dgeqrfp.f:149
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
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