LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dchktz.f
Go to the documentation of this file.
1*> \brief \b DCHKTZ
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 DCHKTZ( DOTYPE, NM, MVAL, NN, NVAL, THRESH, TSTERR, A,
12* COPYA, S, TAU, WORK, NOUT )
13*
14* .. Scalar Arguments ..
15* LOGICAL TSTERR
16* INTEGER NM, NN, NOUT
17* DOUBLE PRECISION THRESH
18* ..
19* .. Array Arguments ..
20* LOGICAL DOTYPE( * )
21* INTEGER MVAL( * ), NVAL( * )
22* DOUBLE PRECISION A( * ), COPYA( * ), S( * ),
23* $ TAU( * ), WORK( * )
24* ..
25*
26*
27*> \par Purpose:
28* =============
29*>
30*> \verbatim
31*>
32*> DCHKTZ tests DTZRZF.
33*> \endverbatim
34*
35* Arguments:
36* ==========
37*
38*> \param[in] DOTYPE
39*> \verbatim
40*> DOTYPE is LOGICAL array, dimension (NTYPES)
41*> The matrix types to be used for testing. Matrices of type j
42*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
43*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
44*> \endverbatim
45*>
46*> \param[in] NM
47*> \verbatim
48*> NM is INTEGER
49*> The number of values of M contained in the vector MVAL.
50*> \endverbatim
51*>
52*> \param[in] MVAL
53*> \verbatim
54*> MVAL is INTEGER array, dimension (NM)
55*> The values of the matrix row dimension M.
56*> \endverbatim
57*>
58*> \param[in] NN
59*> \verbatim
60*> NN is INTEGER
61*> The number of values of N contained in the vector NVAL.
62*> \endverbatim
63*>
64*> \param[in] NVAL
65*> \verbatim
66*> NVAL is INTEGER array, dimension (NN)
67*> The values of the matrix column dimension N.
68*> \endverbatim
69*>
70*> \param[in] THRESH
71*> \verbatim
72*> THRESH is DOUBLE PRECISION
73*> The threshold value for the test ratios. A result is
74*> included in the output file if RESULT >= THRESH. To have
75*> every test ratio printed, use THRESH = 0.
76*> \endverbatim
77*>
78*> \param[in] TSTERR
79*> \verbatim
80*> TSTERR is LOGICAL
81*> Flag that indicates whether error exits are to be tested.
82*> \endverbatim
83*>
84*> \param[out] A
85*> \verbatim
86*> A is DOUBLE PRECISION array, dimension (MMAX*NMAX)
87*> where MMAX is the maximum value of M in MVAL and NMAX is the
88*> maximum value of N in NVAL.
89*> \endverbatim
90*>
91*> \param[out] COPYA
92*> \verbatim
93*> COPYA is DOUBLE PRECISION array, dimension (MMAX*NMAX)
94*> \endverbatim
95*>
96*> \param[out] S
97*> \verbatim
98*> S is DOUBLE PRECISION array, dimension
99*> (min(MMAX,NMAX))
100*> \endverbatim
101*>
102*> \param[out] TAU
103*> \verbatim
104*> TAU is DOUBLE PRECISION array, dimension (MMAX)
105*> \endverbatim
106*>
107*> \param[out] WORK
108*> \verbatim
109*> WORK is DOUBLE PRECISION array, dimension
110*> (MMAX*NMAX + 4*NMAX + MMAX)
111*> \endverbatim
112*>
113*> \param[in] NOUT
114*> \verbatim
115*> NOUT is INTEGER
116*> The unit number for output.
117*> \endverbatim
118*
119* Authors:
120* ========
121*
122*> \author Univ. of Tennessee
123*> \author Univ. of California Berkeley
124*> \author Univ. of Colorado Denver
125*> \author NAG Ltd.
126*
127*> \ingroup double_lin
128*
129* =====================================================================
130 SUBROUTINE dchktz( DOTYPE, NM, MVAL, NN, NVAL, THRESH, TSTERR, A,
131 $ COPYA, S, TAU, WORK, NOUT )
132*
133* -- LAPACK test routine --
134* -- LAPACK is a software package provided by Univ. of Tennessee, --
135* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136*
137* .. Scalar Arguments ..
138 LOGICAL TSTERR
139 INTEGER NM, NN, NOUT
140 DOUBLE PRECISION THRESH
141* ..
142* .. Array Arguments ..
143 LOGICAL DOTYPE( * )
144 INTEGER MVAL( * ), NVAL( * )
145 DOUBLE PRECISION A( * ), COPYA( * ), S( * ),
146 $ tau( * ), work( * )
147* ..
148*
149* =====================================================================
150*
151* .. Parameters ..
152 INTEGER NTYPES
153 parameter( ntypes = 3 )
154 INTEGER NTESTS
155 parameter( ntests = 3 )
156 DOUBLE PRECISION ONE, ZERO
157 parameter( one = 1.0d0, zero = 0.0d0 )
158* ..
159* .. Local Scalars ..
160 CHARACTER*3 PATH
161 INTEGER I, IM, IMODE, IN, INFO, K, LDA, LWORK, M,
162 $ mnmin, mode, n, nerrs, nfail, nrun
163 DOUBLE PRECISION EPS
164* ..
165* .. Local Arrays ..
166 INTEGER ISEED( 4 ), ISEEDY( 4 )
167 DOUBLE PRECISION RESULT( NTESTS )
168* ..
169* .. External Functions ..
170 DOUBLE PRECISION DLAMCH, DQRT12, DRZT01, DRZT02
171 EXTERNAL dlamch, dqrt12, drzt01, drzt02
172* ..
173* .. External Subroutines ..
174 EXTERNAL alahd, alasum, derrtz, dgeqr2, dlacpy, dlaord,
176* ..
177* .. Intrinsic Functions ..
178 INTRINSIC max, min
179* ..
180* .. Scalars in Common ..
181 LOGICAL LERR, OK
182 CHARACTER*32 SRNAMT
183 INTEGER INFOT, IOUNIT
184* ..
185* .. Common blocks ..
186 COMMON / infoc / infot, iounit, ok, lerr
187 COMMON / srnamc / srnamt
188* ..
189* .. Data statements ..
190 DATA iseedy / 1988, 1989, 1990, 1991 /
191* ..
192* .. Executable Statements ..
193*
194* Initialize constants and the random number seed.
195*
196 path( 1: 1 ) = 'Double precision'
197 path( 2: 3 ) = 'TZ'
198 nrun = 0
199 nfail = 0
200 nerrs = 0
201 DO 10 i = 1, 4
202 iseed( i ) = iseedy( i )
203 10 CONTINUE
204 eps = dlamch( 'Epsilon' )
205*
206* Test the error exits
207*
208 IF( tsterr )
209 $ CALL derrtz( path, nout )
210 infot = 0
211*
212 DO 70 im = 1, nm
213*
214* Do for each value of M in MVAL.
215*
216 m = mval( im )
217 lda = max( 1, m )
218*
219 DO 60 in = 1, nn
220*
221* Do for each value of N in NVAL for which M .LE. N.
222*
223 n = nval( in )
224 mnmin = min( m, n )
225 lwork = max( 1, n*n+4*m+n, m*n+2*mnmin+4*n )
226*
227 IF( m.LE.n ) THEN
228 DO 50 imode = 1, ntypes
229 IF( .NOT.dotype( imode ) )
230 $ GO TO 50
231*
232* Do for each type of singular value distribution.
233* 0: zero matrix
234* 1: one small singular value
235* 2: exponential distribution
236*
237 mode = imode - 1
238*
239* Test DTZRQF
240*
241* Generate test matrix of size m by n using
242* singular value distribution indicated by `mode'.
243*
244 IF( mode.EQ.0 ) THEN
245 CALL dlaset( 'Full', m, n, zero, zero, a, lda )
246 DO 30 i = 1, mnmin
247 s( i ) = zero
248 30 CONTINUE
249 ELSE
250 CALL dlatms( m, n, 'Uniform', iseed,
251 $ 'Nonsymmetric', s, imode,
252 $ one / eps, one, m, n, 'No packing', a,
253 $ lda, work, info )
254 CALL dgeqr2( m, n, a, lda, work, work( mnmin+1 ),
255 $ info )
256 CALL dlaset( 'Lower', m-1, n, zero, zero, a( 2 ),
257 $ lda )
258 CALL dlaord( 'Decreasing', mnmin, s, 1 )
259 END IF
260*
261* Save A and its singular values
262*
263 CALL dlacpy( 'All', m, n, a, lda, copya, lda )
264*
265* Call DTZRZF to reduce the upper trapezoidal matrix to
266* upper triangular form.
267*
268 srnamt = 'DTZRZF'
269 CALL dtzrzf( m, n, a, lda, tau, work, lwork, info )
270*
271* Compute norm(svd(a) - svd(r))
272*
273 result( 1 ) = dqrt12( m, m, a, lda, s, work,
274 $ lwork )
275*
276* Compute norm( A - R*Q )
277*
278 result( 2 ) = drzt01( m, n, copya, a, lda, tau, work,
279 $ lwork )
280*
281* Compute norm(Q'*Q - I).
282*
283 result( 3 ) = drzt02( m, n, a, lda, tau, work, lwork )
284*
285* Print information about the tests that did not pass
286* the threshold.
287*
288 DO 40 k = 1, ntests
289 IF( result( k ).GE.thresh ) THEN
290 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
291 $ CALL alahd( nout, path )
292 WRITE( nout, fmt = 9999 )m, n, imode, k,
293 $ result( k )
294 nfail = nfail + 1
295 END IF
296 40 CONTINUE
297 nrun = nrun + 3
298 50 CONTINUE
299 END IF
300 60 CONTINUE
301 70 CONTINUE
302*
303* Print a summary of the results.
304*
305 CALL alasum( path, nout, nfail, nrun, nerrs )
306*
307 9999 FORMAT( ' M =', i5, ', N =', i5, ', type ', i2, ', test ', i2,
308 $ ', ratio =', g12.5 )
309*
310* End if DCHKTZ
311*
312 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine dchktz(dotype, nm, mval, nn, nval, thresh, tsterr, a, copya, s, tau, work, nout)
DCHKTZ
Definition dchktz.f:132
subroutine derrtz(path, nunit)
DERRTZ
Definition derrtz.f:54
subroutine dlaord(job, n, x, incx)
DLAORD
Definition dlaord.f:73
subroutine dlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
DLATMS
Definition dlatms.f:321
subroutine dgeqr2(m, n, a, lda, tau, work, info)
DGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition dgeqr2.f:130
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 dtzrzf(m, n, a, lda, tau, work, lwork, info)
DTZRZF
Definition dtzrzf.f:151