LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zchkqp3rk.f
Go to the documentation of this file.
1*> \brief \b ZCHKQP3RK
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE ZCHKQP3RK( DOTYPE, NM, MVAL, NN, NVAL, NNS, NSVAL,
12* $ NNB, NBVAL, NXVAL, THRESH, A, COPYA,
13* $ B, COPYB, S, TAU,
14* $ WORK, RWORK, IWORK, NOUT )
15* IMPLICIT NONE
16*
17* .. Scalar Arguments ..
18* INTEGER NM, NN, NNB, NOUT
19* DOUBLE PRECISION THRESH
20* ..
21* .. Array Arguments ..
22* LOGICAL DOTYPE( * )
23* INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NVAL( * ),
24* $ NXVAL( * )
25* DOUBLE PRECISION S( * ), RWORK( * )
26* COMPLEX*16 A( * ), COPYA( * ), TAU( * ), WORK( * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> ZCHKQP3RK tests ZGEQP3RK.
36*> \endverbatim
37*
38* Arguments:
39* ==========
40*
41*> \param[in] DOTYPE
42*> \verbatim
43*> DOTYPE is LOGICAL array, dimension (NTYPES)
44*> The matrix types to be used for testing. Matrices of type j
45*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
46*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
47*> \endverbatim
48*>
49*> \param[in] NM
50*> \verbatim
51*> NM is INTEGER
52*> The number of values of M contained in the vector MVAL.
53*> \endverbatim
54*>
55*> \param[in] MVAL
56*> \verbatim
57*> MVAL is INTEGER array, dimension (NM)
58*> The values of the matrix row dimension M.
59*> \endverbatim
60*>
61*> \param[in] NN
62*> \verbatim
63*> NN is INTEGER
64*> The number of values of N contained in the vector NVAL.
65*> \endverbatim
66*>
67*> \param[in] NVAL
68*> \verbatim
69*> NVAL is INTEGER array, dimension (NN)
70*> The values of the matrix column dimension N.
71*> \endverbatim
72*>
73*> \param[in] NNS
74*> \verbatim
75*> NNS is INTEGER
76*> The number of values of NRHS contained in the vector NSVAL.
77*> \endverbatim
78*>
79*> \param[in] NSVAL
80*> \verbatim
81*> NSVAL is INTEGER array, dimension (NNS)
82*> The values of the number of right hand sides NRHS.
83*> \endverbatim
84*> \param[in] NNB
85*> \verbatim
86*> NNB is INTEGER
87*> The number of values of NB and NX contained in the
88*> vectors NBVAL and NXVAL. The blocking parameters are used
89*> in pairs (NB,NX).
90*> \endverbatim
91*>
92*> \param[in] NBVAL
93*> \verbatim
94*> NBVAL is INTEGER array, dimension (NNB)
95*> The values of the blocksize NB.
96*> \endverbatim
97*>
98*> \param[in] NXVAL
99*> \verbatim
100*> NXVAL is INTEGER array, dimension (NNB)
101*> The values of the crossover point NX.
102*> \endverbatim
103*>
104*> \param[in] THRESH
105*> \verbatim
106*> THRESH is DOUBLE PRECISION
107*> The threshold value for the test ratios. A result is
108*> included in the output file if RESULT >= THRESH. To have
109*> every test ratio printed, use THRESH = 0.
110*> \endverbatim
111*>
112*> \param[out] A
113*> \verbatim
114*> A is COMPLEX*16 array, dimension (MMAX*NMAX)
115*> where MMAX is the maximum value of M in MVAL and NMAX is the
116*> maximum value of N in NVAL.
117*> \endverbatim
118*>
119*> \param[out] COPYA
120*> \verbatim
121*> COPYA is COMPLEX*16 array, dimension (MMAX*NMAX)
122*> \endverbatim
123*>
124*> \param[out] B
125*> \verbatim
126*> B is COMPLEX*16 array, dimension (MMAX*NSMAX)
127*> where MMAX is the maximum value of M in MVAL and NSMAX is the
128*> maximum value of NRHS in NSVAL.
129*> \endverbatim
130*>
131*> \param[out] COPYB
132*> \verbatim
133*> COPYB is COMPLEX*16 array, dimension (MMAX*NSMAX)
134*> \endverbatim
135*>
136*> \param[out] S
137*> \verbatim
138*> S is DOUBLE PRECISION array, dimension
139*> (min(MMAX,NMAX))
140*> \endverbatim
141*>
142*> \param[out] TAU
143*> \verbatim
144*> TAU is COMPLEX*16 array, dimension (MMAX)
145*> \endverbatim
146*>
147*> \param[out] WORK
148*> \verbatim
149*> WORK is COMPLEX*16 array, dimension
150*> (max(M*max(M,N) + 4*min(M,N) + max(M,N)))
151*> \endverbatim
152*>
153*> \param[out] RWORK
154*> \verbatim
155*> RWORK is DOUBLE PRECISION array, dimension (4*NMAX)
156*> \endverbatim
157*>
158*> \param[out] IWORK
159*> \verbatim
160*> IWORK is INTEGER array, dimension (2*NMAX)
161*> \endverbatim
162*>
163*> \param[in] NOUT
164*> \verbatim
165*> NOUT is INTEGER
166*> The unit number for output.
167*> \endverbatim
168*
169* Authors:
170* ========
171*
172*> \author Univ. of Tennessee
173*> \author Univ. of California Berkeley
174*> \author Univ. of Colorado Denver
175*> \author NAG Ltd.
176*
177*> \ingroup complex16_lin
178*
179* =====================================================================
180 SUBROUTINE zchkqp3rk( DOTYPE, NM, MVAL, NN, NVAL, NNS, NSVAL,
181 $ NNB, NBVAL, NXVAL, THRESH, A, COPYA,
182 $ B, COPYB, S, TAU,
183 $ WORK, RWORK, IWORK, NOUT )
184 IMPLICIT NONE
185*
186* -- LAPACK test routine --
187* -- LAPACK is a software package provided by Univ. of Tennessee, --
188* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189*
190* .. Scalar Arguments ..
191 INTEGER NM, NN, NNB, NNS, NOUT
192 DOUBLE PRECISION THRESH
193* ..
194* .. Array Arguments ..
195 LOGICAL DOTYPE( * )
196 INTEGER IWORK( * ), NBVAL( * ), MVAL( * ), NVAL( * ),
197 $ NSVAL( * ), NXVAL( * )
198 DOUBLE PRECISION S( * ), RWORK( * )
199 COMPLEX*16 A( * ), COPYA( * ), B( * ), COPYB( * ),
200 $ TAU( * ), WORK( * )
201* ..
202*
203* =====================================================================
204*
205* .. Parameters ..
206 INTEGER NTYPES
207 PARAMETER ( NTYPES = 19 )
208 integer ntests
209 parameter( ntests = 5 )
210 DOUBLE PRECISION ONE, ZERO, BIGNUM
211 COMPLEX*16 CONE, CZERO
212 parameter( one = 1.0d+0, zero = 0.0d+0,
213 $ czero = ( 0.0d+0, 0.0d+0 ),
214 $ cone = ( 1.0d+0, 0.0d+0 ),
215 $ bignum = 1.0d+38 )
216* ..
217* .. Local Scalars ..
218 CHARACTER DIST, TYPE
219 CHARACTER*3 PATH
220 INTEGER I, IHIGH, ILOW, IM, IMAT, IN, INC_ZERO,
221 $ inb, ind_offset_gen,
222 $ ind_in, ind_out, ins, info,
223 $ istep, j, j_inc, j_first_nz, jb_zero,
224 $ kfact, kl, kmax, ku, lda, lw, lwork,
225 $ lwork_mqr, m, minmn, minmnb_gen, mode, n,
226 $ nb, nb_zero, nerrs, nfail, nb_gen, nrhs,
227 $ nrun, nx, t
228 DOUBLE PRECISION ANORM, CNDNUM, EPS, ABSTOL, RELTOL,
229 $ DTEMP, MAXC2NRMK, RELMAXC2NRMK
230* ..
231* .. Local Arrays ..
232 INTEGER ISEED( 4 ), ISEEDY( 4 )
233 DOUBLE PRECISION RESULT( NTESTS ), RDUMMY( 1 )
234* ..
235* .. External Functions ..
236 DOUBLE PRECISION DLAMCH, ZQPT01, ZQRT11, ZQRT12, ZLANGE
237 EXTERNAL DLAMCH, ZQPT01, ZQRT11, ZQRT12, ZLANGE
238* ..
239* .. External Subroutines ..
240 EXTERNAL alaerh, alahd, alasum, dlaord, icopy, zaxpy,
243* ..
244* .. Intrinsic Functions ..
245 INTRINSIC abs, dble, max, min, mod
246* ..
247* .. Scalars in Common ..
248 LOGICAL LERR, OK
249 CHARACTER*32 SRNAMT
250 INTEGER INFOT, IOUNIT, ZUNMQR_LWORK
251* ..
252* .. Common blocks ..
253 COMMON / infoc / infot, iounit, ok, lerr
254 COMMON / srnamc / srnamt
255* ..
256* .. Data statements ..
257 DATA iseedy / 1988, 1989, 1990, 1991 /
258* ..
259* .. Executable Statements ..
260*
261* Initialize constants and the random number seed.
262*
263 path( 1: 1 ) = 'Zomplex precision'
264 path( 2: 3 ) = 'QK'
265 nrun = 0
266 nfail = 0
267 nerrs = 0
268 DO i = 1, 4
269 iseed( i ) = iseedy( i )
270 END DO
271 eps = dlamch( 'Epsilon' )
272 infot = 0
273*
274 DO im = 1, nm
275*
276* Do for each value of M in MVAL.
277*
278 m = mval( im )
279 lda = max( 1, m )
280*
281 DO in = 1, nn
282*
283* Do for each value of N in NVAL.
284*
285 n = nval( in )
286 minmn = min( m, n )
287 lwork = max( 1, m*max( m, n )+4*minmn+max( m, n ),
288 $ m*n + 2*minmn + 4*n )
289*
290 DO ins = 1, nns
291 nrhs = nsval( ins )
292*
293* Set up parameters with ZLATB4 and generate
294* M-by-NRHS B matrix with ZLATMS.
295* IMAT = 14:
296* Random matrix, CNDNUM = 2, NORM = ONE,
297* MODE = 3 (geometric distribution of singular values).
298*
299 CALL zlatb4( path, 14, m, nrhs, TYPE, kl, ku, anorm,
300 $ mode, cndnum, dist )
301*
302 srnamt = 'ZLATMS'
303 CALL zlatms( m, nrhs, dist, iseed, TYPE, s, mode,
304 $ cndnum, anorm, kl, ku, 'No packing',
305 $ copyb, lda, work, info )
306*
307* Check error code from ZLATMS.
308*
309 IF( info.NE.0 ) THEN
310 CALL alaerh( path, 'ZLATMS', info, 0, ' ', m,
311 $ nrhs, -1, -1, -1, 6, nfail, nerrs,
312 $ nout )
313 cycle
314 END IF
315*
316 DO imat = 1, ntypes
317*
318* Do the tests only if DOTYPE( IMAT ) is true.
319*
320 IF( .NOT.dotype( imat ) )
321 $ cycle
322*
323* The type of distribution used to generate the random
324* eigen-/singular values:
325* ( 'S' for symmetric distribution ) => UNIFORM( -1, 1 )
326*
327* Do for each type of NON-SYMMETRIC matrix: CNDNUM NORM MODE
328* 1. Zero matrix
329* 2. Random, Diagonal, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
330* 3. Random, Upper triangular, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
331* 4. Random, Lower triangular, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
332* 5. Random, First column is zero, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
333* 6. Random, Last MINMN column is zero, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
334* 7. Random, Last N column is zero, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
335* 8. Random, Middle column in MINMN is zero, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
336* 9. Random, First half of MINMN columns are zero, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
337* 10. Random, Last columns are zero starting from MINMN/2+1, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
338* 11. Random, Half MINMN columns in the middle are zero starting
339* from MINMN/2-(MINMN/2)/2+1, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
340* 12. Random, Odd columns are ZERO, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
341* 13. Random, Even columns are ZERO, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
342* 14. Random, CNDNUM = 2 CNDNUM = 2 ONE 3 ( geometric distribution of singular values )
343* 15. Random, CNDNUM = sqrt(0.1/EPS) CNDNUM = BADC1 = sqrt(0.1/EPS) ONE 3 ( geometric distribution of singular values )
344* 16. Random, CNDNUM = 0.1/EPS CNDNUM = BADC2 = 0.1/EPS ONE 3 ( geometric distribution of singular values )
345* 17. Random, CNDNUM = 0.1/EPS, CNDNUM = BADC2 = 0.1/EPS ONE 2 ( one small singular value, S(N)=1/CNDNUM )
346* one small singular value S(N)=1/CNDNUM
347* 18. Random, CNDNUM = 2, scaled near underflow CNDNUM = 2 SMALL = SAFMIN
348* 19. Random, CNDNUM = 2, scaled near overflow CNDNUM = 2 LARGE = 1.0/( 0.25 * ( SAFMIN / EPS ) ) 3 ( geometric distribution of singular values )
349*
350 IF( imat.EQ.1 ) THEN
351*
352* Matrix 1: Zero matrix
353*
354 CALL zlaset( 'Full', m, n, czero, czero, copya, lda )
355 DO i = 1, minmn
356 s( i ) = zero
357 END DO
358*
359 ELSE IF( (imat.GE.2 .AND. imat.LE.4 )
360 $ .OR. (imat.GE.14 .AND. imat.LE.19 ) ) THEN
361*
362* Matrices 2-5.
363*
364* Set up parameters with DLATB4 and generate a test
365* matrix with ZLATMS.
366*
367 CALL zlatb4( path, imat, m, n, TYPE, kl, ku, anorm,
368 $ mode, cndnum, dist )
369*
370 srnamt = 'ZLATMS'
371 CALL zlatms( m, n, dist, iseed, TYPE, s, mode,
372 $ cndnum, anorm, kl, ku, 'No packing',
373 $ copya, lda, work, info )
374*
375* Check error code from ZLATMS.
376*
377 IF( info.NE.0 ) THEN
378 CALL alaerh( path, 'ZLATMS', info, 0, ' ', m, n,
379 $ -1, -1, -1, imat, nfail, nerrs,
380 $ nout )
381 cycle
382 END IF
383*
384 CALL dlaord( 'Decreasing', minmn, s, 1 )
385*
386 ELSE IF( minmn.GE.2
387 $ .AND. imat.GE.5 .AND. imat.LE.13 ) THEN
388*
389* Rectangular matrices 5-13 that contain zero columns,
390* only for matrices MINMN >=2.
391*
392* JB_ZERO is the column index of ZERO block.
393* NB_ZERO is the column block size of ZERO block.
394* NB_GEN is the column blcok size of the
395* generated block.
396* J_INC in the non_zero column index increment
397* for matrix 12 and 13.
398* J_FIRS_NZ is the index of the first non-zero
399* column.
400*
401 IF( imat.EQ.5 ) THEN
402*
403* First column is zero.
404*
405 jb_zero = 1
406 nb_zero = 1
407 nb_gen = n - nb_zero
408*
409 ELSE IF( imat.EQ.6 ) THEN
410*
411* Last column MINMN is zero.
412*
413 jb_zero = minmn
414 nb_zero = 1
415 nb_gen = n - nb_zero
416*
417 ELSE IF( imat.EQ.7 ) THEN
418*
419* Last column N is zero.
420*
421 jb_zero = n
422 nb_zero = 1
423 nb_gen = n - nb_zero
424*
425 ELSE IF( imat.EQ.8 ) THEN
426*
427* Middle column in MINMN is zero.
428*
429 jb_zero = minmn / 2 + 1
430 nb_zero = 1
431 nb_gen = n - nb_zero
432*
433 ELSE IF( imat.EQ.9 ) THEN
434*
435* First half of MINMN columns is zero.
436*
437 jb_zero = 1
438 nb_zero = minmn / 2
439 nb_gen = n - nb_zero
440*
441 ELSE IF( imat.EQ.10 ) THEN
442*
443* Last columns are zero columns,
444* starting from (MINMN / 2 + 1) column.
445*
446 jb_zero = minmn / 2 + 1
447 nb_zero = n - jb_zero + 1
448 nb_gen = n - nb_zero
449*
450 ELSE IF( imat.EQ.11 ) THEN
451*
452* Half of the columns in the middle of MINMN
453* columns is zero, starting from
454* MINMN/2 - (MINMN/2)/2 + 1 column.
455*
456 jb_zero = minmn / 2 - (minmn / 2) / 2 + 1
457 nb_zero = minmn / 2
458 nb_gen = n - nb_zero
459*
460 ELSE IF( imat.EQ.12 ) THEN
461*
462* Odd-numbered columns are zero,
463*
464 nb_gen = n / 2
465 nb_zero = n - nb_gen
466 j_inc = 2
467 j_first_nz = 2
468*
469 ELSE IF( imat.EQ.13 ) THEN
470*
471* Even-numbered columns are zero.
472*
473 nb_zero = n / 2
474 nb_gen = n - nb_zero
475 j_inc = 2
476 j_first_nz = 1
477*
478 END IF
479*
480*
481* 1) Set the first NB_ZERO columns in COPYA(1:M,1:N)
482* to zero.
483*
484 CALL zlaset( 'Full', m, nb_zero, czero, czero,
485 $ copya, lda )
486*
487* 2) Generate an M-by-(N-NB_ZERO) matrix with the
488* chosen singular value distribution
489* in COPYA(1:M,NB_ZERO+1:N).
490*
491 CALL zlatb4( path, imat, m, nb_gen, TYPE, kl, ku,
492 $ anorm, mode, cndnum, dist )
493*
494 srnamt = 'ZLATMS'
495*
496 ind_offset_gen = nb_zero * lda
497*
498 CALL zlatms( m, nb_gen, dist, iseed, TYPE, s, mode,
499 $ cndnum, anorm, kl, ku, 'No packing',
500 $ copya( ind_offset_gen + 1 ), lda,
501 $ work, info )
502*
503* Check error code from ZLATMS.
504*
505 IF( info.NE.0 ) THEN
506 CALL alaerh( path, 'ZLATMS', info, 0, ' ', m,
507 $ nb_gen, -1, -1, -1, imat, nfail,
508 $ nerrs, nout )
509 cycle
510 END IF
511*
512* 3) Swap the gererated colums from the right side
513* NB_GEN-size block in COPYA into correct column
514* positions.
515*
516 IF( imat.EQ.6
517 $ .OR. imat.EQ.7
518 $ .OR. imat.EQ.8
519 $ .OR. imat.EQ.10
520 $ .OR. imat.EQ.11 ) THEN
521*
522* Move by swapping the generated columns
523* from the right NB_GEN-size block from
524* (NB_ZERO+1:NB_ZERO+JB_ZERO)
525* into columns (1:JB_ZERO-1).
526*
527 DO j = 1, jb_zero-1, 1
528 CALL zswap( m,
529 $ copya( ( nb_zero+j-1)*lda+1), 1,
530 $ copya( (j-1)*lda + 1 ), 1 )
531 END DO
532*
533 ELSE IF( imat.EQ.12 .OR. imat.EQ.13 ) THEN
534*
535* ( IMAT = 12, Odd-numbered ZERO columns. )
536* Swap the generated columns from the right
537* NB_GEN-size block into the even zero colums in the
538* left NB_ZERO-size block.
539*
540* ( IMAT = 13, Even-numbered ZERO columns. )
541* Swap the generated columns from the right
542* NB_GEN-size block into the odd zero colums in the
543* left NB_ZERO-size block.
544*
545 DO j = 1, nb_gen, 1
546 ind_out = ( nb_zero+j-1 )*lda + 1
547 ind_in = ( j_inc*(j-1)+(j_first_nz-1) )*lda
548 $ + 1
549 CALL zswap( m,
550 $ copya( ind_out ), 1,
551 $ copya( ind_in), 1 )
552 END DO
553*
554 END IF
555*
556* 5) Order the singular values generated by
557* DLAMTS in decreasing order and add trailing zeros
558* that correspond to zero columns.
559* The total number of singular values is MINMN.
560*
561 minmnb_gen = min( m, nb_gen )
562*
563 CALL dlaord( 'Decreasing', minmnb_gen, s, 1 )
564
565 DO i = minmnb_gen+1, minmn
566 s( i ) = zero
567 END DO
568*
569 ELSE
570*
571* IF(MINMN.LT.2) skip this size for this matrix type.
572*
573 cycle
574 END IF
575*
576* Initialize a copy array for a pivot array for DGEQP3RK.
577*
578 DO i = 1, n
579 iwork( i ) = 0
580 END DO
581*
582 DO inb = 1, nnb
583*
584* Do for each pair of values (NB,NX) in NBVAL and NXVAL.
585*
586 nb = nbval( inb )
587 CALL xlaenv( 1, nb )
588 nx = nxval( inb )
589 CALL xlaenv( 3, nx )
590*
591* We do MIN(M,N)+1 because we need a test for KMAX > N,
592* when KMAX is larger than MIN(M,N), KMAX should be
593* KMAX = MIN(M,N)
594*
595 DO kmax = 0, min(m,n)+1
596*
597* Get a working copy of COPYA into A( 1:M,1:N ).
598* Get a working copy of COPYB into A( 1:M, (N+1):NRHS ).
599* Get a working copy of COPYB into into B( 1:M, 1:NRHS ).
600* Get a working copy of IWORK(1:N) awith zeroes into
601* which is going to be used as pivot array IWORK( N+1:2N ).
602* NOTE: IWORK(2N+1:3N) is going to be used as a WORK array
603* for the routine.
604*
605 CALL zlacpy( 'All', m, n, copya, lda, a, lda )
606 CALL zlacpy( 'All', m, nrhs, copyb, lda,
607 $ a( lda*n + 1 ), lda )
608 CALL zlacpy( 'All', m, nrhs, copyb, lda,
609 $ b, lda )
610 CALL icopy( n, iwork( 1 ), 1, iwork( n+1 ), 1 )
611*
612 abstol = -1.0
613 reltol = -1.0
614*
615* Compute the QR factorization with pivoting of A
616*
617 lw = max( 1, max( 2*n + nb*( n+nrhs+1 ),
618 $ 3*n + nrhs - 1 ) )
619*
620* Compute ZGEQP3RK factorization of A.
621*
622 srnamt = 'ZGEQP3RK'
623 CALL zgeqp3rk( m, n, nrhs, kmax, abstol, reltol,
624 $ a, lda, kfact, maxc2nrmk,
625 $ relmaxc2nrmk, iwork( n+1 ), tau,
626 $ work, lw, rwork, iwork( 2*n+1 ),
627 $ info )
628*
629* Check error code from ZGEQP3RK.
630*
631 IF( info.LT.0 )
632 $ CALL alaerh( path, 'ZGEQP3RK', info, 0, ' ',
633 $ m, n, nx, -1, nb, imat,
634 $ nfail, nerrs, nout )
635*
636 IF( kfact.EQ.minmn ) THEN
637*
638* Compute test 1:
639*
640* This test in only for the full rank factorization of
641* the matrix A.
642*
643* Array S(1:min(M,N)) contains svd(A) the sigular values
644* of the original matrix A in decreasing absolute value
645* order. The test computes svd(R), the vector sigular
646* values of the upper trapezoid of A(1:M,1:N) that
647* contains the factor R, in decreasing order. The test
648* returns the ratio:
649*
650* 2-norm(svd(R) - svd(A)) / ( max(M,N) * 2-norm(svd(A)) * EPS )
651*
652 result( 1 ) = zqrt12( m, n, a, lda, s, work,
653 $ lwork , rwork )
654*
655 DO t = 1, 1
656 IF( result( t ).GE.thresh ) THEN
657 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
658 $ CALL alahd( nout, path )
659 WRITE( nout, fmt = 9999 ) 'ZGEQP3RK', m, n,
660 $ nrhs, kmax, abstol, reltol, nb, nx,
661 $ imat, t, result( t )
662 nfail = nfail + 1
663 END IF
664 END DO
665 nrun = nrun + 1
666*
667* End test 1
668*
669 END IF
670
671* Compute test 2:
672*
673* The test returns the ratio:
674*
675* 1-norm( A*P - Q*R ) / ( max(M,N) * 1-norm(A) * EPS )
676*
677 result( 2 ) = zqpt01( m, n, kfact, copya, a, lda, tau,
678 $ iwork( n+1 ), work, lwork )
679*
680* Compute test 3:
681*
682* The test returns the ratio:
683*
684* 1-norm( Q**T * Q - I ) / ( M * EPS )
685*
686 result( 3 ) = zqrt11( m, kfact, a, lda, tau, work,
687 $ lwork )
688*
689* Print information about the tests that did not pass
690* the threshold.
691*
692 DO t = 2, 3
693 IF( result( t ).GE.thresh ) THEN
694 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
695 $ CALL alahd( nout, path )
696 WRITE( nout, fmt = 9999 ) 'ZGEQP3RK', m, n,
697 $ nrhs, kmax, abstol, reltol,
698 $ nb, nx, imat, t, result( t )
699 nfail = nfail + 1
700 END IF
701 END DO
702 nrun = nrun + 2
703*
704* Compute test 4:
705*
706* This test is only for the factorizations with the
707* rank greater than 2.
708* The elements on the diagonal of R should be non-
709* increasing.
710*
711* The test returns the ratio:
712*
713* Returns 1.0D+100 if abs(R(K+1,K+1)) > abs(R(K,K)),
714* K=1:KFACT-1
715*
716 IF( min(kfact, minmn).GE.2 ) THEN
717*
718 DO j = 1, kfact-1, 1
719*
720 dtemp = (( abs( a( (j-1)*m+j ) ) -
721 $ abs( a( (j)*m+j+1 ) ) ) /
722 $ abs( a(1) ) )
723*
724 IF( dtemp.LT.zero ) THEN
725 result( 4 ) = bignum
726 END IF
727*
728 END DO
729*
730* Print information about the tests that did not
731* pass the threshold.
732*
733 DO t = 4, 4
734 IF( result( t ).GE.thresh ) THEN
735 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
736 $ CALL alahd( nout, path )
737 WRITE( nout, fmt = 9999 ) 'ZGEQP3RK',
738 $ m, n, nrhs, kmax, abstol, reltol,
739 $ nb, nx, imat, t,
740 $ result( t )
741 nfail = nfail + 1
742 END IF
743 END DO
744 nrun = nrun + 1
745*
746* End test 4.
747*
748 END IF
749*
750* Compute test 5:
751*
752* This test in only for matrix A with min(M,N) > 0.
753*
754* The test returns the ratio:
755*
756* 1-norm(Q**T * B - Q**T * B ) /
757* ( M * EPS )
758*
759* (1) Compute B:=Q**T * B in the matrix B.
760*
761 IF( minmn.GT.0 ) THEN
762*
763 lwork_mqr = max(1, nrhs)
764 CALL zunmqr( 'Left', 'Conjugate transpose',
765 $ m, nrhs, kfact, a, lda, tau, b, lda,
766 $ work, lwork_mqr, info )
767*
768 DO i = 1, nrhs
769*
770* Compare N+J-th column of A and J-column of B.
771*
772 CALL zaxpy( m, -cone, a( ( n+i-1 )*lda+1 ), 1,
773 $ b( ( i-1 )*lda+1 ), 1 )
774 END DO
775*
776 result( 5 ) =
777 $ abs(
778 $ zlange( 'One-norm', m, nrhs, b, lda, rdummy ) /
779 $ ( dble( m )*dlamch( 'Epsilon' ) )
780 $ )
781*
782* Print information about the tests that did not pass
783* the threshold.
784*
785 DO t = 5, 5
786 IF( result( t ).GE.thresh ) THEN
787 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
788 $ CALL alahd( nout, path )
789 WRITE( nout, fmt = 9999 ) 'ZGEQP3RK', m, n,
790 $ nrhs, kmax, abstol, reltol,
791 $ nb, nx, imat, t, result( t )
792 nfail = nfail + 1
793 END IF
794 END DO
795 nrun = nrun + 1
796*
797* End compute test 5.
798*
799 END IF
800*
801* END DO KMAX = 1, MIN(M,N)+1
802*
803 END DO
804*
805* END DO for INB = 1, NNB
806*
807 END DO
808*
809* END DO for IMAT = 1, NTYPES
810*
811 END DO
812*
813* END DO for INS = 1, NNS
814*
815 END DO
816*
817* END DO for IN = 1, NN
818*
819 END DO
820*
821* END DO for IM = 1, NM
822*
823 END DO
824*
825* Print a summary of the results.
826*
827 CALL alasum( path, nout, nfail, nrun, nerrs )
828*
829 9999 FORMAT( 1x, a, ' M =', i5, ', N =', i5, ', NRHS =', i5,
830 $ ', KMAX =', i5, ', ABSTOL =', g12.5,
831 $ ', RELTOL =', g12.5, ', NB =', i4, ', NX =', i4,
832 $ ', type ', i2, ', test ', i2, ', ratio =', g12.5 )
833*
834* End of ZCHKQP3RK
835*
836 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
subroutine dlaord(job, n, x, incx)
DLAORD
Definition dlaord.f:73
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
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
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 zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
subroutine zunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
ZUNMQR
Definition zunmqr.f:167
subroutine icopy(n, sx, incx, sy, incy)
ICOPY
Definition icopy.f:75
subroutine zchkqp3rk(dotype, nm, mval, nn, nval, nns, nsval, nnb, nbval, nxval, thresh, a, copya, b, copyb, s, tau, work, rwork, iwork, nout)
ZCHKQP3RK
Definition zchkqp3rk.f:184
subroutine zgeqp3rk(m, n, nrhs, kmax, abstol, reltol, a, lda, k, maxc2nrmk, relmaxc2nrmk, jpiv, tau, work, lwork, rwork, iwork, info)
ZGEQP3RK computes a truncated Householder QR factorization with column pivoting of a complex m-by-n m...
Definition zgeqp3rk.f:591
subroutine zlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB4
Definition zlatb4.f:121
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332