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

SDRVSY_ROOK

Purpose:
 SDRVSY_ROOK tests the driver routines SSYSV_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 REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]AINV
          AINV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NRHS)
[out]X
          X is REAL array, dimension (NMAX*NRHS)
[out]XACT
          XACT is REAL array, dimension (NMAX*NRHS)
[out]WORK
          WORK is REAL array, dimension
                      (NMAX*max(2,NRHS))
[out]RWORK
          RWORK is REAL 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 155 of file sdrvsy_rook.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: