LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zchkrq.f
Go to the documentation of this file.
1*> \brief \b ZCHKRQ
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 ZCHKRQ( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
12* NRHS, THRESH, TSTERR, NMAX, A, AF, AQ, AR, AC,
13* B, X, XACT, TAU, WORK, RWORK, IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NM, NMAX, NN, NNB, NOUT, NRHS
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NVAL( * ),
23* $ NXVAL( * )
24* DOUBLE PRECISION RWORK( * )
25* COMPLEX*16 A( * ), AC( * ), AF( * ), AQ( * ), AR( * ),
26* $ B( * ), TAU( * ), WORK( * ), X( * ), XACT( * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> ZCHKRQ tests ZGERQF, ZUNGRQ and ZUNMRQ.
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] NNB
74*> \verbatim
75*> NNB is INTEGER
76*> The number of values of NB and NX contained in the
77*> vectors NBVAL and NXVAL. The blocking parameters are used
78*> in pairs (NB,NX).
79*> \endverbatim
80*>
81*> \param[in] NBVAL
82*> \verbatim
83*> NBVAL is INTEGER array, dimension (NNB)
84*> The values of the blocksize NB.
85*> \endverbatim
86*>
87*> \param[in] NXVAL
88*> \verbatim
89*> NXVAL is INTEGER array, dimension (NNB)
90*> The values of the crossover point NX.
91*> \endverbatim
92*>
93*> \param[in] NRHS
94*> \verbatim
95*> NRHS is INTEGER
96*> The number of right hand side vectors to be generated for
97*> each linear system.
98*> \endverbatim
99*>
100*> \param[in] THRESH
101*> \verbatim
102*> THRESH is DOUBLE PRECISION
103*> The threshold value for the test ratios. A result is
104*> included in the output file if RESULT >= THRESH. To have
105*> every test ratio printed, use THRESH = 0.
106*> \endverbatim
107*>
108*> \param[in] TSTERR
109*> \verbatim
110*> TSTERR is LOGICAL
111*> Flag that indicates whether error exits are to be tested.
112*> \endverbatim
113*>
114*> \param[in] NMAX
115*> \verbatim
116*> NMAX is INTEGER
117*> The maximum value permitted for M or N, used in dimensioning
118*> the work arrays.
119*> \endverbatim
120*>
121*> \param[out] A
122*> \verbatim
123*> A is COMPLEX*16 array, dimension (NMAX*NMAX)
124*> \endverbatim
125*>
126*> \param[out] AF
127*> \verbatim
128*> AF is COMPLEX*16 array, dimension (NMAX*NMAX)
129*> \endverbatim
130*>
131*> \param[out] AQ
132*> \verbatim
133*> AQ is COMPLEX*16 array, dimension (NMAX*NMAX)
134*> \endverbatim
135*>
136*> \param[out] AR
137*> \verbatim
138*> AR is COMPLEX*16 array, dimension (NMAX*NMAX)
139*> \endverbatim
140*>
141*> \param[out] AC
142*> \verbatim
143*> AC is COMPLEX*16 array, dimension (NMAX*NMAX)
144*> \endverbatim
145*>
146*> \param[out] B
147*> \verbatim
148*> B is COMPLEX*16 array, dimension (NMAX*NRHS)
149*> \endverbatim
150*>
151*> \param[out] X
152*> \verbatim
153*> X is COMPLEX*16 array, dimension (NMAX*NRHS)
154*> \endverbatim
155*>
156*> \param[out] XACT
157*> \verbatim
158*> XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
159*> \endverbatim
160*>
161*> \param[out] TAU
162*> \verbatim
163*> TAU is COMPLEX*16 array, dimension (NMAX)
164*> \endverbatim
165*>
166*> \param[out] WORK
167*> \verbatim
168*> WORK is COMPLEX*16 array, dimension (NMAX*NMAX)
169*> \endverbatim
170*>
171*> \param[out] RWORK
172*> \verbatim
173*> RWORK is DOUBLE PRECISION array, dimension (NMAX)
174*> \endverbatim
175*>
176*> \param[out] IWORK
177*> \verbatim
178*> IWORK is INTEGER array, dimension (NMAX)
179*> \endverbatim
180*>
181*> \param[in] NOUT
182*> \verbatim
183*> NOUT is INTEGER
184*> The unit number for output.
185*> \endverbatim
186*
187* Authors:
188* ========
189*
190*> \author Univ. of Tennessee
191*> \author Univ. of California Berkeley
192*> \author Univ. of Colorado Denver
193*> \author NAG Ltd.
194*
195*> \ingroup complex16_lin
196*
197* =====================================================================
198 SUBROUTINE zchkrq( DOTYPE, NM, MVAL, NN, NVAL, NNB, NBVAL, NXVAL,
199 $ NRHS, THRESH, TSTERR, NMAX, A, AF, AQ, AR, AC,
200 $ B, X, XACT, TAU, WORK, RWORK, IWORK, NOUT )
201*
202* -- LAPACK test routine --
203* -- LAPACK is a software package provided by Univ. of Tennessee, --
204* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
205*
206* .. Scalar Arguments ..
207 LOGICAL TSTERR
208 INTEGER NM, NMAX, NN, NNB, NOUT, NRHS
209 DOUBLE PRECISION THRESH
210* ..
211* .. Array Arguments ..
212 LOGICAL DOTYPE( * )
213 INTEGER IWORK( * ), MVAL( * ), NBVAL( * ), NVAL( * ),
214 $ nxval( * )
215 DOUBLE PRECISION RWORK( * )
216 COMPLEX*16 A( * ), AC( * ), AF( * ), AQ( * ), AR( * ),
217 $ b( * ), tau( * ), work( * ), x( * ), xact( * )
218* ..
219*
220* =====================================================================
221*
222* .. Parameters ..
223 INTEGER NTESTS
224 PARAMETER ( NTESTS = 7 )
225 INTEGER NTYPES
226 parameter( ntypes = 8 )
227 DOUBLE PRECISION ZERO
228 parameter( zero = 0.0d0 )
229* ..
230* .. Local Scalars ..
231 CHARACTER DIST, TYPE
232 CHARACTER*3 PATH
233 INTEGER I, IK, IM, IMAT, IN, INB, INFO, K, KL, KU, LDA,
234 $ lwork, m, minmn, mode, n, nb, nerrs, nfail, nk,
235 $ nrun, nt, nx
236 DOUBLE PRECISION ANORM, CNDNUM
237* ..
238* .. Local Arrays ..
239 INTEGER ISEED( 4 ), ISEEDY( 4 ), KVAL( 4 )
240 DOUBLE PRECISION RESULT( NTESTS )
241* ..
242* .. External Subroutines ..
243 EXTERNAL alaerh, alahd, alasum, xlaenv, zerrrq, zgerqs,
245 $ zrqt02, zrqt03
246* ..
247* .. Intrinsic Functions ..
248 INTRINSIC max, min
249* ..
250* .. Scalars in Common ..
251 LOGICAL LERR, OK
252 CHARACTER*32 SRNAMT
253 INTEGER INFOT, NUNIT
254* ..
255* .. Common blocks ..
256 COMMON / infoc / infot, nunit, ok, lerr
257 COMMON / srnamc / srnamt
258* ..
259* .. Data statements ..
260 DATA iseedy / 1988, 1989, 1990, 1991 /
261* ..
262* .. Executable Statements ..
263*
264* Initialize constants and the random number seed.
265*
266 path( 1: 1 ) = 'Zomplex precision'
267 path( 2: 3 ) = 'RQ'
268 nrun = 0
269 nfail = 0
270 nerrs = 0
271 DO 10 i = 1, 4
272 iseed( i ) = iseedy( i )
273 10 CONTINUE
274*
275* Test the error exits
276*
277 IF( tsterr )
278 $ CALL zerrrq( path, nout )
279 infot = 0
280 CALL xlaenv( 2, 2 )
281*
282 lda = nmax
283 lwork = nmax*max( nmax, nrhs )
284*
285* Do for each value of M in MVAL.
286*
287 DO 70 im = 1, nm
288 m = mval( im )
289*
290* Do for each value of N in NVAL.
291*
292 DO 60 in = 1, nn
293 n = nval( in )
294 minmn = min( m, n )
295 DO 50 imat = 1, ntypes
296*
297* Do the tests only if DOTYPE( IMAT ) is true.
298*
299 IF( .NOT.dotype( imat ) )
300 $ GO TO 50
301*
302* Set up parameters with ZLATB4 and generate a test matrix
303* with ZLATMS.
304*
305 CALL zlatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
306 $ cndnum, dist )
307*
308 srnamt = 'ZLATMS'
309 CALL zlatms( m, n, dist, iseed, TYPE, rwork, mode,
310 $ cndnum, anorm, kl, ku, 'No packing', a, lda,
311 $ work, info )
312*
313* Check error code from ZLATMS.
314*
315 IF( info.NE.0 ) THEN
316 CALL alaerh( path, 'ZLATMS', info, 0, ' ', m, n, -1,
317 $ -1, -1, imat, nfail, nerrs, nout )
318 GO TO 50
319 END IF
320*
321* Set some values for K: the first value must be MINMN,
322* corresponding to the call of ZRQT01; other values are
323* used in the calls of ZRQT02, and must not exceed MINMN.
324*
325 kval( 1 ) = minmn
326 kval( 2 ) = 0
327 kval( 3 ) = 1
328 kval( 4 ) = minmn / 2
329 IF( minmn.EQ.0 ) THEN
330 nk = 1
331 ELSE IF( minmn.EQ.1 ) THEN
332 nk = 2
333 ELSE IF( minmn.LE.3 ) THEN
334 nk = 3
335 ELSE
336 nk = 4
337 END IF
338*
339* Do for each value of K in KVAL
340*
341 DO 40 ik = 1, nk
342 k = kval( ik )
343*
344* Do for each pair of values (NB,NX) in NBVAL and NXVAL.
345*
346 DO 30 inb = 1, nnb
347 nb = nbval( inb )
348 CALL xlaenv( 1, nb )
349 nx = nxval( inb )
350 CALL xlaenv( 3, nx )
351 DO i = 1, ntests
352 result( i ) = zero
353 END DO
354 nt = 2
355 IF( ik.EQ.1 ) THEN
356*
357* Test ZGERQF
358*
359 CALL zrqt01( m, n, a, af, aq, ar, lda, tau,
360 $ work, lwork, rwork, result( 1 ) )
361 ELSE IF( m.LE.n ) THEN
362*
363* Test ZUNGRQ, using factorization
364* returned by ZRQT01
365*
366 CALL zrqt02( m, n, k, a, af, aq, ar, lda, tau,
367 $ work, lwork, rwork, result( 1 ) )
368 END IF
369 IF( m.GE.k ) THEN
370*
371* Test ZUNMRQ, using factorization returned
372* by ZRQT01
373*
374 CALL zrqt03( m, n, k, af, ac, ar, aq, lda, tau,
375 $ work, lwork, rwork, result( 3 ) )
376 nt = nt + 4
377*
378* If M>=N and K=N, call ZGERQS to solve a system
379* with NRHS right hand sides and compute the
380* residual.
381*
382 IF( k.EQ.m .AND. inb.EQ.1 ) THEN
383*
384* Generate a solution and set the right
385* hand side.
386*
387 srnamt = 'ZLARHS'
388 CALL zlarhs( path, 'New', 'Full',
389 $ 'No transpose', m, n, 0, 0,
390 $ nrhs, a, lda, xact, lda, b, lda,
391 $ iseed, info )
392*
393 CALL zlacpy( 'Full', m, nrhs, b, lda,
394 $ x( n-m+1 ), lda )
395 srnamt = 'ZGERQS'
396 CALL zgerqs( m, n, nrhs, af, lda, tau, x,
397 $ lda, work, lwork, info )
398*
399* Check error code from ZGERQS.
400*
401 IF( info.NE.0 )
402 $ CALL alaerh( path, 'ZGERQS', info, 0, ' ',
403 $ m, n, nrhs, -1, nb, imat,
404 $ nfail, nerrs, nout )
405*
406 CALL zget02( 'No transpose', m, n, nrhs, a,
407 $ lda, x, lda, b, lda, rwork,
408 $ result( 7 ) )
409 nt = nt + 1
410 END IF
411 END IF
412*
413* Print information about the tests that did not
414* pass the threshold.
415*
416 DO 20 i = 1, nt
417 IF( result( i ).GE.thresh ) THEN
418 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
419 $ CALL alahd( nout, path )
420 WRITE( nout, fmt = 9999 )m, n, k, nb, nx,
421 $ imat, i, result( i )
422 nfail = nfail + 1
423 END IF
424 20 CONTINUE
425 nrun = nrun + nt
426 30 CONTINUE
427 40 CONTINUE
428 50 CONTINUE
429 60 CONTINUE
430 70 CONTINUE
431*
432* Print a summary of the results.
433*
434 CALL alasum( path, nout, nfail, nrun, nerrs )
435*
436 9999 FORMAT( ' M=', i5, ', N=', i5, ', K=', i5, ', NB=', i4, ', NX=',
437 $ i5, ', type ', i2, ', test(', i2, ')=', g12.5 )
438 RETURN
439*
440* End of ZCHKRQ
441*
442 END
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine zget02(trans, m, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZGET02
Definition zget02.f:134
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.f:208
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 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 zchkrq(dotype, nm, mval, nn, nval, nnb, nbval, nxval, nrhs, thresh, tsterr, nmax, a, af, aq, ar, ac, b, x, xact, tau, work, rwork, iwork, nout)
ZCHKRQ
Definition zchkrq.f:201
subroutine zerrrq(path, nunit)
ZERRRQ
Definition zerrrq.f:55
subroutine zgerqs(m, n, nrhs, a, lda, tau, b, ldb, work, lwork, info)
ZGERQS
Definition zgerqs.f:122
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
subroutine zrqt01(m, n, a, af, q, r, lda, tau, work, lwork, rwork, result)
ZRQT01
Definition zrqt01.f:126
subroutine zrqt02(m, n, k, a, af, q, r, lda, tau, work, lwork, rwork, result)
ZRQT02
Definition zrqt02.f:136
subroutine zrqt03(m, n, k, af, c, cc, q, lda, tau, work, lwork, rwork, result)
ZRQT03
Definition zrqt03.f:136