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

◆ zchkq3()

subroutine zchkq3 ( logical, dimension( * ) dotype,
integer nm,
integer, dimension( * ) mval,
integer nn,
integer, dimension( * ) nval,
integer nnb,
integer, dimension( * ) nbval,
integer, dimension( * ) nxval,
double precision thresh,
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, dimension( * ) iwork,
integer nout )

ZCHKQ3

Purpose:
!>
!> ZCHKQ3 tests ZGEQP3.
!> 
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]NNB
!>          NNB is INTEGER
!>          The number of values of NB and NX contained in the
!>          vectors NBVAL and NXVAL.  The blocking parameters are used
!>          in pairs (NB,NX).
!> 
[in]NBVAL
!>          NBVAL is INTEGER array, dimension (NNB)
!>          The values of the blocksize NB.
!> 
[in]NXVAL
!>          NXVAL is INTEGER array, dimension (NNB)
!>          The values of the crossover point NX.
!> 
[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.
!> 
[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
!>                      (max(M*max(M,N) + 4*min(M,N) + max(M,N)))
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (4*NMAX)
!> 
[out]IWORK
!>          IWORK is INTEGER 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 155 of file zchkq3.f.

158*
159* -- LAPACK test routine --
160* -- LAPACK is a software package provided by Univ. of Tennessee, --
161* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
162*
163* .. Scalar Arguments ..
164 INTEGER NM, NN, NNB, NOUT
165 DOUBLE PRECISION THRESH
166* ..
167* .. Array Arguments ..
168 LOGICAL DOTYPE( * )
169 INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NVAL( * ),
170 $ NXVAL( * )
171 DOUBLE PRECISION S( * ), RWORK( * )
172 COMPLEX*16 A( * ), COPYA( * ), TAU( * ), WORK( * )
173* ..
174*
175* =====================================================================
176*
177* .. Parameters ..
178 INTEGER NTYPES
179 parameter( ntypes = 6 )
180 INTEGER NTESTS
181 parameter( ntests = 3 )
182 DOUBLE PRECISION ONE, ZERO
183 COMPLEX*16 CZERO
184 parameter( one = 1.0d+0, zero = 0.0d+0,
185 $ czero = ( 0.0d+0, 0.0d+0 ) )
186* ..
187* .. Local Scalars ..
188 CHARACTER*3 PATH
189 INTEGER I, IHIGH, ILOW, IM, IMODE, IN, INB, INFO,
190 $ ISTEP, K, LDA, LW, LWORK, M, MNMIN, MODE, N,
191 $ NB, NERRS, NFAIL, NRUN, NX
192 DOUBLE PRECISION EPS
193* ..
194* .. Local Arrays ..
195 INTEGER ISEED( 4 ), ISEEDY( 4 )
196 DOUBLE PRECISION RESULT( NTESTS )
197* ..
198* .. External Functions ..
199 DOUBLE PRECISION DLAMCH, ZQPT01, ZQRT11, ZQRT12
200 EXTERNAL dlamch, zqpt01, zqrt11, zqrt12
201* ..
202* .. External Subroutines ..
203 EXTERNAL alahd, alasum, dlaord, icopy, xlaenv, zgeqp3,
205* ..
206* .. Intrinsic Functions ..
207 INTRINSIC max, min
208* ..
209* .. Scalars in Common ..
210 LOGICAL LERR, OK
211 CHARACTER*32 SRNAMT
212 INTEGER INFOT, IOUNIT
213* ..
214* .. Common blocks ..
215 COMMON / infoc / infot, iounit, ok, lerr
216 COMMON / srnamc / srnamt
217* ..
218* .. Data statements ..
219 DATA iseedy / 1988, 1989, 1990, 1991 /
220* ..
221* .. Executable Statements ..
222*
223* Initialize constants and the random number seed.
224*
225 path( 1: 1 ) = 'Zomplex precision'
226 path( 2: 3 ) = 'Q3'
227 nrun = 0
228 nfail = 0
229 nerrs = 0
230 DO 10 i = 1, 4
231 iseed( i ) = iseedy( i )
232 10 CONTINUE
233 eps = dlamch( 'Epsilon' )
234 infot = 0
235*
236 DO 90 im = 1, nm
237*
238* Do for each value of M in MVAL.
239*
240 m = mval( im )
241 lda = max( 1, m )
242*
243 DO 80 in = 1, nn
244*
245* Do for each value of N in NVAL.
246*
247 n = nval( in )
248 mnmin = min( m, n )
249 lwork = max( 1, m*max( m, n )+4*mnmin+max( m, n ) )
250*
251 DO 70 imode = 1, ntypes
252 IF( .NOT.dotype( imode ) )
253 $ GO TO 70
254*
255* Do for each type of matrix
256* 1: zero matrix
257* 2: one small singular value
258* 3: geometric distribution of singular values
259* 4: first n/2 columns fixed
260* 5: last n/2 columns fixed
261* 6: every second column fixed
262*
263 mode = imode
264 IF( imode.GT.3 )
265 $ mode = 1
266*
267* Generate test matrix of size m by n using
268* singular value distribution indicated by `mode'.
269*
270 DO 20 i = 1, n
271 iwork( i ) = 0
272 20 CONTINUE
273 IF( imode.EQ.1 ) THEN
274 CALL zlaset( 'Full', m, n, czero, czero, copya, lda )
275 DO 30 i = 1, mnmin
276 s( i ) = zero
277 30 CONTINUE
278 ELSE
279 CALL zlatms( m, n, 'Uniform', iseed, 'Nonsymm', s,
280 $ mode, one / eps, one, m, n, 'No packing',
281 $ copya, lda, work, info )
282 IF( imode.GE.4 ) THEN
283 IF( imode.EQ.4 ) THEN
284 ilow = 1
285 istep = 1
286 ihigh = max( 1, n / 2 )
287 ELSE IF( imode.EQ.5 ) THEN
288 ilow = max( 1, n / 2 )
289 istep = 1
290 ihigh = n
291 ELSE IF( imode.EQ.6 ) THEN
292 ilow = 1
293 istep = 2
294 ihigh = n
295 END IF
296 DO 40 i = ilow, ihigh, istep
297 iwork( i ) = 1
298 40 CONTINUE
299 END IF
300 CALL dlaord( 'Decreasing', mnmin, s, 1 )
301 END IF
302*
303 DO 60 inb = 1, nnb
304*
305* Do for each pair of values (NB,NX) in NBVAL and NXVAL.
306*
307 nb = nbval( inb )
308 CALL xlaenv( 1, nb )
309 nx = nxval( inb )
310 CALL xlaenv( 3, nx )
311*
312* Save A and its singular values and a copy of
313* vector IWORK.
314*
315 CALL zlacpy( 'All', m, n, copya, lda, a, lda )
316 CALL icopy( n, iwork( 1 ), 1, iwork( n+1 ), 1 )
317*
318* Workspace needed.
319*
320 lw = nb*( n+1 )
321*
322 srnamt = 'ZGEQP3'
323 CALL zgeqp3( m, n, a, lda, iwork( n+1 ), tau, work,
324 $ lw, rwork, info )
325*
326* Compute norm(svd(a) - svd(r))
327*
328 result( 1 ) = zqrt12( m, n, a, lda, s, work,
329 $ lwork, rwork )
330*
331* Compute norm( A*P - Q*R )
332*
333 result( 2 ) = zqpt01( m, n, mnmin, copya, a, lda, tau,
334 $ iwork( n+1 ), work, lwork )
335*
336* Compute Q'*Q
337*
338 result( 3 ) = zqrt11( m, mnmin, a, lda, tau, work,
339 $ lwork )
340*
341* Print information about the tests that did not pass
342* the threshold.
343*
344 DO 50 k = 1, ntests
345 IF( result( k ).GE.thresh ) THEN
346 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
347 $ CALL alahd( nout, path )
348 WRITE( nout, fmt = 9999 )'ZGEQP3', m, n, nb,
349 $ imode, k, result( k )
350 nfail = nfail + 1
351 END IF
352 50 CONTINUE
353 nrun = nrun + ntests
354*
355 60 CONTINUE
356 70 CONTINUE
357 80 CONTINUE
358 90 CONTINUE
359*
360* Print a summary of the results.
361*
362 CALL alasum( path, nout, nfail, nrun, nerrs )
363*
364 9999 FORMAT( 1x, a, ' M =', i5, ', N =', i5, ', NB =', i4, ', type ',
365 $ i2, ', test ', i2, ', ratio =', g12.5 )
366*
367* End of ZCHKQ3
368*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine dlaord(job, n, x, incx)
DLAORD
Definition dlaord.f:73
subroutine zgeqp3(m, n, a, lda, jpvt, tau, work, lwork, rwork, info)
ZGEQP3
Definition zgeqp3.f:157
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
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:104
subroutine icopy(n, sx, incx, sy, incy)
ICOPY
Definition icopy.f:75
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 zqpt01(m, n, k, a, af, lda, tau, jpvt, work, lwork)
ZQPT01
Definition zqpt01.f:120
double precision function zqrt11(m, k, a, lda, tau, work, lwork)
ZQRT11
Definition zqrt11.f:98
double precision function zqrt12(m, n, a, lda, s, work, lwork, rwork)
ZQRT12
Definition zqrt12.f:97
Here is the call graph for this function:
Here is the caller graph for this function: