LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine schktz ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
integer  NN,
integer, dimension( * )  NVAL,
real  THRESH,
logical  TSTERR,
real, dimension( * )  A,
real, dimension( * )  COPYA,
real, dimension( * )  S,
real, dimension( * )  TAU,
real, dimension( * )  WORK,
integer  NOUT 
)

SCHKTZ

Purpose:
 SCHKTZ tests STZRZF.
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 REAL
          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 REAL 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 REAL array, dimension (MMAX*NMAX)
[out]S
          S is REAL array, dimension
                      (min(MMAX,NMAX))
[out]TAU
          TAU is REAL array, dimension (MMAX)
[out]WORK
          WORK is REAL array, dimension
                      (MMAX*NMAX + 4*NMAX + MMAX)
[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 134 of file schktz.f.

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