LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 CUNMQR.
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.
Date
November 2015

Definition at line 203 of file zchkqr.f.

203 *
204 * -- LAPACK test routine (version 3.6.0) --
205 * -- LAPACK is a software package provided by Univ. of Tennessee, --
206 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
207 * November 2015
208 *
209 * .. Scalar Arguments ..
210  LOGICAL tsterr
211  INTEGER nm, nmax, nn, nnb, nout, nrhs
212  DOUBLE PRECISION thresh
213 * ..
214 * .. Array Arguments ..
215  LOGICAL dotype( * )
216  INTEGER iwork( * ), mval( * ), nbval( * ), nval( * ),
217  $ nxval( * )
218  DOUBLE PRECISION rwork( * )
219  COMPLEX*16 a( * ), ac( * ), af( * ), aq( * ), ar( * ),
220  $ b( * ), tau( * ), work( * ), x( * ), xact( * )
221 * ..
222 *
223 * =====================================================================
224 *
225 * .. Parameters ..
226  INTEGER ntests
227  parameter ( ntests = 9 )
228  INTEGER ntypes
229  parameter ( ntypes = 8 )
230  DOUBLE PRECISION zero
231  parameter ( zero = 0.0d0 )
232 * ..
233 * .. Local Scalars ..
234  CHARACTER dist, type
235  CHARACTER*3 path
236  INTEGER i, ik, im, imat, in, inb, info, k, kl, ku, lda,
237  $ lwork, m, minmn, mode, n, nb, nerrs, nfail, nk,
238  $ nrun, nt, nx
239  DOUBLE PRECISION anorm, cndnum
240 * ..
241 * .. Local Arrays ..
242  INTEGER iseed( 4 ), iseedy( 4 ), kval( 4 )
243  DOUBLE PRECISION result( ntests )
244 * ..
245 * .. External Functions ..
246  LOGICAL zgennd
247  EXTERNAL zgennd
248 * ..
249 * .. External Subroutines ..
250  EXTERNAL alaerh, alahd, alasum, xlaenv, zerrqr, zgeqrs,
252  $ zqrt01p, zqrt02, zqrt03
253 * ..
254 * .. Intrinsic Functions ..
255  INTRINSIC max, min
256 * ..
257 * .. Scalars in Common ..
258  LOGICAL lerr, ok
259  CHARACTER*32 srnamt
260  INTEGER infot, nunit
261 * ..
262 * .. Common blocks ..
263  COMMON / infoc / infot, nunit, ok, lerr
264  COMMON / srnamc / srnamt
265 * ..
266 * .. Data statements ..
267  DATA iseedy / 1988, 1989, 1990, 1991 /
268 * ..
269 * .. Executable Statements ..
270 *
271 * Initialize constants and the random number seed.
272 *
273  path( 1: 1 ) = 'Zomplex precision'
274  path( 2: 3 ) = 'QR'
275  nrun = 0
276  nfail = 0
277  nerrs = 0
278  DO 10 i = 1, 4
279  iseed( i ) = iseedy( i )
280  10 CONTINUE
281 *
282 * Test the error exits
283 *
284  IF( tsterr )
285  $ CALL zerrqr( path, nout )
286  infot = 0
287  CALL xlaenv( 2, 2 )
288 *
289  lda = nmax
290  lwork = nmax*max( nmax, nrhs )
291 *
292 * Do for each value of M in MVAL.
293 *
294  DO 70 im = 1, nm
295  m = mval( im )
296 *
297 * Do for each value of N in NVAL.
298 *
299  DO 60 in = 1, nn
300  n = nval( in )
301  minmn = min( m, n )
302  DO 50 imat = 1, ntypes
303 *
304 * Do the tests only if DOTYPE( IMAT ) is true.
305 *
306  IF( .NOT.dotype( imat ) )
307  $ GO TO 50
308 *
309 * Set up parameters with ZLATB4 and generate a test matrix
310 * with ZLATMS.
311 *
312  CALL zlatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
313  $ cndnum, dist )
314 *
315  srnamt = 'ZLATMS'
316  CALL zlatms( m, n, dist, iseed, TYPE, rwork, mode,
317  $ cndnum, anorm, kl, ku, 'No packing', a, lda,
318  $ work, info )
319 *
320 * Check error code from ZLATMS.
321 *
322  IF( info.NE.0 ) THEN
323  CALL alaerh( path, 'ZLATMS', info, 0, ' ', m, n, -1,
324  $ -1, -1, imat, nfail, nerrs, nout )
325  GO TO 50
326  END IF
327 *
328 * Set some values for K: the first value must be MINMN,
329 * corresponding to the call of ZQRT01; other values are
330 * used in the calls of ZQRT02, and must not exceed MINMN.
331 *
332  kval( 1 ) = minmn
333  kval( 2 ) = 0
334  kval( 3 ) = 1
335  kval( 4 ) = minmn / 2
336  IF( minmn.EQ.0 ) THEN
337  nk = 1
338  ELSE IF( minmn.EQ.1 ) THEN
339  nk = 2
340  ELSE IF( minmn.LE.3 ) THEN
341  nk = 3
342  ELSE
343  nk = 4
344  END IF
345 *
346 * Do for each value of K in KVAL
347 *
348  DO 40 ik = 1, nk
349  k = kval( ik )
350 *
351 * Do for each pair of values (NB,NX) in NBVAL and NXVAL.
352 *
353  DO 30 inb = 1, nnb
354  nb = nbval( inb )
355  CALL xlaenv( 1, nb )
356  nx = nxval( inb )
357  CALL xlaenv( 3, nx )
358  DO i = 1, ntests
359  result( i ) = zero
360  END DO
361  nt = 2
362  IF( ik.EQ.1 ) THEN
363 *
364 * Test ZGEQRF
365 *
366  CALL zqrt01( m, n, a, af, aq, ar, lda, tau,
367  $ work, lwork, rwork, result( 1 ) )
368 *
369 * Test ZGEQRFP
370 *
371  CALL zqrt01p( m, n, a, af, aq, ar, lda, tau,
372  $ work, lwork, rwork, result( 8 ) )
373 
374  IF( .NOT. zgennd( m, n, af, lda ) )
375  $ result( 9 ) = 2*thresh
376  nt = nt + 1
377  ELSE IF( m.GE.n ) THEN
378 *
379 * Test ZUNGQR, using factorization
380 * returned by ZQRT01
381 *
382  CALL zqrt02( m, n, k, a, af, aq, ar, lda, tau,
383  $ work, lwork, rwork, result( 1 ) )
384  END IF
385  IF( m.GE.k ) THEN
386 *
387 * Test ZUNMQR, using factorization returned
388 * by ZQRT01
389 *
390  CALL zqrt03( m, n, k, af, ac, ar, aq, lda, tau,
391  $ work, lwork, rwork, result( 3 ) )
392  nt = nt + 4
393 *
394 * If M>=N and K=N, call ZGEQRS to solve a system
395 * with NRHS right hand sides and compute the
396 * residual.
397 *
398  IF( k.EQ.n .AND. inb.EQ.1 ) THEN
399 *
400 * Generate a solution and set the right
401 * hand side.
402 *
403  srnamt = 'ZLARHS'
404  CALL zlarhs( path, 'New', 'Full',
405  $ 'No transpose', m, n, 0, 0,
406  $ nrhs, a, lda, xact, lda, b, lda,
407  $ iseed, info )
408 *
409  CALL zlacpy( 'Full', m, nrhs, b, lda, x,
410  $ lda )
411  srnamt = 'ZGEQRS'
412  CALL zgeqrs( m, n, nrhs, af, lda, tau, x,
413  $ lda, work, lwork, info )
414 *
415 * Check error code from ZGEQRS.
416 *
417  IF( info.NE.0 )
418  $ CALL alaerh( path, 'ZGEQRS', info, 0, ' ',
419  $ m, n, nrhs, -1, nb, imat,
420  $ nfail, nerrs, nout )
421 *
422  CALL zget02( 'No transpose', m, n, nrhs, a,
423  $ lda, x, lda, b, lda, rwork,
424  $ result( 7 ) )
425  nt = nt + 1
426  END IF
427  END IF
428 *
429 * Print information about the tests that did not
430 * pass the threshold.
431 *
432  DO 20 i = 1, ntests
433  IF( result( i ).GE.thresh ) THEN
434  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
435  $ CALL alahd( nout, path )
436  WRITE( nout, fmt = 9999 )m, n, k, nb, nx,
437  $ imat, i, result( i )
438  nfail = nfail + 1
439  END IF
440  20 CONTINUE
441  nrun = nrun + ntests
442  30 CONTINUE
443  40 CONTINUE
444  50 CONTINUE
445  60 CONTINUE
446  70 CONTINUE
447 *
448 * Print a summary of the results.
449 *
450  CALL alasum( path, nout, nfail, nrun, nerrs )
451 *
452  9999 FORMAT( ' M=', i5, ', N=', i5, ', K=', i5, ', NB=', i4, ', NX=',
453  $ i5, ', type ', i2, ', test(', i2, ')=', g12.5 )
454  RETURN
455 *
456 * End of ZCHKQR
457 *
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zqrt02(M, N, K, A, AF, Q, R, LDA, TAU, WORK, LWORK, RWORK, RESULT)
ZQRT02
Definition: zqrt02.f:137
subroutine zerrqr(PATH, NUNIT)
ZERRQR
Definition: zerrqr.f:57
subroutine zlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
ZLARHS
Definition: zlarhs.f:211
subroutine zqrt01p(M, N, A, AF, Q, R, LDA, TAU, WORK, LWORK, RWORK, RESULT)
ZQRT01P
Definition: zqrt01p.f:128
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
logical function zgennd(M, N, A, LDA)
ZGENND
Definition: zgennd.f:70
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zqrt03(M, N, K, AF, C, CC, Q, LDA, TAU, WORK, LWORK, RWORK, RESULT)
ZQRT03
Definition: zqrt03.f:138
subroutine zqrt01(M, N, A, AF, Q, R, LDA, TAU, WORK, LWORK, RWORK, RESULT)
ZQRT01
Definition: zqrt01.f:128
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
subroutine zget02(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZGET02
Definition: zget02.f:135
subroutine zgeqrs(M, N, NRHS, A, LDA, TAU, B, LDB, WORK, LWORK, INFO)
ZGEQRS
Definition: zgeqrs.f:123
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75

Here is the call graph for this function:

Here is the caller graph for this function: