LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zchktz()

subroutine zchktz ( logical, dimension( * )  dotype,
integer  nm,
integer, dimension( * )  mval,
integer  nn,
integer, dimension( * )  nval,
double precision  thresh,
logical  tsterr,
complex*16, dimension( * )  a,
complex*16, dimension( * )  copya,
double precision, dimension( * )  s,
complex*16, dimension( * )  tau,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer  nout 
)

ZCHKTZ

Purpose:
 ZCHKTZ tests ZTZRZF.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NM
          NM is INTEGER
          The number of values of M contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix column dimension N.
[in]THRESH
          THRESH is DOUBLE PRECISION
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[out]A
          A is COMPLEX*16 array, dimension (MMAX*NMAX)
          where MMAX is the maximum value of M in MVAL and NMAX is the
          maximum value of N in NVAL.
[out]COPYA
          COPYA is COMPLEX*16 array, dimension (MMAX*NMAX)
[out]S
          S is DOUBLE PRECISION array, dimension
                      (min(MMAX,NMAX))
[out]TAU
          TAU is COMPLEX*16 array, dimension (MMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (MMAX*NMAX + 4*NMAX + MMAX)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (2*NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 135 of file zchktz.f.

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*
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
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
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 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
double precision function zqrt12(m, n, a, lda, s, work, lwork, rwork)
ZQRT12
Definition zqrt12.f:97
double precision function zrzt01(m, n, a, af, lda, tau, work, lwork)
ZRZT01
Definition zrzt01.f:98
double precision function zrzt02(m, n, af, lda, tau, work, lwork)
ZRZT02
Definition zrzt02.f:91
Here is the call graph for this function:
Here is the caller graph for this function: