LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ zchkhe_aa()

subroutine zchkhe_aa ( logical, dimension( * )  DOTYPE,
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 
)

ZCHKHE_AA

Purpose:
 ZCHKHE_AA tests ZHETRF_AA, -TRS_AA.
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]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 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 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(NMAX,2*NSMAX))
[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 2017

Definition at line 174 of file zchkhe_aa.f.

174 *
175 * -- LAPACK test routine (version 3.8.0) --
176 * -- LAPACK is a software package provided by Univ. of Tennessee, --
177 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
178 * November 2017
179 *
180  IMPLICIT NONE
181 *
182 * .. Scalar Arguments ..
183  LOGICAL TSTERR
184  INTEGER NMAX, NN, NNB, NNS, NOUT
185  DOUBLE PRECISION THRESH
186 * ..
187 * .. Array Arguments ..
188  LOGICAL DOTYPE( * )
189  INTEGER IWORK( * ), NBVAL( * ), NSVAL( * ), NVAL( * )
190  DOUBLE PRECISION RWORK( * )
191  COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
192  $ WORK( * ), X( * ), XACT( * )
193 * ..
194 *
195 * =====================================================================
196 *
197 * .. Parameters ..
198  DOUBLE PRECISION ZERO
199  parameter( zero = 0.0d+0 )
200  COMPLEX*16 CZERO
201  parameter( czero = ( 0.0d+0, 0.0d+0 ) )
202  INTEGER NTYPES
203  parameter( ntypes = 10 )
204  INTEGER NTESTS
205  parameter( ntests = 9 )
206 * ..
207 * .. Local Scalars ..
208  LOGICAL ZEROT
209  CHARACTER DIST, TYPE, UPLO, XTYPE
210  CHARACTER*3 PATH, MATPATH
211  INTEGER I, I1, I2, IMAT, IN, INB, INFO, IOFF, IRHS,
212  $ IUPLO, IZERO, J, K, KL, KU, LDA, LWORK, MODE,
213  $ N, NB, NERRS, NFAIL, NIMAT, NRHS, NRUN, NT
214  DOUBLE PRECISION ANORM, CNDNUM
215 * ..
216 * .. Local Arrays ..
217  CHARACTER UPLOS( 2 )
218  INTEGER ISEED( 4 ), ISEEDY( 4 )
219  DOUBLE PRECISION RESULT( NTESTS )
220 * ..
221 * .. External Subroutines ..
222  EXTERNAL alaerh, alahd, alasum, xlaenv, zerrhe,
225 * ..
226 * .. Intrinsic Functions ..
227  INTRINSIC max, min
228 * ..
229 * .. Scalars in Common ..
230  LOGICAL LERR, OK
231  CHARACTER*32 SRNAMT
232  INTEGER INFOT, NUNIT
233 * ..
234 * .. Common blocks ..
235  COMMON / infoc / infot, nunit, ok, lerr
236  COMMON / srnamc / srnamt
237 * ..
238 * .. Data statements ..
239  DATA iseedy / 1988, 1989, 1990, 1991 /
240  DATA uplos / 'U', 'L' /
241 * ..
242 * .. Executable Statements ..
243 *
244 * Initialize constants and the random number seed.
245 *
246 * Test path
247 *
248  path( 1: 1 ) = 'Zomplex precision'
249  path( 2: 3 ) = 'HA'
250 *
251 * Path to generate matrices
252 *
253  matpath( 1: 1 ) = 'Zomplex precision'
254  matpath( 2: 3 ) = 'HE'
255  nrun = 0
256  nfail = 0
257  nerrs = 0
258  DO 10 i = 1, 4
259  iseed( i ) = iseedy( i )
260  10 CONTINUE
261 *
262 * Test the error exits
263 *
264  IF( tsterr )
265  $ CALL zerrhe( path, nout )
266  infot = 0
267 *
268 * Set the minimum block size for which the block routine should
269 * be used, which will be later returned by ILAENV
270 *
271  CALL xlaenv( 2, 2 )
272 *
273 * Do for each value of N in NVAL
274 *
275  DO 180 in = 1, nn
276  n = nval( in )
277  IF( n .GT. nmax ) THEN
278  nfail = nfail + 1
279  WRITE(nout, 9995) 'M ', n, nmax
280  GO TO 180
281  END IF
282  lda = max( n, 1 )
283  xtype = 'N'
284  nimat = ntypes
285  IF( n.LE.0 )
286  $ nimat = 1
287 *
288  izero = 0
289  DO 170 imat = 1, nimat
290 *
291 * Do the tests only if DOTYPE( IMAT ) is true.
292 *
293  IF( .NOT.dotype( imat ) )
294  $ GO TO 170
295 *
296 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
297 *
298  zerot = imat.GE.3 .AND. imat.LE.6
299  IF( zerot .AND. n.LT.imat-2 )
300  $ GO TO 170
301 *
302 * Do first for UPLO = 'U', then for UPLO = 'L'
303 *
304  DO 160 iuplo = 1, 2
305  uplo = uplos( iuplo )
306 *
307 * Set up parameters with ZLATB4 for the matrix generator
308 * based on the type of matrix to be generated.
309 *
310  CALL zlatb4( matpath, imat, n, n, TYPE, KL, KU,
311  $ ANORM, MODE, CNDNUM, DIST )
312 *
313 * Generate a matrix with ZLATMS.
314 *
315  srnamt = 'ZLATMS'
316  CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
317  $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
318  $ INFO )
319 *
320 * Check error code from ZLATMS and handle error.
321 *
322  IF( info.NE.0 ) THEN
323  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
324  $ -1, -1, imat, nfail, nerrs, nout )
325 *
326 * Skip all tests for this generated matrix
327 *
328  GO TO 160
329  END IF
330 *
331 * For types 3-6, zero one or more rows and columns of
332 * the matrix to test that INFO is returned correctly.
333 *
334  IF( zerot ) THEN
335  IF( imat.EQ.3 ) THEN
336  izero = 1
337  ELSE IF( imat.EQ.4 ) THEN
338  izero = n
339  ELSE
340  izero = n / 2 + 1
341  END IF
342 *
343  IF( imat.LT.6 ) THEN
344 *
345 * Set row and column IZERO to zero.
346 *
347  IF( iuplo.EQ.1 ) THEN
348  ioff = ( izero-1 )*lda
349  DO 20 i = 1, izero - 1
350  a( ioff+i ) = czero
351  20 CONTINUE
352  ioff = ioff + izero
353  DO 30 i = izero, n
354  a( ioff ) = czero
355  ioff = ioff + lda
356  30 CONTINUE
357  ELSE
358  ioff = izero
359  DO 40 i = 1, izero - 1
360  a( ioff ) = czero
361  ioff = ioff + lda
362  40 CONTINUE
363  ioff = ioff - izero
364  DO 50 i = izero, n
365  a( ioff+i ) = czero
366  50 CONTINUE
367  END IF
368  ELSE
369  IF( iuplo.EQ.1 ) THEN
370 *
371 * Set the first IZERO rows and columns to zero.
372 *
373  ioff = 0
374  DO 70 j = 1, n
375  i2 = min( j, izero )
376  DO 60 i = 1, i2
377  a( ioff+i ) = czero
378  60 CONTINUE
379  ioff = ioff + lda
380  70 CONTINUE
381  izero = 1
382  ELSE
383 *
384 * Set the last IZERO rows and columns to zero.
385 *
386  ioff = 0
387  DO 90 j = 1, n
388  i1 = max( j, izero )
389  DO 80 i = i1, n
390  a( ioff+i ) = czero
391  80 CONTINUE
392  ioff = ioff + lda
393  90 CONTINUE
394  END IF
395  END IF
396  ELSE
397  izero = 0
398  END IF
399 *
400 * End generate test matrix A.
401 *
402 *
403 * Set the imaginary part of the diagonals.
404 *
405  CALL zlaipd( n, a, lda+1, 0 )
406 *
407 * Do for each value of NB in NBVAL
408 *
409  DO 150 inb = 1, nnb
410 *
411 * Set the optimal blocksize, which will be later
412 * returned by ILAENV.
413 *
414  nb = nbval( inb )
415  CALL xlaenv( 1, nb )
416 *
417 * Copy the test matrix A into matrix AFAC which
418 * will be factorized in place. This is needed to
419 * preserve the test matrix A for subsequent tests.
420 *
421  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
422 *
423 * Compute the L*D*L**T or U*D*U**T factorization of the
424 * matrix. IWORK stores details of the interchanges and
425 * the block structure of D. AINV is a work array for
426 * block factorization, LWORK is the length of AINV.
427 *
428  lwork = max( 1, ( nb+1 )*lda )
429  srnamt = 'ZHETRF_AA'
430  CALL zhetrf_aa( uplo, n, afac, lda, iwork, ainv,
431  $ lwork, info )
432 *
433 * Adjust the expected value of INFO to account for
434 * pivoting.
435 *
436 c IF( IZERO.GT.0 ) THEN
437 c J = 1
438 c K = IZERO
439 c 100 CONTINUE
440 c IF( J.EQ.K ) THEN
441 c K = IWORK( J )
442 c ELSE IF( IWORK( J ).EQ.K ) THEN
443 c K = J
444 c END IF
445 c IF( J.LT.K ) THEN
446 c J = J + 1
447 c GO TO 100
448 c END IF
449 c ELSE
450  k = 0
451 c END IF
452 *
453 * Check error code from ZHETRF and handle error.
454 *
455  IF( info.NE.k ) THEN
456  CALL alaerh( path, 'ZHETRF_AA', info, k, uplo,
457  $ n, n, -1, -1, nb, imat, nfail, nerrs,
458  $ nout )
459  END IF
460 *
461 *+ TEST 1
462 * Reconstruct matrix from factors and compute residual.
463 *
464  CALL zhet01_aa( uplo, n, a, lda, afac, lda, iwork,
465  $ ainv, lda, rwork, result( 1 ) )
466  nt = 1
467 *
468 *
469 * Print information about the tests that did not pass
470 * the threshold.
471 *
472  DO 110 k = 1, nt
473  IF( result( k ).GE.thresh ) THEN
474  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
475  $ CALL alahd( nout, path )
476  WRITE( nout, fmt = 9999 )uplo, n, nb, imat, k,
477  $ result( k )
478  nfail = nfail + 1
479  END IF
480  110 CONTINUE
481  nrun = nrun + nt
482 *
483 * Skip solver test if INFO is not 0.
484 *
485  IF( info.NE.0 ) THEN
486  GO TO 140
487  END IF
488 *
489 * Do for each value of NRHS in NSVAL.
490 *
491  DO 130 irhs = 1, nns
492  nrhs = nsval( irhs )
493 *
494 *+ TEST 2 (Using TRS)
495 * Solve and compute residual for A * X = B.
496 *
497 * Choose a set of NRHS random solution vectors
498 * stored in XACT and set up the right hand side B
499 *
500  srnamt = 'ZLARHS'
501  CALL zlarhs( matpath, xtype, uplo, ' ', n, n,
502  $ kl, ku, nrhs, a, lda, xact, lda,
503  $ b, lda, iseed, info )
504  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
505 *
506  srnamt = 'ZHETRS_AA'
507  lwork = max( 1, 3*n-2 )
508  CALL zhetrs_aa( uplo, n, nrhs, afac, lda, iwork,
509  $ x, lda, work, lwork, info )
510 *
511 * Check error code from ZHETRS and handle error.
512 *
513  IF( info.NE.0 ) THEN
514  IF( izero.EQ.0 ) THEN
515  CALL alaerh( path, 'ZHETRS_AA', info, 0,
516  $ uplo, n, n, -1, -1, nrhs, imat,
517  $ nfail, nerrs, nout )
518  END IF
519  ELSE
520 *
521  CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda
522  $ )
523 *
524 * Compute the residual for the solution
525 *
526  CALL zpot02( uplo, n, nrhs, a, lda, x, lda,
527  $ work, lda, rwork, result( 2 ) )
528 *
529 * Print information about the tests that did not pass
530 * the threshold.
531 *
532  DO 120 k = 2, 2
533  IF( result( k ).GE.thresh ) THEN
534  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
535  $ CALL alahd( nout, path )
536  WRITE( nout, fmt = 9998 )uplo, n, nrhs,
537  $ imat, k, result( k )
538  nfail = nfail + 1
539  END IF
540  120 CONTINUE
541  END IF
542  nrun = nrun + 1
543 *
544 * End do for each value of NRHS in NSVAL.
545 *
546  130 CONTINUE
547  140 CONTINUE
548  150 CONTINUE
549  160 CONTINUE
550  170 CONTINUE
551  180 CONTINUE
552 *
553 * Print a summary of the results.
554 *
555  CALL alasum( path, nout, nfail, nrun, nerrs )
556 *
557  9999 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NB =', i4, ', type ',
558  $ i2, ', test ', i2, ', ratio =', g12.5 )
559  9998 FORMAT( ' UPLO = ''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
560  $ i2, ', test(', i2, ') =', g12.5 )
561 c 9997 FORMAT( ' UPLO = ''', A1, ''', N =', I5, ',', 10X, ' type ', I2,
562 c $ ', test(', I2, ') =', G12.5 )
563  9995 FORMAT( ' Invalid input value: ', a4, '=', i6, '; must be <=',
564  $ i6 )
565  RETURN
566 *
567 * End of ZCHKHE_AA
568 *
Here is the call graph for this function:
Here is the caller graph for this function:
zhet01_aa
subroutine zhet01_aa(UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C, LDC, RWORK, RESID)
ZHET01_AA
Definition: zhet01_aa.f:127
zpot02
subroutine zpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZPOT02
Definition: zpot02.f:129
zlarhs
subroutine zlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
ZLARHS
Definition: zlarhs.f:211
alahd
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:109
zlaipd
subroutine zlaipd(N, A, INDA, VINDA)
ZLAIPD
Definition: zlaipd.f:85
zhetrf_aa
subroutine zhetrf_aa(UPLO, N, A, LDA, IPIV, WORK, LWORK, INFO)
ZHETRF_AA
Definition: zhetrf_aa.f:134
zerrhe
subroutine zerrhe(PATH, NUNIT)
ZERRHE
Definition: zerrhe.f:57
zlacpy
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
alaerh
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
alasum
subroutine alasum(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASUM
Definition: alasum.f:75
zhetrs_aa
subroutine zhetrs_aa(UPLO, N, NRHS, A, LDA, IPIV, B, LDB, WORK, LWORK, INFO)
ZHETRS_AA
Definition: zhetrs_aa.f:134
zlatb4
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
xlaenv
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
zlatms
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334