LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
November 2015

Definition at line 139 of file zchktz.f.

139 *
140 * -- LAPACK test routine (version 3.6.0) --
141 * -- LAPACK is a software package provided by Univ. of Tennessee, --
142 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
143 * November 2015
144 *
145 * .. Scalar Arguments ..
146  LOGICAL tsterr
147  INTEGER nm, nn, nout
148  DOUBLE PRECISION thresh
149 * ..
150 * .. Array Arguments ..
151  LOGICAL dotype( * )
152  INTEGER mval( * ), nval( * )
153  DOUBLE PRECISION s( * ), rwork( * )
154  COMPLEX*16 a( * ), copya( * ), tau( * ), work( * )
155 * ..
156 *
157 * =====================================================================
158 *
159 * .. Parameters ..
160  INTEGER ntypes
161  parameter ( ntypes = 3 )
162  INTEGER ntests
163  parameter ( ntests = 3 )
164  DOUBLE PRECISION one, zero
165  parameter ( one = 1.0d0, zero = 0.0d0 )
166 * ..
167 * .. Local Scalars ..
168  CHARACTER*3 path
169  INTEGER i, im, imode, in, info, k, lda, lwork, m,
170  $ mnmin, mode, n, nerrs, nfail, nrun
171  DOUBLE PRECISION eps
172 * ..
173 * .. Local Arrays ..
174  INTEGER iseed( 4 ), iseedy( 4 )
175  DOUBLE PRECISION result( ntests )
176 * ..
177 * .. External Functions ..
178  DOUBLE PRECISION dlamch, zqrt12, zrzt01, zrzt02
179  EXTERNAL dlamch, zqrt12, zrzt01, zrzt02
180 * ..
181 * .. External Subroutines ..
182  EXTERNAL alahd, alasum, dlaord, zerrtz, zgeqr2, zlacpy,
183  $ zlaset, zlatms, ztzrzf
184 * ..
185 * .. Intrinsic Functions ..
186  INTRINSIC dcmplx, max, min
187 * ..
188 * .. Scalars in Common ..
189  LOGICAL lerr, ok
190  CHARACTER*32 srnamt
191  INTEGER infot, iounit
192 * ..
193 * .. Common blocks ..
194  COMMON / infoc / infot, iounit, ok, lerr
195  COMMON / srnamc / srnamt
196 * ..
197 * .. Data statements ..
198  DATA iseedy / 1988, 1989, 1990, 1991 /
199 * ..
200 * .. Executable Statements ..
201 *
202 * Initialize constants and the random number seed.
203 *
204  path( 1: 1 ) = 'Zomplex precision'
205  path( 2: 3 ) = 'TZ'
206  nrun = 0
207  nfail = 0
208  nerrs = 0
209  DO 10 i = 1, 4
210  iseed( i ) = iseedy( i )
211  10 CONTINUE
212  eps = dlamch( 'Epsilon' )
213 *
214 * Test the error exits
215 *
216  IF( tsterr )
217  $ CALL zerrtz( path, nout )
218  infot = 0
219 *
220  DO 70 im = 1, nm
221 *
222 * Do for each value of M in MVAL.
223 *
224  m = mval( im )
225  lda = max( 1, m )
226 *
227  DO 60 in = 1, nn
228 *
229 * Do for each value of N in NVAL for which M .LE. N.
230 *
231  n = nval( in )
232  mnmin = min( m, n )
233  lwork = max( 1, n*n+4*m+n )
234 *
235  IF( m.LE.n ) THEN
236  DO 50 imode = 1, ntypes
237  IF( .NOT.dotype( imode ) )
238  $ GO TO 50
239 *
240 * Do for each type of singular value distribution.
241 * 0: zero matrix
242 * 1: one small singular value
243 * 2: exponential distribution
244 *
245  mode = imode - 1
246 *
247 * Test ZTZRQF
248 *
249 * Generate test matrix of size m by n using
250 * singular value distribution indicated by `mode'.
251 *
252  IF( mode.EQ.0 ) THEN
253  CALL zlaset( 'Full', m, n, dcmplx( zero ),
254  $ dcmplx( zero ), a, lda )
255  DO 30 i = 1, mnmin
256  s( i ) = zero
257  30 CONTINUE
258  ELSE
259  CALL zlatms( m, n, 'Uniform', iseed,
260  $ 'Nonsymmetric', s, imode,
261  $ one / eps, one, m, n, 'No packing', a,
262  $ lda, work, info )
263  CALL zgeqr2( m, n, a, lda, work, work( mnmin+1 ),
264  $ info )
265  CALL zlaset( 'Lower', m-1, n, dcmplx( zero ),
266  $ dcmplx( zero ), a( 2 ), lda )
267  CALL dlaord( 'Decreasing', mnmin, s, 1 )
268  END IF
269 *
270 * Save A and its singular values
271 *
272  CALL zlacpy( 'All', m, n, a, lda, copya, lda )
273 *
274 * Call ZTZRZF to reduce the upper trapezoidal matrix to
275 * upper triangular form.
276 *
277  srnamt = 'ZTZRZF'
278  CALL ztzrzf( m, n, a, lda, tau, work, lwork, info )
279 *
280 * Compute norm(svd(a) - svd(r))
281 *
282  result( 1 ) = zqrt12( m, m, a, lda, s, work,
283  $ lwork, rwork )
284 *
285 * Compute norm( A - R*Q )
286 *
287  result( 2 ) = zrzt01( m, n, copya, a, lda, tau, work,
288  $ lwork )
289 *
290 * Compute norm(Q'*Q - I).
291 *
292  result( 3 ) = zrzt02( m, n, a, lda, tau, work, lwork )
293 *
294 * Print information about the tests that did not pass
295 * the threshold.
296 *
297  DO 40 k = 1, ntests
298  IF( result( k ).GE.thresh ) THEN
299  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
300  $ CALL alahd( nout, path )
301  WRITE( nout, fmt = 9999 )m, n, imode, k,
302  $ result( k )
303  nfail = nfail + 1
304  END IF
305  40 CONTINUE
306  nrun = nrun + 3
307  50 CONTINUE
308  END IF
309  60 CONTINUE
310  70 CONTINUE
311 *
312 * Print a summary of the results.
313 *
314  CALL alasum( path, nout, nfail, nrun, nerrs )
315 *
316  9999 FORMAT( ' M =', i5, ', N =', i5, ', type ', i2, ', test ', i2,
317  $ ', ratio =', g12.5 )
318 *
319 * End if ZCHKTZ
320 *
subroutine ztzrzf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZTZRZF
Definition: ztzrzf.f:153
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
double precision function zqrt12(M, N, A, LDA, S, WORK, LWORK, RWORK)
ZQRT12
Definition: zqrt12.f:99
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
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:108
subroutine dlaord(JOB, N, X, INCX)
DLAORD
Definition: dlaord.f:75
double precision function zrzt02(M, N, AF, LDA, TAU, WORK, LWORK)
ZRZT02
Definition: zrzt02.f:93
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:123
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
subroutine zerrtz(PATH, NUNIT)
ZERRTZ
Definition: zerrtz.f:56
double precision function zrzt01(M, N, A, AF, LDA, TAU, WORK, LWORK)
ZRZT01
Definition: zrzt01.f:100
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75

Here is the call graph for this function:

Here is the caller graph for this function: