LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine cdrvsy_rk ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
real  THRESH,
logical  TSTERR,
integer  NMAX,
complex, dimension( * )  A,
complex, dimension( * )  AFAC,
complex, dimension( * )  E,
complex, dimension( * )  AINV,
complex, dimension( * )  B,
complex, dimension( * )  X,
complex, dimension( * )  XACT,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

CDRVSY_RK

Purpose:
 CDRVSY_RK tests the driver routines CSYSV_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 REAL
          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 array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is COMPLEX array, dimension (NMAX*NMAX)
[out]E
          E is COMPLEX array, dimension (NMAX)
 \param[out] AINV

 \verbatim
          AINV is COMPLEX array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX array, dimension (NMAX*NRHS)
[out]X
          X is COMPLEX array, dimension (NMAX*NRHS)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NRHS)
[out]WORK
 
[out]RWORK
          RWORK is REAL 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 158 of file cdrvsy_rk.f.

158 *
159 * -- LAPACK test routine (version 3.7.0) --
160 * -- LAPACK is a software package provided by Univ. of Tennessee, --
161 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
162 * December 2016
163 *
164 * .. Scalar Arguments ..
165  LOGICAL tsterr
166  INTEGER nmax, nn, nout, nrhs
167  REAL thresh
168 * ..
169 * .. Array Arguments ..
170  LOGICAL dotype( * )
171  INTEGER iwork( * ), nval( * )
172  REAL rwork( * )
173  COMPLEX a( * ), afac( * ), ainv( * ), b( * ), e( * ),
174  $ work( * ), x( * ), xact( * )
175 * ..
176 *
177 * =====================================================================
178 *
179 * .. Parameters ..
180  REAL one, zero
181  parameter ( one = 1.0e+0, zero = 0.0e+0 )
182  INTEGER ntypes, ntests
183  parameter ( ntypes = 11, ntests = 3 )
184  INTEGER nfact
185  parameter ( nfact = 2 )
186 * ..
187 * .. Local Scalars ..
188  LOGICAL zerot
189  CHARACTER dist, fact, TYPE, uplo, xtype
190  CHARACTER*3 matpath, path
191  INTEGER i, i1, i2, ifact, imat, in, info, ioff, iuplo,
192  $ izero, j, k, kl, ku, lda, lwork, mode, n,
193  $ nb, nbmin, nerrs, nfail, nimat, nrun, nt
194  REAL ainvnm, anorm, cndnum, rcondc
195 * ..
196 * .. Local Arrays ..
197  CHARACTER facts( nfact ), uplos( 2 )
198  INTEGER iseed( 4 ), iseedy( 4 )
199  REAL result( ntests )
200 
201 * ..
202 * .. External Functions ..
203  REAL clansy
204  EXTERNAL clansy
205 * ..
206 * .. External Subroutines ..
207  EXTERNAL aladhd, alaerh, alasvm, xlaenv, cerrvx, cget04,
210 * ..
211 * .. Scalars in Common ..
212  LOGICAL lerr, ok
213  CHARACTER*32 srnamt
214  INTEGER infot, nunit
215 * ..
216 * .. Common blocks ..
217  COMMON / infoc / infot, nunit, ok, lerr
218  COMMON / srnamc / srnamt
219 * ..
220 * .. Intrinsic Functions ..
221  INTRINSIC max, min
222 * ..
223 * .. Data statements ..
224  DATA iseedy / 1988, 1989, 1990, 1991 /
225  DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
226 * ..
227 * .. Executable Statements ..
228 *
229 * Initialize constants and the random number seed.
230 *
231 * Test path
232 *
233  path( 1: 1 ) = 'Complex precision'
234  path( 2: 3 ) = 'SK'
235 *
236 * Path to generate matrices
237 *
238  matpath( 1: 1 ) = 'Complex precision'
239  matpath( 2: 3 ) = 'SY'
240 *
241  nrun = 0
242  nfail = 0
243  nerrs = 0
244  DO 10 i = 1, 4
245  iseed( i ) = iseedy( i )
246  10 CONTINUE
247  lwork = max( 2*nmax, nmax*nrhs )
248 *
249 * Test the error exits
250 *
251  IF( tsterr )
252  $ CALL cerrvx( path, nout )
253  infot = 0
254 *
255 * Set the block size and minimum block size for which the block
256 * routine should be used, which will be later returned by ILAENV.
257 *
258  nb = 1
259  nbmin = 2
260  CALL xlaenv( 1, nb )
261  CALL xlaenv( 2, nbmin )
262 *
263 * Do for each value of N in NVAL
264 *
265  DO 180 in = 1, nn
266  n = nval( in )
267  lda = max( n, 1 )
268  xtype = 'N'
269  nimat = ntypes
270  IF( n.LE.0 )
271  $ nimat = 1
272 *
273  DO 170 imat = 1, nimat
274 *
275 * Do the tests only if DOTYPE( IMAT ) is true.
276 *
277  IF( .NOT.dotype( imat ) )
278  $ GO TO 170
279 *
280 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
281 *
282  zerot = imat.GE.3 .AND. imat.LE.6
283  IF( zerot .AND. n.LT.imat-2 )
284  $ GO TO 170
285 *
286 * Do first for UPLO = 'U', then for UPLO = 'L'
287 *
288  DO 160 iuplo = 1, 2
289  uplo = uplos( iuplo )
290 *
291  IF( imat.NE.ntypes ) THEN
292 *
293 * Begin generate the test matrix A.
294 *
295 * Set up parameters with CLATB4 for the matrix generator
296 * based on the type of matrix to be generated.
297 *
298  CALL clatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
299  $ mode, cndnum, dist )
300 *
301 * Generate a matrix with CLATMS.
302 *
303  srnamt = 'CLATMS'
304  CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
305  $ cndnum, anorm, kl, ku, uplo, a, lda,
306  $ work, info )
307 *
308 * Check error code from CLATMS and handle error.
309 *
310  IF( info.NE.0 ) THEN
311  CALL alaerh( path, 'CLATMS', 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  ELSE
387 *
388 * IMAT = NTYPES: Use a special block diagonal matrix to
389 * test alternate code for the 2-by-2 blocks.
390 *
391  CALL clatsy( uplo, n, a, lda, iseed )
392  END IF
393 *
394  DO 150 ifact = 1, nfact
395 *
396 * Do first for FACT = 'F', then for other values.
397 *
398  fact = facts( ifact )
399 *
400 * Compute the condition number
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 = clansy( '1', uplo, n, a, lda, rwork )
412 *
413 * Factor the matrix A.
414 *
415 
416  CALL clacpy( uplo, n, n, a, lda, afac, lda )
417  CALL csytrf_rk( uplo, n, afac, lda, e, iwork, work,
418  $ lwork, info )
419 *
420 * Compute inv(A) and take its norm.
421 *
422  CALL clacpy( 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 csytri_3( uplo, n, ainv, lda, e, iwork,
429  $ work, lwork, info )
430  ainvnm = clansy( '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 = 'CLARHS'
444  CALL clarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
445  $ nrhs, a, lda, xact, lda, b, lda, iseed,
446  $ info )
447  xtype = 'C'
448 *
449 * --- Test CSYSV_RK ---
450 *
451  IF( ifact.EQ.2 ) THEN
452  CALL clacpy( uplo, n, n, a, lda, afac, lda )
453  CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
454 *
455 * Factor the matrix and solve the system using
456 * CSYSV_RK.
457 *
458  srnamt = 'CSYSV_RK'
459  CALL csysv_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 CSYSV_RK and handle error.
480 *
481  IF( info.NE.k ) THEN
482  CALL alaerh( path, 'CSYSV_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 csyt01_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 clacpy( 'Full', n, nrhs, b, lda, work, lda )
500  CALL csyt02( 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 cget04( 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 )'CSYSV_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 CDRVSY_RK
541 *
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 clatsy(UPLO, N, X, LDX, ISEED)
CLATSY
Definition: clatsy.f:91
subroutine csysv_rk(UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB, WORK, LWORK, INFO)
CSYSV_RK computes the solution to system of linear equations A * X = B for SY matrices ...
Definition: csysv_rk.f:230
subroutine clarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
CLARHS
Definition: clarhs.f:211
subroutine csyt02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
CSYT02
Definition: csyt02.f:129
subroutine cerrvx(PATH, NUNIT)
CERRVX
Definition: cerrvx.f:57
subroutine csytrf_rk(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
CSYTRF_RK computes the factorization of a complex symmetric indefinite matrix using the bounded Bunch...
Definition: csytrf_rk.f:261
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine clatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
CLATMS
Definition: clatms.f:334
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:92
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine csyt01_3(UPLO, N, A, LDA, AFAC, LDAFAC, E, IPIV, C, LDC, RWORK, RESID)
CSYT01_3
Definition: csyt01_3.f:143
subroutine csytri_3(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
CSYTRI_3
Definition: csytri_3.f:172
real function clansy(NORM, UPLO, N, A, LDA, WORK)
CLANSY 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: clansy.f:125
subroutine clatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
CLATB4
Definition: clatb4.f:123
subroutine cget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
CGET04
Definition: cget04.f:104

Here is the call graph for this function:

Here is the caller graph for this function: