LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine ddrvsy_rook ( 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( * )  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_ROOK

Purpose:
 DDRVSY_ROOK tests the driver routines DSYSV_ROOK.
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]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
November 2013

Definition at line 154 of file ddrvsy_rook.f.

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