LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine zdrvhe_rk ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
complex*16, dimension( * )  A,
complex*16, dimension( * )  AFAC,
complex*16, dimension( * )  E,
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 
)

ZDRVHE_RK

Purpose:
 ZDRVHE_RK tests the driver routines ZHESV_RK.
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]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 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]E
          E is COMPLEX*16 array, dimension (NMAX)

 \param[out] AINV
 \verbatim
          AINV 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]WORK
          WORK is COMPLEX*16 array, dimension (NMAX*max(2,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[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
December 2016

Definition at line 159 of file zdrvhe_rk.f.

159 *
160 * -- LAPACK test routine (version 3.7.0) --
161 * -- LAPACK is a software package provided by Univ. of Tennessee, --
162 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163 * December 2016
164 *
165 * .. Scalar Arguments ..
166  LOGICAL tsterr
167  INTEGER nmax, nn, nout, nrhs
168  DOUBLE PRECISION thresh
169 * ..
170 * .. Array Arguments ..
171  LOGICAL dotype( * )
172  INTEGER iwork( * ), nval( * )
173  DOUBLE PRECISION rwork( * )
174  COMPLEX*16 a( * ), afac( * ), ainv( * ), b( * ), e( * ),
175  $ work( * ), x( * ), xact( * )
176 * ..
177 *
178 * =====================================================================
179 *
180 * .. Parameters ..
181  DOUBLE PRECISION one, zero
182  parameter ( one = 1.0d+0, zero = 0.0d+0 )
183  INTEGER ntypes, ntests
184  parameter ( ntypes = 10, ntests = 3 )
185  INTEGER nfact
186  parameter ( nfact = 2 )
187 * ..
188 * .. Local Scalars ..
189  LOGICAL zerot
190  CHARACTER dist, fact, TYPE, uplo, xtype
191  CHARACTER*3 matpath, path
192  INTEGER i, i1, i2, ifact, imat, in, info, ioff, iuplo,
193  $ izero, j, k, kl, ku, lda, lwork, mode, n,
194  $ nb, nbmin, nerrs, nfail, nimat, nrun, nt
195  DOUBLE PRECISION ainvnm, anorm, cndnum, rcondc
196 * ..
197 * .. Local Arrays ..
198  CHARACTER facts( nfact ), uplos( 2 )
199  INTEGER iseed( 4 ), iseedy( 4 )
200  DOUBLE PRECISION result( ntests )
201 
202 * ..
203 * .. External Functions ..
204  DOUBLE PRECISION zlanhe
205  EXTERNAL zlanhe
206 * ..
207 * .. External Subroutines ..
208  EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx,
211 * ..
212 * .. Scalars in Common ..
213  LOGICAL lerr, ok
214  CHARACTER*32 srnamt
215  INTEGER infot, nunit
216 * ..
217 * .. Common blocks ..
218  COMMON / infoc / infot, nunit, ok, lerr
219  COMMON / srnamc / srnamt
220 * ..
221 * .. Intrinsic Functions ..
222  INTRINSIC max, min
223 * ..
224 * .. Data statements ..
225  DATA iseedy / 1988, 1989, 1990, 1991 /
226  DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
227 * ..
228 * .. Executable Statements ..
229 *
230 * Initialize constants and the random number seed.
231 *
232 * Test path
233 *
234  path( 1: 1 ) = 'Zomplex precision'
235  path( 2: 3 ) = 'HK'
236 *
237 * Path to generate matrices
238 *
239  matpath( 1: 1 ) = 'Zomplex precision'
240  matpath( 2: 3 ) = 'HE'
241 *
242  nrun = 0
243  nfail = 0
244  nerrs = 0
245  DO 10 i = 1, 4
246  iseed( i ) = iseedy( i )
247  10 CONTINUE
248  lwork = max( 2*nmax, nmax*nrhs )
249 *
250 * Test the error exits
251 *
252  IF( tsterr )
253  $ CALL zerrvx( path, nout )
254  infot = 0
255 *
256 * Set the block size and minimum block size for which the block
257 * routine should be used, which will be later returned by ILAENV.
258 *
259  nb = 1
260  nbmin = 2
261  CALL xlaenv( 1, nb )
262  CALL xlaenv( 2, nbmin )
263 *
264 * Do for each value of N in NVAL
265 *
266  DO 180 in = 1, nn
267  n = nval( in )
268  lda = max( n, 1 )
269  xtype = 'N'
270  nimat = ntypes
271  IF( n.LE.0 )
272  $ nimat = 1
273 *
274  DO 170 imat = 1, nimat
275 *
276 * Do the tests only if DOTYPE( IMAT ) is true.
277 *
278  IF( .NOT.dotype( imat ) )
279  $ GO TO 170
280 *
281 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
282 *
283  zerot = imat.GE.3 .AND. imat.LE.6
284  IF( zerot .AND. n.LT.imat-2 )
285  $ GO TO 170
286 *
287 * Do first for UPLO = 'U', then for UPLO = 'L'
288 *
289  DO 160 iuplo = 1, 2
290  uplo = uplos( iuplo )
291 *
292 * Begin generate the test matrix A.
293 *
294 * Set up parameters with ZLATB4 for the matrix generator
295 * based on the type of matrix to be generated.
296 *
297  CALL zlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
298  $ mode, cndnum, dist )
299 *
300 * Generate a matrix with ZLATMS.
301 *
302  srnamt = 'ZLATMS'
303  CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
304  $ cndnum, anorm, kl, ku, uplo, a, lda,
305  $ work, info )
306 *
307 * Check error code from ZLATMS and handle error.
308 *
309  IF( info.NE.0 ) THEN
310  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
311  $ -1, -1, -1, imat, nfail, nerrs, nout )
312  GO TO 160
313  END IF
314 *
315 * For types 3-6, zero one or more rows and columns of
316 * the matrix to test that INFO is returned correctly.
317 *
318  IF( zerot ) THEN
319  IF( imat.EQ.3 ) THEN
320  izero = 1
321  ELSE IF( imat.EQ.4 ) THEN
322  izero = n
323  ELSE
324  izero = n / 2 + 1
325  END IF
326 *
327  IF( imat.LT.6 ) THEN
328 *
329 * Set row and column IZERO to zero.
330 *
331  IF( iuplo.EQ.1 ) THEN
332  ioff = ( izero-1 )*lda
333  DO 20 i = 1, izero - 1
334  a( ioff+i ) = zero
335  20 CONTINUE
336  ioff = ioff + izero
337  DO 30 i = izero, n
338  a( ioff ) = zero
339  ioff = ioff + lda
340  30 CONTINUE
341  ELSE
342  ioff = izero
343  DO 40 i = 1, izero - 1
344  a( ioff ) = zero
345  ioff = ioff + lda
346  40 CONTINUE
347  ioff = ioff - izero
348  DO 50 i = izero, n
349  a( ioff+i ) = zero
350  50 CONTINUE
351  END IF
352  ELSE
353  IF( iuplo.EQ.1 ) THEN
354 *
355 * Set the first IZERO rows and columns to zero.
356 *
357  ioff = 0
358  DO 70 j = 1, n
359  i2 = min( j, izero )
360  DO 60 i = 1, i2
361  a( ioff+i ) = zero
362  60 CONTINUE
363  ioff = ioff + lda
364  70 CONTINUE
365  ELSE
366 *
367 * Set the first IZERO rows and columns to zero.
368 *
369  ioff = 0
370  DO 90 j = 1, n
371  i1 = max( j, izero )
372  DO 80 i = i1, n
373  a( ioff+i ) = zero
374  80 CONTINUE
375  ioff = ioff + lda
376  90 CONTINUE
377  END IF
378  END IF
379  ELSE
380  izero = 0
381  END IF
382 *
383 * End generate the test matrix A.
384 *
385 *
386  DO 150 ifact = 1, nfact
387 *
388 * Do first for FACT = 'F', then for other values.
389 *
390  fact = facts( ifact )
391 *
392 * Compute the condition number
393 *
394  IF( zerot ) THEN
395  IF( ifact.EQ.1 )
396  $ GO TO 150
397  rcondc = zero
398 *
399  ELSE IF( ifact.EQ.1 ) THEN
400 *
401 * Compute the 1-norm of A.
402 *
403  anorm = zlanhe( '1', uplo, n, a, lda, rwork )
404 *
405 * Factor the matrix A.
406 *
407 
408  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
409  CALL zhetrf_rk( uplo, n, afac, lda, e, iwork, work,
410  $ lwork, info )
411 *
412 * Compute inv(A) and take its norm.
413 *
414  CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
415  lwork = (n+nb+1)*(nb+3)
416 *
417 * We need to copute the invesrse to compute
418 * RCONDC that is used later in TEST3.
419 *
420  CALL zhetri_3( uplo, n, ainv, lda, e, iwork,
421  $ work, lwork, info )
422  ainvnm = zlanhe( '1', uplo, n, ainv, lda, rwork )
423 *
424 * Compute the 1-norm condition number of A.
425 *
426  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
427  rcondc = one
428  ELSE
429  rcondc = ( one / anorm ) / ainvnm
430  END IF
431  END IF
432 *
433 * Form an exact solution and set the right hand side.
434 *
435  srnamt = 'ZLARHS'
436  CALL zlarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
437  $ nrhs, a, lda, xact, lda, b, lda, iseed,
438  $ info )
439  xtype = 'C'
440 *
441 * --- Test ZHESV_RK ---
442 *
443  IF( ifact.EQ.2 ) THEN
444  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
445  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
446 *
447 * Factor the matrix and solve the system using
448 * ZHESV_RK.
449 *
450  srnamt = 'ZHESV_RK'
451  CALL zhesv_rk( uplo, n, nrhs, afac, lda, e, iwork,
452  $ x, lda, work, lwork, info )
453 *
454 * Adjust the expected value of INFO to account for
455 * pivoting.
456 *
457  k = izero
458  IF( k.GT.0 ) THEN
459  100 CONTINUE
460  IF( iwork( k ).LT.0 ) THEN
461  IF( iwork( k ).NE.-k ) THEN
462  k = -iwork( k )
463  GO TO 100
464  END IF
465  ELSE IF( iwork( k ).NE.k ) THEN
466  k = iwork( k )
467  GO TO 100
468  END IF
469  END IF
470 *
471 * Check error code from ZHESV_RK and handle error.
472 *
473  IF( info.NE.k ) THEN
474  CALL alaerh( path, 'ZHESV_RK', info, k, uplo,
475  $ n, n, -1, -1, nrhs, imat, nfail,
476  $ nerrs, nout )
477  GO TO 120
478  ELSE IF( info.NE.0 ) THEN
479  GO TO 120
480  END IF
481 *
482 *+ TEST 1 Reconstruct matrix from factors and compute
483 * residual.
484 *
485  CALL zhet01_3( uplo, n, a, lda, afac, lda, e,
486  $ iwork, ainv, lda, rwork,
487  $ result( 1 ) )
488 *
489 *+ TEST 2 Compute residual of the computed solution.
490 *
491  CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
492  CALL zpot02( uplo, n, nrhs, a, lda, x, lda, work,
493  $ lda, rwork, result( 2 ) )
494 *
495 *+ TEST 3
496 * Check solution from generated exact solution.
497 *
498  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
499  $ result( 3 ) )
500  nt = 3
501 *
502 * Print information about the tests that did not pass
503 * the threshold.
504 *
505  DO 110 k = 1, nt
506  IF( result( k ).GE.thresh ) THEN
507  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
508  $ CALL aladhd( nout, path )
509  WRITE( nout, fmt = 9999 )'ZHESV_RK', uplo,
510  $ n, imat, k, result( k )
511  nfail = nfail + 1
512  END IF
513  110 CONTINUE
514  nrun = nrun + nt
515  120 CONTINUE
516  END IF
517 *
518  150 CONTINUE
519 *
520  160 CONTINUE
521  170 CONTINUE
522  180 CONTINUE
523 *
524 * Print a summary of the results.
525 *
526  CALL alasvm( path, nout, nfail, nrun, nerrs )
527 *
528  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
529  $ ', test ', i2, ', ratio =', g12.5 )
530  RETURN
531 *
532 * End of ZDRVHE_RK
533 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine zhetrf_rk(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
ZHETRF_RK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bunch...
Definition: zhetrf_rk.f:261
subroutine zhesv_rk(UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB, WORK, LWORK, INFO)
ZHESV_RK computes the solution to system of linear equations A * X = B for SY matrices ...
Definition: zhesv_rk.f:230
double precision function zlanhe(NORM, UPLO, N, A, LDA, WORK)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix.
Definition: zlanhe.f:126
subroutine zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
subroutine zhet01_3(UPLO, N, A, LDA, AFAC, LDAFAC, E, IPIV, C, LDC, RWORK, RESID)
ZHET01_3
Definition: zhet01_3.f:143
subroutine zpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZPOT02
Definition: zpot02.f:129
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 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 xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zhetri_3(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
ZHETRI_3
Definition: zhetri_3.f:172
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:92
subroutine zerrvx(PATH, NUNIT)
ZERRVX
Definition: zerrvx.f:57
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334

Here is the call graph for this function:

Here is the caller graph for this function: