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

◆ schktz()

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.

Definition at line 130 of file schktz.f.

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 REAL THRESH
141* ..
142* .. Array Arguments ..
143 LOGICAL DOTYPE( * )
144 INTEGER MVAL( * ), NVAL( * )
145 REAL 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 REAL ONE, ZERO
157 parameter( one = 1.0e0, zero = 0.0e0 )
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 REAL EPS
164* ..
165* .. Local Arrays ..
166 INTEGER ISEED( 4 ), ISEEDY( 4 )
167 REAL RESULT( NTESTS )
168* ..
169* .. External Functions ..
170 REAL SLAMCH, SQRT12, SRZT01, SRZT02
171 EXTERNAL slamch, sqrt12, srzt01, srzt02
172* ..
173* .. External Subroutines ..
174 EXTERNAL alahd, alasum, serrtz, sgeqr2, slacpy, slaord,
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 ) = 'Single 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 = slamch( 'Epsilon' )
205*
206* Test the error exits
207*
208 IF( tsterr )
209 $ CALL serrtz( 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 STZRQF
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 slaset( 'Full', m, n, zero, zero, a, lda )
246 DO 30 i = 1, mnmin
247 s( i ) = zero
248 30 CONTINUE
249 ELSE
250 CALL slatms( m, n, 'Uniform', iseed,
251 $ 'Nonsymmetric', s, imode,
252 $ one / eps, one, m, n, 'No packing', a,
253 $ lda, work, info )
254 CALL sgeqr2( m, n, a, lda, work, work( mnmin+1 ),
255 $ info )
256 CALL slaset( 'Lower', m-1, n, zero, zero, a( 2 ),
257 $ lda )
258 CALL slaord( 'Decreasing', mnmin, s, 1 )
259 END IF
260*
261* Save A and its singular values
262*
263 CALL slacpy( 'All', m, n, a, lda, copya, lda )
264*
265* Call STZRZF to reduce the upper trapezoidal matrix to
266* upper triangular form.
267*
268 srnamt = 'STZRZF'
269 CALL stzrzf( m, n, a, lda, tau, work, lwork, info )
270*
271* Compute norm(svd(a) - svd(r))
272*
273 result( 1 ) = sqrt12( m, m, a, lda, s, work,
274 $ lwork )
275*
276* Compute norm( A - R*Q )
277*
278 result( 2 ) = srzt01( m, n, copya, a, lda, tau, work,
279 $ lwork )
280*
281* Compute norm(Q'*Q - I).
282*
283 result( 3 ) = srzt02( 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 SCHKTZ
311*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
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:128
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
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:108
subroutine stzrzf(m, n, a, lda, tau, work, lwork, info)
STZRZF
Definition stzrzf.f:149
subroutine serrtz(path, nunit)
SERRTZ
Definition serrtz.f:54
subroutine slaord(job, n, x, incx)
SLAORD
Definition slaord.f:73
subroutine slatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
SLATMS
Definition slatms.f:321
real function sqrt12(m, n, a, lda, s, work, lwork)
SQRT12
Definition sqrt12.f:89
real function srzt01(m, n, a, af, lda, tau, work, lwork)
SRZT01
Definition srzt01.f:98
real function srzt02(m, n, af, lda, tau, work, lwork)
SRZT02
Definition srzt02.f:91
Here is the call graph for this function:
Here is the caller graph for this function: