LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zchkge ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
complex*16, dimension( * )  A,
complex*16, dimension( * )  AFAC,
complex*16, dimension( * )  AINV,
complex*16, dimension( * )  B,
complex*16, dimension( * )  X,
complex*16, dimension( * )  XACT,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

ZCHKGE

Purpose:
 ZCHKGE tests ZGETRF, -TRI, -TRS, -RFS, and -CON.
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 contained in the vector NBVAL.
[in]NBVAL
          NBVAL is INTEGER array, dimension (NBVAL)
          The values of the blocksize NB.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[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]AFAC
          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension
                      (max(2*NMAX,2*NSMAX+NWORK))
[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 2011

Definition at line 188 of file zchkge.f.

188 *
189 * -- LAPACK test routine (version 3.4.0) --
190 * -- LAPACK is a software package provided by Univ. of Tennessee, --
191 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
192 * November 2011
193 *
194 * .. Scalar Arguments ..
195  LOGICAL tsterr
196  INTEGER nm, nmax, nn, nnb, nns, nout
197  DOUBLE PRECISION thresh
198 * ..
199 * .. Array Arguments ..
200  LOGICAL dotype( * )
201  INTEGER iwork( * ), mval( * ), nbval( * ), nsval( * ),
202  $ nval( * )
203  DOUBLE PRECISION rwork( * )
204  COMPLEX*16 a( * ), afac( * ), ainv( * ), b( * ),
205  $ work( * ), x( * ), xact( * )
206 * ..
207 *
208 * =====================================================================
209 *
210 * .. Parameters ..
211  DOUBLE PRECISION one, zero
212  parameter ( one = 1.0d+0, zero = 0.0d+0 )
213  INTEGER ntypes
214  parameter ( ntypes = 11 )
215  INTEGER ntests
216  parameter ( ntests = 8 )
217  INTEGER ntran
218  parameter ( ntran = 3 )
219 * ..
220 * .. Local Scalars ..
221  LOGICAL trfcon, zerot
222  CHARACTER dist, norm, trans, TYPE, xtype
223  CHARACTER*3 path
224  INTEGER i, im, imat, in, inb, info, ioff, irhs, itran,
225  $ izero, k, kl, ku, lda, lwork, m, mode, n, nb,
226  $ nerrs, nfail, nimat, nrhs, nrun, nt
227  DOUBLE PRECISION ainvnm, anorm, anormi, anormo, cndnum, dummy,
228  $ rcond, rcondc, rcondi, rcondo
229 * ..
230 * .. Local Arrays ..
231  CHARACTER transs( ntran )
232  INTEGER iseed( 4 ), iseedy( 4 )
233  DOUBLE PRECISION result( ntests )
234 * ..
235 * .. External Functions ..
236  DOUBLE PRECISION dget06, zlange
237  EXTERNAL dget06, zlange
238 * ..
239 * .. External Subroutines ..
240  EXTERNAL alaerh, alahd, alasum, xlaenv, zerrge, zgecon,
243  $ zlatb4, zlatms
244 * ..
245 * .. Intrinsic Functions ..
246  INTRINSIC dcmplx, max, min
247 * ..
248 * .. Scalars in Common ..
249  LOGICAL lerr, ok
250  CHARACTER*32 srnamt
251  INTEGER infot, nunit
252 * ..
253 * .. Common blocks ..
254  COMMON / infoc / infot, nunit, ok, lerr
255  COMMON / srnamc / srnamt
256 * ..
257 * .. Data statements ..
258  DATA iseedy / 1988, 1989, 1990, 1991 / ,
259  $ transs / 'N', 'T', 'C' /
260 * ..
261 * .. Executable Statements ..
262 *
263 * Initialize constants and the random number seed.
264 *
265  path( 1: 1 ) = 'Zomplex precision'
266  path( 2: 3 ) = 'GE'
267  nrun = 0
268  nfail = 0
269  nerrs = 0
270  DO 10 i = 1, 4
271  iseed( i ) = iseedy( i )
272  10 CONTINUE
273 *
274 * Test the error exits
275 *
276  CALL xlaenv( 1, 1 )
277  IF( tsterr )
278  $ CALL zerrge( path, nout )
279  infot = 0
280  CALL xlaenv( 2, 2 )
281 *
282 * Do for each value of M in MVAL
283 *
284  DO 120 im = 1, nm
285  m = mval( im )
286  lda = max( 1, m )
287 *
288 * Do for each value of N in NVAL
289 *
290  DO 110 in = 1, nn
291  n = nval( in )
292  xtype = 'N'
293  nimat = ntypes
294  IF( m.LE.0 .OR. n.LE.0 )
295  $ nimat = 1
296 *
297  DO 100 imat = 1, nimat
298 *
299 * Do the tests only if DOTYPE( IMAT ) is true.
300 *
301  IF( .NOT.dotype( imat ) )
302  $ GO TO 100
303 *
304 * Skip types 5, 6, or 7 if the matrix size is too small.
305 *
306  zerot = imat.GE.5 .AND. imat.LE.7
307  IF( zerot .AND. n.LT.imat-4 )
308  $ GO TO 100
309 *
310 * Set up parameters with ZLATB4 and generate a test matrix
311 * with ZLATMS.
312 *
313  CALL zlatb4( path, imat, m, n, TYPE, kl, ku, anorm, mode,
314  $ cndnum, dist )
315 *
316  srnamt = 'ZLATMS'
317  CALL zlatms( m, n, dist, iseed, TYPE, rwork, mode,
318  $ cndnum, anorm, kl, ku, 'No packing', a, lda,
319  $ work, info )
320 *
321 * Check error code from ZLATMS.
322 *
323  IF( info.NE.0 ) THEN
324  CALL alaerh( path, 'ZLATMS', info, 0, ' ', m, n, -1,
325  $ -1, -1, imat, nfail, nerrs, nout )
326  GO TO 100
327  END IF
328 *
329 * For types 5-7, zero one or more columns of the matrix to
330 * test that INFO is returned correctly.
331 *
332  IF( zerot ) THEN
333  IF( imat.EQ.5 ) THEN
334  izero = 1
335  ELSE IF( imat.EQ.6 ) THEN
336  izero = min( m, n )
337  ELSE
338  izero = min( m, n ) / 2 + 1
339  END IF
340  ioff = ( izero-1 )*lda
341  IF( imat.LT.7 ) THEN
342  DO 20 i = 1, m
343  a( ioff+i ) = zero
344  20 CONTINUE
345  ELSE
346  CALL zlaset( 'Full', m, n-izero+1, dcmplx( zero ),
347  $ dcmplx( zero ), a( ioff+1 ), lda )
348  END IF
349  ELSE
350  izero = 0
351  END IF
352 *
353 * These lines, if used in place of the calls in the DO 60
354 * loop, cause the code to bomb on a Sun SPARCstation.
355 *
356 * ANORMO = ZLANGE( 'O', M, N, A, LDA, RWORK )
357 * ANORMI = ZLANGE( 'I', M, N, A, LDA, RWORK )
358 *
359 * Do for each blocksize in NBVAL
360 *
361  DO 90 inb = 1, nnb
362  nb = nbval( inb )
363  CALL xlaenv( 1, nb )
364 *
365 * Compute the LU factorization of the matrix.
366 *
367  CALL zlacpy( 'Full', m, n, a, lda, afac, lda )
368  srnamt = 'ZGETRF'
369  CALL zgetrf( m, n, afac, lda, iwork, info )
370 *
371 * Check error code from ZGETRF.
372 *
373  IF( info.NE.izero )
374  $ CALL alaerh( path, 'ZGETRF', info, izero, ' ', m,
375  $ n, -1, -1, nb, imat, nfail, nerrs,
376  $ nout )
377  trfcon = .false.
378 *
379 *+ TEST 1
380 * Reconstruct matrix from factors and compute residual.
381 *
382  CALL zlacpy( 'Full', m, n, afac, lda, ainv, lda )
383  CALL zget01( m, n, a, lda, ainv, lda, iwork, rwork,
384  $ result( 1 ) )
385  nt = 1
386 *
387 *+ TEST 2
388 * Form the inverse if the factorization was successful
389 * and compute the residual.
390 *
391  IF( m.EQ.n .AND. info.EQ.0 ) THEN
392  CALL zlacpy( 'Full', n, n, afac, lda, ainv, lda )
393  srnamt = 'ZGETRI'
394  nrhs = nsval( 1 )
395  lwork = nmax*max( 3, nrhs )
396  CALL zgetri( n, ainv, lda, iwork, work, lwork,
397  $ info )
398 *
399 * Check error code from ZGETRI.
400 *
401  IF( info.NE.0 )
402  $ CALL alaerh( path, 'ZGETRI', info, 0, ' ', n, n,
403  $ -1, -1, nb, imat, nfail, nerrs,
404  $ nout )
405 *
406 * Compute the residual for the matrix times its
407 * inverse. Also compute the 1-norm condition number
408 * of A.
409 *
410  CALL zget03( n, a, lda, ainv, lda, work, lda,
411  $ rwork, rcondo, result( 2 ) )
412  anormo = zlange( 'O', m, n, a, lda, rwork )
413 *
414 * Compute the infinity-norm condition number of A.
415 *
416  anormi = zlange( 'I', m, n, a, lda, rwork )
417  ainvnm = zlange( 'I', n, n, ainv, lda, rwork )
418  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
419  rcondi = one
420  ELSE
421  rcondi = ( one / anormi ) / ainvnm
422  END IF
423  nt = 2
424  ELSE
425 *
426 * Do only the condition estimate if INFO > 0.
427 *
428  trfcon = .true.
429  anormo = zlange( 'O', m, n, a, lda, rwork )
430  anormi = zlange( 'I', m, n, a, lda, rwork )
431  rcondo = zero
432  rcondi = zero
433  END IF
434 *
435 * Print information about the tests so far that did not
436 * pass the threshold.
437 *
438  DO 30 k = 1, nt
439  IF( result( k ).GE.thresh ) THEN
440  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
441  $ CALL alahd( nout, path )
442  WRITE( nout, fmt = 9999 )m, n, nb, imat, k,
443  $ result( k )
444  nfail = nfail + 1
445  END IF
446  30 CONTINUE
447  nrun = nrun + nt
448 *
449 * Skip the remaining tests if this is not the first
450 * block size or if M .ne. N. Skip the solve tests if
451 * the matrix is singular.
452 *
453  IF( inb.GT.1 .OR. m.NE.n )
454  $ GO TO 90
455  IF( trfcon )
456  $ GO TO 70
457 *
458  DO 60 irhs = 1, nns
459  nrhs = nsval( irhs )
460  xtype = 'N'
461 *
462  DO 50 itran = 1, ntran
463  trans = transs( itran )
464  IF( itran.EQ.1 ) THEN
465  rcondc = rcondo
466  ELSE
467  rcondc = rcondi
468  END IF
469 *
470 *+ TEST 3
471 * Solve and compute residual for A * X = B.
472 *
473  srnamt = 'ZLARHS'
474  CALL zlarhs( path, xtype, ' ', trans, n, n, kl,
475  $ ku, nrhs, a, lda, xact, lda, b,
476  $ lda, iseed, info )
477  xtype = 'C'
478 *
479  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
480  srnamt = 'ZGETRS'
481  CALL zgetrs( trans, n, nrhs, afac, lda, iwork,
482  $ x, lda, info )
483 *
484 * Check error code from ZGETRS.
485 *
486  IF( info.NE.0 )
487  $ CALL alaerh( path, 'ZGETRS', info, 0, trans,
488  $ n, n, -1, -1, nrhs, imat, nfail,
489  $ nerrs, nout )
490 *
491  CALL zlacpy( 'Full', n, nrhs, b, lda, work,
492  $ lda )
493  CALL zget02( trans, n, n, nrhs, a, lda, x, lda,
494  $ work, lda, rwork, result( 3 ) )
495 *
496 *+ TEST 4
497 * Check solution from generated exact solution.
498 *
499  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
500  $ result( 4 ) )
501 *
502 *+ TESTS 5, 6, and 7
503 * Use iterative refinement to improve the
504 * solution.
505 *
506  srnamt = 'ZGERFS'
507  CALL zgerfs( trans, n, nrhs, a, lda, afac, lda,
508  $ iwork, b, lda, x, lda, rwork,
509  $ rwork( nrhs+1 ), work,
510  $ rwork( 2*nrhs+1 ), info )
511 *
512 * Check error code from ZGERFS.
513 *
514  IF( info.NE.0 )
515  $ CALL alaerh( path, 'ZGERFS', info, 0, trans,
516  $ n, n, -1, -1, nrhs, imat, nfail,
517  $ nerrs, nout )
518 *
519  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
520  $ result( 5 ) )
521  CALL zget07( trans, n, nrhs, a, lda, b, lda, x,
522  $ lda, xact, lda, rwork, .true.,
523  $ rwork( nrhs+1 ), result( 6 ) )
524 *
525 * Print information about the tests that did not
526 * pass the threshold.
527 *
528  DO 40 k = 3, 7
529  IF( result( k ).GE.thresh ) THEN
530  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
531  $ CALL alahd( nout, path )
532  WRITE( nout, fmt = 9998 )trans, n, nrhs,
533  $ imat, k, result( k )
534  nfail = nfail + 1
535  END IF
536  40 CONTINUE
537  nrun = nrun + 5
538  50 CONTINUE
539  60 CONTINUE
540 *
541 *+ TEST 8
542 * Get an estimate of RCOND = 1/CNDNUM.
543 *
544  70 CONTINUE
545  DO 80 itran = 1, 2
546  IF( itran.EQ.1 ) THEN
547  anorm = anormo
548  rcondc = rcondo
549  norm = 'O'
550  ELSE
551  anorm = anormi
552  rcondc = rcondi
553  norm = 'I'
554  END IF
555  srnamt = 'ZGECON'
556  CALL zgecon( norm, n, afac, lda, anorm, rcond,
557  $ work, rwork, info )
558 *
559 * Check error code from ZGECON.
560 *
561  IF( info.NE.0 )
562  $ CALL alaerh( path, 'ZGECON', info, 0, norm, n,
563  $ n, -1, -1, -1, imat, nfail, nerrs,
564  $ nout )
565 *
566 * This line is needed on a Sun SPARCstation.
567 *
568  dummy = rcond
569 *
570  result( 8 ) = dget06( rcond, rcondc )
571 *
572 * Print information about the tests that did not pass
573 * the threshold.
574 *
575  IF( result( 8 ).GE.thresh ) THEN
576  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
577  $ CALL alahd( nout, path )
578  WRITE( nout, fmt = 9997 )norm, n, imat, 8,
579  $ result( 8 )
580  nfail = nfail + 1
581  END IF
582  nrun = nrun + 1
583  80 CONTINUE
584  90 CONTINUE
585  100 CONTINUE
586 *
587  110 CONTINUE
588  120 CONTINUE
589 *
590 * Print a summary of the results.
591 *
592  CALL alasum( path, nout, nfail, nrun, nerrs )
593 *
594  9999 FORMAT( ' M = ', i5, ', N =', i5, ', NB =', i4, ', type ', i2,
595  $ ', test(', i2, ') =', g12.5 )
596  9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
597  $ i2, ', test(', i2, ') =', g12.5 )
598  9997 FORMAT( ' NORM =''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
599  $ ', test(', i2, ') =', g12.5 )
600  RETURN
601 *
602 * End of ZCHKGE
603 *
subroutine zgetrs(TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
ZGETRS
Definition: zgetrs.f:123
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 zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
subroutine zgetrf(M, N, A, LDA, IPIV, INFO)
ZGETRF VARIANT: Crout Level 3 BLAS version of the algorithm.
Definition: zgetrf.f:102
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 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:108
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine zget07(TRANS, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, CHKFERR, BERR, RESLTS)
ZGET07
Definition: zget07.f:168
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zget03(N, A, LDA, AINV, LDAINV, WORK, LDWORK, RWORK, RCOND, RESID)
ZGET03
Definition: zget03.f:112
double precision function zlange(NORM, M, N, A, LDA, WORK)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: zlange.f:117
subroutine zgecon(NORM, N, A, LDA, ANORM, RCOND, WORK, RWORK, INFO)
ZGECON
Definition: zgecon.f:126
subroutine zerrge(PATH, NUNIT)
ZERRGE
Definition: zerrge.f:57
subroutine zgetri(N, A, LDA, IPIV, WORK, LWORK, INFO)
ZGETRI
Definition: zgetri.f:116
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine zgerfs(TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
ZGERFS
Definition: zgerfs.f:188
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
subroutine zget01(M, N, A, LDA, AFAC, LDAFAC, IPIV, RWORK, RESID)
ZGET01
Definition: zget01.f:110
subroutine zget02(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZGET02
Definition: zget02.f:135
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: