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

DDRVSY_RK

Purpose:
 DDRVSY_RK tests the driver routines DSYSV_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 DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]E
          E is DOUBLE PRECISION array, dimension (NMAX)
[out]AINV
          AINV is DOUBLE PRECISION array, dimension (NMAX*NMAX)
[out]B
          B is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]X
          X is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]XACT
          XACT is DOUBLE PRECISION array, dimension (NMAX*NRHS)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (NMAX*max(2,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[out]IWORK
          IWORK is INTEGER array, dimension (2*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 ddrvsy_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  DOUBLE PRECISION thresh
168 * ..
169 * .. Array Arguments ..
170  LOGICAL dotype( * )
171  INTEGER iwork( * ), nval( * )
172  DOUBLE PRECISION a( * ), afac( * ), ainv( * ), b( * ), e( * ),
173  $ rwork( * ), work( * ), x( * ), xact( * )
174 * ..
175 *
176 * =====================================================================
177 *
178 * .. Parameters ..
179  DOUBLE PRECISION one, zero
180  parameter ( one = 1.0d+0, zero = 0.0d+0 )
181  INTEGER ntypes, ntests
182  parameter ( ntypes = 10, ntests = 3 )
183  INTEGER nfact
184  parameter ( nfact = 2 )
185 * ..
186 * .. Local Scalars ..
187  LOGICAL zerot
188  CHARACTER dist, fact, TYPE, uplo, xtype
189  CHARACTER*3 path, matpath
190  INTEGER i, i1, i2, ifact, imat, in, info, ioff, iuplo,
191  $ izero, j, k, kl, ku, lda, lwork, mode, n,
192  $ nb, nbmin, nerrs, nfail, nimat, nrun, nt
193  DOUBLE PRECISION ainvnm, anorm, cndnum, rcondc
194 * ..
195 * .. Local Arrays ..
196  CHARACTER facts( nfact ), uplos( 2 )
197  INTEGER iseed( 4 ), iseedy( 4 )
198  DOUBLE PRECISION result( ntests )
199 * ..
200 * .. External Functions ..
201  DOUBLE PRECISION dlansy
202  EXTERNAL dlansy
203 * ..
204 * .. External Subroutines ..
205  EXTERNAL aladhd, alaerh, alasvm, derrvx, dget04, dlacpy,
208 * ..
209 * .. Scalars in Common ..
210  LOGICAL lerr, ok
211  CHARACTER*32 srnamt
212  INTEGER infot, nunit
213 * ..
214 * .. Common blocks ..
215  COMMON / infoc / infot, nunit, ok, lerr
216  COMMON / srnamc / srnamt
217 * ..
218 * .. Intrinsic Functions ..
219  INTRINSIC max, min
220 * ..
221 * .. Data statements ..
222  DATA iseedy / 1988, 1989, 1990, 1991 /
223  DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
224 * ..
225 * .. Executable Statements ..
226 *
227 * Initialize constants and the random number seed.
228 *
229 * Test path
230 *
231  path( 1: 1 ) = 'Double precision'
232  path( 2: 3 ) = 'SK'
233 *
234 * Path to generate matrices
235 *
236  matpath( 1: 1 ) = 'Double precision'
237  matpath( 2: 3 ) = 'SY'
238 *
239  nrun = 0
240  nfail = 0
241  nerrs = 0
242  DO 10 i = 1, 4
243  iseed( i ) = iseedy( i )
244  10 CONTINUE
245  lwork = max( 2*nmax, nmax*nrhs )
246 *
247 * Test the error exits
248 *
249  IF( tsterr )
250  $ CALL derrvx( path, nout )
251  infot = 0
252 *
253 * Set the block size and minimum block size for which the block
254 * routine should be used, which will be later returned by ILAENV.
255 *
256  nb = 1
257  nbmin = 2
258  CALL xlaenv( 1, nb )
259  CALL xlaenv( 2, nbmin )
260 *
261 * Do for each value of N in NVAL
262 *
263  DO 180 in = 1, nn
264  n = nval( in )
265  lda = max( n, 1 )
266  xtype = 'N'
267  nimat = ntypes
268  IF( n.LE.0 )
269  $ nimat = 1
270 *
271  DO 170 imat = 1, nimat
272 *
273 * Do the tests only if DOTYPE( IMAT ) is true.
274 *
275  IF( .NOT.dotype( imat ) )
276  $ GO TO 170
277 *
278 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
279 *
280  zerot = imat.GE.3 .AND. imat.LE.6
281  IF( zerot .AND. n.LT.imat-2 )
282  $ GO TO 170
283 *
284 * Do first for UPLO = 'U', then for UPLO = 'L'
285 *
286  DO 160 iuplo = 1, 2
287  uplo = uplos( iuplo )
288 *
289 * Begin generate the test matrix A.
290 *
291 * Set up parameters with DLATB4 for the matrix generator
292 * based on the type of matrix to be generated.
293 *
294  CALL dlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
295  $ mode, cndnum, dist )
296 *
297 * Generate a matrix with DLATMS.
298 *
299  srnamt = 'DLATMS'
300  CALL dlatms( n, n, dist, iseed, TYPE, rwork, mode,
301  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
302  $ info )
303 *
304 * Check error code from DLATMS and handle error.
305 *
306  IF( info.NE.0 ) THEN
307  CALL alaerh( path, 'DLATMS', info, 0, uplo, n, n, -1,
308  $ -1, -1, imat, nfail, nerrs, nout )
309 *
310 * Skip all tests for this generated matrix
311 *
312  GO TO 160
313  END IF
314 *
315 * For types 3-6, zero one or more rows and columns of the
316 * 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  ioff = 0
354  IF( iuplo.EQ.1 ) THEN
355 *
356 * Set the first IZERO rows and columns to zero.
357 *
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 last IZERO rows and columns to zero.
368 *
369  DO 90 j = 1, n
370  i1 = max( j, izero )
371  DO 80 i = i1, n
372  a( ioff+i ) = zero
373  80 CONTINUE
374  ioff = ioff + lda
375  90 CONTINUE
376  END IF
377  END IF
378  ELSE
379  izero = 0
380  END IF
381 *
382 * End generate the test matrix A.
383 *
384  DO 150 ifact = 1, nfact
385 *
386 * Do first for FACT = 'F', then for other values.
387 *
388  fact = facts( ifact )
389 *
390 * Compute the condition number
391 *
392  IF( zerot ) THEN
393  IF( ifact.EQ.1 )
394  $ GO TO 150
395  rcondc = zero
396 *
397  ELSE IF( ifact.EQ.1 ) THEN
398 *
399 * Compute the 1-norm of A.
400 *
401  anorm = dlansy( '1', uplo, n, a, lda, rwork )
402 *
403 * Factor the matrix A.
404 *
405  CALL dlacpy( uplo, n, n, a, lda, afac, lda )
406  CALL dsytrf_rk( uplo, n, afac, lda, e, iwork, work,
407  $ lwork, info )
408 *
409 * Compute inv(A) and take its norm.
410 *
411  CALL dlacpy( uplo, n, n, afac, lda, ainv, lda )
412  lwork = (n+nb+1)*(nb+3)
413 *
414 * We need to copute the invesrse to compute
415 * RCONDC that is used later in TEST3.
416 *
417  CALL dsytri_3( uplo, n, ainv, lda, e, iwork,
418  $ work, lwork, info )
419  ainvnm = dlansy( '1', uplo, n, ainv, lda, rwork )
420 *
421 * Compute the 1-norm condition number of A.
422 *
423  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
424  rcondc = one
425  ELSE
426  rcondc = ( one / anorm ) / ainvnm
427  END IF
428  END IF
429 *
430 * Form an exact solution and set the right hand side.
431 *
432  srnamt = 'DLARHS'
433  CALL dlarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
434  $ nrhs, a, lda, xact, lda, b, lda, iseed,
435  $ info )
436  xtype = 'C'
437 *
438 * --- Test DSYSV_RK ---
439 *
440  IF( ifact.EQ.2 ) THEN
441  CALL dlacpy( uplo, n, n, a, lda, afac, lda )
442  CALL dlacpy( 'Full', n, nrhs, b, lda, x, lda )
443 *
444 * Factor the matrix and solve the system using
445 * DSYSV_RK.
446 *
447  srnamt = 'DSYSV_RK'
448  CALL dsysv_rk( uplo, n, nrhs, afac, lda, e, iwork,
449  $ x, lda, work, lwork, info )
450 *
451 * Adjust the expected value of INFO to account for
452 * pivoting.
453 *
454  k = izero
455  IF( k.GT.0 ) THEN
456  100 CONTINUE
457  IF( iwork( k ).LT.0 ) THEN
458  IF( iwork( k ).NE.-k ) THEN
459  k = -iwork( k )
460  GO TO 100
461  END IF
462  ELSE IF( iwork( k ).NE.k ) THEN
463  k = iwork( k )
464  GO TO 100
465  END IF
466  END IF
467 *
468 * Check error code from DSYSV_RK and handle error.
469 *
470  IF( info.NE.k ) THEN
471  CALL alaerh( path, 'DSYSV_RK', info, k, uplo,
472  $ n, n, -1, -1, nrhs, imat, nfail,
473  $ nerrs, nout )
474  GO TO 120
475  ELSE IF( info.NE.0 ) THEN
476  GO TO 120
477  END IF
478 *
479 *+ TEST 1 Reconstruct matrix from factors and compute
480 * residual.
481 *
482  CALL dsyt01_3( uplo, n, a, lda, afac, lda, e,
483  $ iwork, ainv, lda, rwork,
484  $ result( 1 ) )
485 *
486 *+ TEST 2 Compute residual of the computed solution.
487 *
488  CALL dlacpy( 'Full', n, nrhs, b, lda, work, lda )
489  CALL dpot02( uplo, n, nrhs, a, lda, x, lda, work,
490  $ lda, rwork, result( 2 ) )
491 *
492 *+ TEST 3
493 * Check solution from generated exact solution.
494 *
495  CALL dget04( n, nrhs, x, lda, xact, lda, rcondc,
496  $ result( 3 ) )
497  nt = 3
498 *
499 * Print information about the tests that did not pass
500 * the threshold.
501 *
502  DO 110 k = 1, nt
503  IF( result( k ).GE.thresh ) THEN
504  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
505  $ CALL aladhd( nout, path )
506  WRITE( nout, fmt = 9999 )'DSYSV_RK', uplo,
507  $ n, imat, k, result( k )
508  nfail = nfail + 1
509  END IF
510  110 CONTINUE
511  nrun = nrun + nt
512  120 CONTINUE
513  END IF
514 *
515  150 CONTINUE
516 *
517  160 CONTINUE
518  170 CONTINUE
519  180 CONTINUE
520 *
521 * Print a summary of the results.
522 *
523  CALL alasvm( path, nout, nfail, nrun, nerrs )
524 *
525  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
526  $ ', test ', i2, ', ratio =', g12.5 )
527  RETURN
528 *
529 * End of DDRVSY_RK
530 *
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric matrix.
Definition: dlansy.f:124
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine dlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
DLARHS
Definition: dlarhs.f:206
subroutine dsytri_3(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
DSYTRI_3
Definition: dsytri_3.f:172
subroutine dsyt01_3(UPLO, N, A, LDA, AFAC, LDAFAC, E, IPIV, C, LDC, RWORK, RESID)
DSYT01_3
Definition: dsyt01_3.f:142
subroutine dsytrf_rk(UPLO, N, A, LDA, E, IPIV, WORK, LWORK, INFO)
DSYTRF_RK computes the factorization of a real symmetric indefinite matrix using the bounded Bunch-Ka...
Definition: dsytrf_rk.f:261
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine dlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
DLATB4
Definition: dlatb4.f:122
subroutine dget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
DGET04
Definition: dget04.f:104
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:92
subroutine derrvx(PATH, NUNIT)
DERRVX
Definition: derrvx.f:57
subroutine dsysv_rk(UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB, WORK, LWORK, INFO)
DSYSV_RK computes the solution to system of linear equations A * X = B for SY matrices ...
Definition: dsysv_rk.f:230
subroutine dlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
DLATMS
Definition: dlatms.f:323
subroutine dpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
DPOT02
Definition: dpot02.f:129

Here is the call graph for this function:

Here is the caller graph for this function: