LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zchktz.f
Go to the documentation of this file.
1*> \brief \b ZCHKTZ
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 ZCHKTZ( DOTYPE, NM, MVAL, NN, NVAL, THRESH, TSTERR, A,
12* COPYA, S, TAU, WORK, RWORK, 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 S( * ), RWORK( * )
23* COMPLEX*16 A( * ), COPYA( * ), TAU( * ), WORK( * )
24* ..
25*
26*
27*> \par Purpose:
28* =============
29*>
30*> \verbatim
31*>
32*> ZCHKTZ tests ZTZRZF.
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 COMPLEX*16 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 COMPLEX*16 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 COMPLEX*16 array, dimension (MMAX)
105*> \endverbatim
106*>
107*> \param[out] WORK
108*> \verbatim
109*> WORK is COMPLEX*16 array, dimension
110*> (MMAX*NMAX + 4*NMAX + MMAX)
111*> \endverbatim
112*>
113*> \param[out] RWORK
114*> \verbatim
115*> RWORK is DOUBLE PRECISION array, dimension (2*NMAX)
116*> \endverbatim
117*>
118*> \param[in] NOUT
119*> \verbatim
120*> NOUT is INTEGER
121*> The unit number for output.
122*> \endverbatim
123*
124* Authors:
125* ========
126*
127*> \author Univ. of Tennessee
128*> \author Univ. of California Berkeley
129*> \author Univ. of Colorado Denver
130*> \author NAG Ltd.
131*
132*> \ingroup complex16_lin
133*
134* =====================================================================
135 SUBROUTINE zchktz( DOTYPE, NM, MVAL, NN, NVAL, THRESH, TSTERR, A,
136 $ COPYA, S, TAU, WORK, RWORK, NOUT )
137*
138* -- LAPACK test routine --
139* -- LAPACK is a software package provided by Univ. of Tennessee, --
140* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
141*
142* .. Scalar Arguments ..
143 LOGICAL TSTERR
144 INTEGER NM, NN, NOUT
145 DOUBLE PRECISION THRESH
146* ..
147* .. Array Arguments ..
148 LOGICAL DOTYPE( * )
149 INTEGER MVAL( * ), NVAL( * )
150 DOUBLE PRECISION S( * ), RWORK( * )
151 COMPLEX*16 A( * ), COPYA( * ), TAU( * ), WORK( * )
152* ..
153*
154* =====================================================================
155*
156* .. Parameters ..
157 INTEGER NTYPES
158 parameter( ntypes = 3 )
159 INTEGER NTESTS
160 parameter( ntests = 3 )
161 DOUBLE PRECISION ONE, ZERO
162 parameter( one = 1.0d0, zero = 0.0d0 )
163* ..
164* .. Local Scalars ..
165 CHARACTER*3 PATH
166 INTEGER I, IM, IMODE, IN, INFO, K, LDA, LWORK, M,
167 $ mnmin, mode, n, nerrs, nfail, nrun
168 DOUBLE PRECISION EPS
169* ..
170* .. Local Arrays ..
171 INTEGER ISEED( 4 ), ISEEDY( 4 )
172 DOUBLE PRECISION RESULT( NTESTS )
173* ..
174* .. External Functions ..
175 DOUBLE PRECISION DLAMCH, ZQRT12, ZRZT01, ZRZT02
176 EXTERNAL dlamch, zqrt12, zrzt01, zrzt02
177* ..
178* .. External Subroutines ..
179 EXTERNAL alahd, alasum, dlaord, zerrtz, zgeqr2, zlacpy,
181* ..
182* .. Intrinsic Functions ..
183 INTRINSIC dcmplx, max, min
184* ..
185* .. Scalars in Common ..
186 LOGICAL LERR, OK
187 CHARACTER*32 SRNAMT
188 INTEGER INFOT, IOUNIT
189* ..
190* .. Common blocks ..
191 COMMON / infoc / infot, iounit, ok, lerr
192 COMMON / srnamc / srnamt
193* ..
194* .. Data statements ..
195 DATA iseedy / 1988, 1989, 1990, 1991 /
196* ..
197* .. Executable Statements ..
198*
199* Initialize constants and the random number seed.
200*
201 path( 1: 1 ) = 'Zomplex precision'
202 path( 2: 3 ) = 'TZ'
203 nrun = 0
204 nfail = 0
205 nerrs = 0
206 DO 10 i = 1, 4
207 iseed( i ) = iseedy( i )
208 10 CONTINUE
209 eps = dlamch( 'Epsilon' )
210*
211* Test the error exits
212*
213 IF( tsterr )
214 $ CALL zerrtz( path, nout )
215 infot = 0
216*
217 DO 70 im = 1, nm
218*
219* Do for each value of M in MVAL.
220*
221 m = mval( im )
222 lda = max( 1, m )
223*
224 DO 60 in = 1, nn
225*
226* Do for each value of N in NVAL for which M .LE. N.
227*
228 n = nval( in )
229 mnmin = min( m, n )
230 lwork = max( 1, n*n+4*m+n )
231*
232 IF( m.LE.n ) THEN
233 DO 50 imode = 1, ntypes
234 IF( .NOT.dotype( imode ) )
235 $ GO TO 50
236*
237* Do for each type of singular value distribution.
238* 0: zero matrix
239* 1: one small singular value
240* 2: exponential distribution
241*
242 mode = imode - 1
243*
244* Test ZTZRQF
245*
246* Generate test matrix of size m by n using
247* singular value distribution indicated by `mode'.
248*
249 IF( mode.EQ.0 ) THEN
250 CALL zlaset( 'Full', m, n, dcmplx( zero ),
251 $ dcmplx( zero ), a, lda )
252 DO 30 i = 1, mnmin
253 s( i ) = zero
254 30 CONTINUE
255 ELSE
256 CALL zlatms( m, n, 'Uniform', iseed,
257 $ 'Nonsymmetric', s, imode,
258 $ one / eps, one, m, n, 'No packing', a,
259 $ lda, work, info )
260 CALL zgeqr2( m, n, a, lda, work, work( mnmin+1 ),
261 $ info )
262 CALL zlaset( 'Lower', m-1, n, dcmplx( zero ),
263 $ dcmplx( zero ), a( 2 ), lda )
264 CALL dlaord( 'Decreasing', mnmin, s, 1 )
265 END IF
266*
267* Save A and its singular values
268*
269 CALL zlacpy( 'All', m, n, a, lda, copya, lda )
270*
271* Call ZTZRZF to reduce the upper trapezoidal matrix to
272* upper triangular form.
273*
274 srnamt = 'ZTZRZF'
275 CALL ztzrzf( m, n, a, lda, tau, work, lwork, info )
276*
277* Compute norm(svd(a) - svd(r))
278*
279 result( 1 ) = zqrt12( m, m, a, lda, s, work,
280 $ lwork, rwork )
281*
282* Compute norm( A - R*Q )
283*
284 result( 2 ) = zrzt01( m, n, copya, a, lda, tau, work,
285 $ lwork )
286*
287* Compute norm(Q'*Q - I).
288*
289 result( 3 ) = zrzt02( m, n, a, lda, tau, work, lwork )
290*
291* Print information about the tests that did not pass
292* the threshold.
293*
294 DO 40 k = 1, ntests
295 IF( result( k ).GE.thresh ) THEN
296 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
297 $ CALL alahd( nout, path )
298 WRITE( nout, fmt = 9999 )m, n, imode, k,
299 $ result( k )
300 nfail = nfail + 1
301 END IF
302 40 CONTINUE
303 nrun = nrun + 3
304 50 CONTINUE
305 END IF
306 60 CONTINUE
307 70 CONTINUE
308*
309* Print a summary of the results.
310*
311 CALL alasum( path, nout, nfail, nrun, nerrs )
312*
313 9999 FORMAT( ' M =', i5, ', N =', i5, ', type ', i2, ', test ', i2,
314 $ ', ratio =', g12.5 )
315*
316* End if ZCHKTZ
317*
318 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine dlaord(job, n, x, incx)
DLAORD
Definition dlaord.f:73
subroutine zgeqr2(m, n, a, lda, tau, work, info)
ZGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition zgeqr2.f:130
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine ztzrzf(m, n, a, lda, tau, work, lwork, info)
ZTZRZF
Definition ztzrzf.f:151
subroutine zchktz(dotype, nm, mval, nn, nval, thresh, tsterr, a, copya, s, tau, work, rwork, nout)
ZCHKTZ
Definition zchktz.f:137
subroutine zerrtz(path, nunit)
ZERRTZ
Definition zerrtz.f:54
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332