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

◆ zchkqr()

subroutine zchkqr ( logical, dimension( * )  dotype,
integer  nm,
integer, dimension( * )  mval,
integer  nn,
integer, dimension( * )  nval,
integer  nnb,
integer, dimension( * )  nbval,
integer, dimension( * )  nxval,
integer  nrhs,
double precision  thresh,
logical  tsterr,
integer  nmax,
complex*16, dimension( * )  a,
complex*16, dimension( * )  af,
complex*16, dimension( * )  aq,
complex*16, dimension( * )  ar,
complex*16, dimension( * )  ac,
complex*16, dimension( * )  b,
complex*16, dimension( * )  x,
complex*16, dimension( * )  xact,
complex*16, dimension( * )  tau,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

ZCHKQR

Purpose:
 ZCHKQR tests ZGEQRF, ZUNGQR and ZUNMQR.
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]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[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.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for M or N, used in dimensioning
          the work arrays.
[out]A
          A is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AF
          AF is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AQ
          AQ is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AR
          AR is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AC
          AC is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]TAU
          TAU is COMPLEX*16 array, dimension (NMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (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 198 of file zchkqr.f.

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 = 9 )
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 Functions ..
243 LOGICAL ZGENND
244 EXTERNAL zgennd
245* ..
246* .. External Subroutines ..
247 EXTERNAL alaerh, alahd, alasum, xlaenv, zerrqr, zgels,
250* ..
251* .. Intrinsic Functions ..
252 INTRINSIC max, min
253* ..
254* .. Scalars in Common ..
255 LOGICAL LERR, OK
256 CHARACTER*32 SRNAMT
257 INTEGER INFOT, NUNIT
258* ..
259* .. Common blocks ..
260 COMMON / infoc / infot, nunit, ok, lerr
261 COMMON / srnamc / srnamt
262* ..
263* .. Data statements ..
264 DATA iseedy / 1988, 1989, 1990, 1991 /
265* ..
266* .. Executable Statements ..
267*
268* Initialize constants and the random number seed.
269*
270 path( 1: 1 ) = 'Zomplex precision'
271 path( 2: 3 ) = 'QR'
272 nrun = 0
273 nfail = 0
274 nerrs = 0
275 DO 10 i = 1, 4
276 iseed( i ) = iseedy( i )
277 10 CONTINUE
278*
279* Test the error exits
280*
281 IF( tsterr )
282 $ CALL zerrqr( path, nout )
283 infot = 0
284 CALL xlaenv( 2, 2 )
285*
286 lda = nmax
287 lwork = nmax*max( nmax, nrhs )
288*
289* Do for each value of M in MVAL.
290*
291 DO 70 im = 1, nm
292 m = mval( im )
293*
294* Do for each value of N in NVAL.
295*
296 DO 60 in = 1, nn
297 n = nval( in )
298 minmn = min( m, n )
299 DO 50 imat = 1, ntypes
300*
301* Do the tests only if DOTYPE( IMAT ) is true.
302*
303 IF( .NOT.dotype( imat ) )
304 $ GO TO 50
305*
306* Set up parameters with ZLATB4 and generate a test matrix
307* with ZLATMS.
308*
309 CALL zlatb4( path, imat, m, n, TYPE, KL, KU, ANORM, MODE,
310 $ CNDNUM, DIST )
311*
312 srnamt = 'ZLATMS'
313 CALL zlatms( m, n, dist, iseed, TYPE, RWORK, MODE,
314 $ CNDNUM, ANORM, KL, KU, 'No packing', A, LDA,
315 $ WORK, INFO )
316*
317* Check error code from ZLATMS.
318*
319 IF( info.NE.0 ) THEN
320 CALL alaerh( path, 'ZLATMS', info, 0, ' ', m, n, -1,
321 $ -1, -1, imat, nfail, nerrs, nout )
322 GO TO 50
323 END IF
324*
325* Set some values for K: the first value must be MINMN,
326* corresponding to the call of ZQRT01; other values are
327* used in the calls of ZQRT02, and must not exceed MINMN.
328*
329 kval( 1 ) = minmn
330 kval( 2 ) = 0
331 kval( 3 ) = 1
332 kval( 4 ) = minmn / 2
333 IF( minmn.EQ.0 ) THEN
334 nk = 1
335 ELSE IF( minmn.EQ.1 ) THEN
336 nk = 2
337 ELSE IF( minmn.LE.3 ) THEN
338 nk = 3
339 ELSE
340 nk = 4
341 END IF
342*
343* Do for each value of K in KVAL
344*
345 DO 40 ik = 1, nk
346 k = kval( ik )
347*
348* Do for each pair of values (NB,NX) in NBVAL and NXVAL.
349*
350 DO 30 inb = 1, nnb
351 nb = nbval( inb )
352 CALL xlaenv( 1, nb )
353 nx = nxval( inb )
354 CALL xlaenv( 3, nx )
355 DO i = 1, ntests
356 result( i ) = zero
357 END DO
358 nt = 2
359 IF( ik.EQ.1 ) THEN
360*
361* Test ZGEQRF
362*
363 CALL zqrt01( m, n, a, af, aq, ar, lda, tau,
364 $ work, lwork, rwork, result( 1 ) )
365*
366* Test ZGEQRFP
367*
368 CALL zqrt01p( m, n, a, af, aq, ar, lda, tau,
369 $ work, lwork, rwork, result( 8 ) )
370
371 IF( .NOT. zgennd( m, n, af, lda ) )
372 $ result( 9 ) = 2*thresh
373 nt = nt + 1
374 ELSE IF( m.GE.n ) THEN
375*
376* Test ZUNGQR, using factorization
377* returned by ZQRT01
378*
379 CALL zqrt02( m, n, k, a, af, aq, ar, lda, tau,
380 $ work, lwork, rwork, result( 1 ) )
381 END IF
382 IF( m.GE.k ) THEN
383*
384* Test ZUNMQR, using factorization returned
385* by ZQRT01
386*
387 CALL zqrt03( m, n, k, af, ac, ar, aq, lda, tau,
388 $ work, lwork, rwork, result( 3 ) )
389 nt = nt + 4
390*
391* If M>=N and K=N, call ZGELS to solve a system
392* with NRHS right hand sides and compute the
393* residual.
394*
395 IF( k.EQ.n .AND. inb.EQ.1 ) THEN
396*
397* Generate a solution and set the right
398* hand side.
399*
400 srnamt = 'ZLARHS'
401 CALL zlarhs( path, 'New', 'Full',
402 $ 'No transpose', m, n, 0, 0,
403 $ nrhs, a, lda, xact, lda, b, lda,
404 $ iseed, info )
405*
406 CALL zlacpy( 'Full', m, nrhs, b, lda, x,
407 $ lda )
408*
409* Reset AF to the original matrix. ZGELS
410* factors the matrix before solving the system.
411*
412 CALL zlacpy( 'Full', m, n, a, lda, af, lda )
413*
414 srnamt = 'ZGELS'
415 CALL zgels( 'No transpose', m, n, nrhs, af,
416 $ lda, x, lda, work, lwork, info )
417*
418* Check error code from ZGELS.
419*
420 IF( info.NE.0 )
421 $ CALL alaerh( path, 'ZGELS', info, 0, 'N',
422 $ m, n, nrhs, -1, nb, imat,
423 $ nfail, nerrs, nout )
424*
425 CALL zget02( 'No transpose', m, n, nrhs, a,
426 $ lda, x, lda, b, lda, rwork,
427 $ result( 7 ) )
428 nt = nt + 1
429 END IF
430 END IF
431*
432* Print information about the tests that did not
433* pass the threshold.
434*
435 DO 20 i = 1, ntests
436 IF( result( i ).GE.thresh ) THEN
437 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
438 $ CALL alahd( nout, path )
439 WRITE( nout, fmt = 9999 )m, n, k, nb, nx,
440 $ imat, i, result( i )
441 nfail = nfail + 1
442 END IF
443 20 CONTINUE
444 nrun = nrun + ntests
445 30 CONTINUE
446 40 CONTINUE
447 50 CONTINUE
448 60 CONTINUE
449 70 CONTINUE
450*
451* Print a summary of the results.
452*
453 CALL alasum( path, nout, nfail, nrun, nerrs )
454*
455 9999 FORMAT( ' M=', i5, ', N=', i5, ', K=', i5, ', NB=', i4, ', NX=',
456 $ i5, ', type ', i2, ', test(', i2, ')=', g12.5 )
457 RETURN
458*
459* End of ZCHKQR
460*
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 zgels(trans, m, n, nrhs, a, lda, b, ldb, work, lwork, info)
ZGELS solves overdetermined or underdetermined systems for GE matrices
Definition zgels.f:182
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 zerrqr(path, nunit)
ZERRQR
Definition zerrqr.f:55
logical function zgennd(m, n, a, lda)
ZGENND
Definition zgennd.f:68
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 zqrt01(m, n, a, af, q, r, lda, tau, work, lwork, rwork, result)
ZQRT01
Definition zqrt01.f:126
subroutine zqrt01p(m, n, a, af, q, r, lda, tau, work, lwork, rwork, result)
ZQRT01P
Definition zqrt01p.f:126
subroutine zqrt02(m, n, k, a, af, q, r, lda, tau, work, lwork, rwork, result)
ZQRT02
Definition zqrt02.f:135
subroutine zqrt03(m, n, k, af, c, cc, q, lda, tau, work, lwork, rwork, result)
ZQRT03
Definition zqrt03.f:136
Here is the call graph for this function:
Here is the caller graph for this function: