LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ zdrvhe_rk()

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)
[out]AINV
          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 160 of file zdrvhe_rk.f.

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