LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine zdrvsy_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 
)

ZDRVSY_RK

Purpose:
 ZDRVSY_RK tests the driver routines ZSYSV_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 zdrvsy_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 = 11, 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 zlansy
206  EXTERNAL zlansy
207 * ..
208 * .. External Subroutines ..
209  EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx, zget04,
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 ) = 'SK'
237 *
238 * Path to generate matrices
239 *
240  matpath( 1: 1 ) = 'Zomplex precision'
241  matpath( 2: 3 ) = 'SY'
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  IF( imat.NE.ntypes ) THEN
294 *
295 * Begin generate the test matrix A.
296 *
297 * Set up parameters with ZLATB4 for the matrix generator
298 * based on the type of matrix to be generated.
299 *
300  CALL zlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
301  $ mode, cndnum, dist )
302 *
303 * Generate a matrix with ZLATMS.
304 *
305  srnamt = 'ZLATMS'
306  CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
307  $ cndnum, anorm, kl, ku, uplo, a, lda,
308  $ work, info )
309 *
310 * Check error code from DLATMS and handle error.
311 *
312  IF( info.NE.0 ) THEN
313  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
314  $ -1, -1, -1, imat, nfail, nerrs, nout )
315  GO TO 160
316  END IF
317 *
318 * For types 3-6, zero one or more rows and columns of
319 * the matrix to test that INFO is returned correctly.
320 *
321  IF( zerot ) THEN
322  IF( imat.EQ.3 ) THEN
323  izero = 1
324  ELSE IF( imat.EQ.4 ) THEN
325  izero = n
326  ELSE
327  izero = n / 2 + 1
328  END IF
329 *
330  IF( imat.LT.6 ) THEN
331 *
332 * Set row and column IZERO to zero.
333 *
334  IF( iuplo.EQ.1 ) THEN
335  ioff = ( izero-1 )*lda
336  DO 20 i = 1, izero - 1
337  a( ioff+i ) = zero
338  20 CONTINUE
339  ioff = ioff + izero
340  DO 30 i = izero, n
341  a( ioff ) = zero
342  ioff = ioff + lda
343  30 CONTINUE
344  ELSE
345  ioff = izero
346  DO 40 i = 1, izero - 1
347  a( ioff ) = zero
348  ioff = ioff + lda
349  40 CONTINUE
350  ioff = ioff - izero
351  DO 50 i = izero, n
352  a( ioff+i ) = zero
353  50 CONTINUE
354  END IF
355  ELSE
356  IF( iuplo.EQ.1 ) THEN
357 *
358 * Set the first IZERO rows and columns to zero.
359 *
360  ioff = 0
361  DO 70 j = 1, n
362  i2 = min( j, izero )
363  DO 60 i = 1, i2
364  a( ioff+i ) = zero
365  60 CONTINUE
366  ioff = ioff + lda
367  70 CONTINUE
368  ELSE
369 *
370 * Set the first IZERO rows and columns to zero.
371 *
372  ioff = 0
373  DO 90 j = 1, n
374  i1 = max( j, izero )
375  DO 80 i = i1, n
376  a( ioff+i ) = zero
377  80 CONTINUE
378  ioff = ioff + lda
379  90 CONTINUE
380  END IF
381  END IF
382  ELSE
383  izero = 0
384  END IF
385  ELSE
386 *
387 * IMAT = NTYPES: Use a special block diagonal matrix to
388 * test alternate code for the 2-by-2 blocks.
389 *
390  CALL zlatsy( uplo, n, a, lda, iseed )
391  END IF
392 *
393  DO 150 ifact = 1, nfact
394 *
395 * Do first for FACT = 'F', then for other values.
396 *
397  fact = facts( ifact )
398 *
399 * Compute the condition number for comparison with
400 * the value returned by ZSYSVX_ROOK.
401 *
402  IF( zerot ) THEN
403  IF( ifact.EQ.1 )
404  $ GO TO 150
405  rcondc = zero
406 *
407  ELSE IF( ifact.EQ.1 ) THEN
408 *
409 * Compute the 1-norm of A.
410 *
411  anorm = zlansy( '1', uplo, n, a, lda, rwork )
412 *
413 * Factor the matrix A.
414 *
415 
416  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
417  CALL zsytrf_rk( uplo, n, afac, lda, e, iwork, ainv,
418  $ lwork, info )
419 *
420 * Compute inv(A) and take its norm.
421 *
422  CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
423  lwork = (n+nb+1)*(nb+3)
424 *
425 * We need to copute the invesrse to compute
426 * RCONDC that is used later in TEST3.
427 *
428  CALL zsytri_3( uplo, n, ainv, lda, e, iwork,
429  $ work, lwork, info )
430  ainvnm = zlansy( '1', uplo, n, ainv, lda, rwork )
431 *
432 * Compute the 1-norm condition number of A.
433 *
434  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
435  rcondc = one
436  ELSE
437  rcondc = ( one / anorm ) / ainvnm
438  END IF
439  END IF
440 *
441 * Form an exact solution and set the right hand side.
442 *
443  srnamt = 'ZLARHS'
444  CALL zlarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
445  $ nrhs, a, lda, xact, lda, b, lda, iseed,
446  $ info )
447  xtype = 'C'
448 *
449 * --- Test ZSYSV_RK ---
450 *
451  IF( ifact.EQ.2 ) THEN
452  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
453  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
454 *
455 * Factor the matrix and solve the system using
456 * ZSYSV_RK.
457 *
458  srnamt = 'ZSYSV_RK'
459  CALL zsysv_rk( uplo, n, nrhs, afac, lda, e, iwork,
460  $ x, lda, work, lwork, info )
461 *
462 * Adjust the expected value of INFO to account for
463 * pivoting.
464 *
465  k = izero
466  IF( k.GT.0 ) THEN
467  100 CONTINUE
468  IF( iwork( k ).LT.0 ) THEN
469  IF( iwork( k ).NE.-k ) THEN
470  k = -iwork( k )
471  GO TO 100
472  END IF
473  ELSE IF( iwork( k ).NE.k ) THEN
474  k = iwork( k )
475  GO TO 100
476  END IF
477  END IF
478 *
479 * Check error code from ZSYSV_RK and handle error.
480 *
481  IF( info.NE.k ) THEN
482  CALL alaerh( path, 'ZSYSV_RK', info, k, uplo,
483  $ n, n, -1, -1, nrhs, imat, nfail,
484  $ nerrs, nout )
485  GO TO 120
486  ELSE IF( info.NE.0 ) THEN
487  GO TO 120
488  END IF
489 *
490 *+ TEST 1 Reconstruct matrix from factors and compute
491 * residual.
492 *
493  CALL zsyt01_3( uplo, n, a, lda, afac, lda, e,
494  $ iwork, ainv, lda, rwork,
495  $ result( 1 ) )
496 *
497 *+ TEST 2 Compute residual of the computed solution.
498 *
499  CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
500  CALL zsyt02( uplo, n, nrhs, a, lda, x, lda, work,
501  $ lda, rwork, result( 2 ) )
502 *
503 *+ TEST 3
504 * Check solution from generated exact solution.
505 *
506  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
507  $ result( 3 ) )
508  nt = 3
509 *
510 * Print information about the tests that did not pass
511 * the threshold.
512 *
513  DO 110 k = 1, nt
514  IF( result( k ).GE.thresh ) THEN
515  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
516  $ CALL aladhd( nout, path )
517  WRITE( nout, fmt = 9999 )'ZSYSV_RK', uplo,
518  $ n, imat, k, result( k )
519  nfail = nfail + 1
520  END IF
521  110 CONTINUE
522  nrun = nrun + nt
523  120 CONTINUE
524  END IF
525 *
526  150 CONTINUE
527 *
528  160 CONTINUE
529  170 CONTINUE
530  180 CONTINUE
531 *
532 * Print a summary of the results.
533 *
534  CALL alasvm( path, nout, nfail, nrun, nerrs )
535 *
536  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
537  $ ', test ', i2, ', ratio =', g12.5 )
538  RETURN
539 *
540 * End of ZDRVSY_RK
541 *
subroutine zsysv_rk(UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB, WORK, LWORK, INFO)
ZSYSV_RK computes the solution to system of linear equations A * X = B for SY matrices ...
Definition: zsysv_rk.f:230
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 zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
subroutine zsytri_3(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
ZSYTRI_3
Definition: zsytri_3.f:172
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 zsyt01_3(UPLO, N, A, LDA, AFAC, LDAFAC, E, IPIV, C, LDC, RWORK, RESID)
ZSYT01_3
Definition: zsyt01_3.f:143
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine zlatsy(UPLO, N, X, LDX, ISEED)
ZLATSY
Definition: zlatsy.f:91
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zsytrf_rk(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
ZSYTRF_RK computes the factorization of a complex symmetric indefinite matrix using the bounded Bunch...
Definition: zsytrf_rk.f:261
double precision function zlansy(NORM, UPLO, N, A, LDA, WORK)
ZLANSY 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 symmetric matrix.
Definition: zlansy.f:125
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:92
subroutine zsyt02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZSYT02
Definition: zsyt02.f:129
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: