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

SDRVGT

Purpose:
 SDRVGT tests SGTSV and -SVX.
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 sides, NRHS >= 0.
[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.
[out]A
          A is REAL array, dimension (NMAX*4)
[out]AF
          AF is REAL array, dimension (NMAX*4)
[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(3,NRHS))
[out]RWORK
          RWORK is REAL array, dimension
                      (max(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 2011

Definition at line 141 of file sdrvgt.f.

141 *
142 * -- LAPACK test routine (version 3.4.0) --
143 * -- LAPACK is a software package provided by Univ. of Tennessee, --
144 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
145 * November 2011
146 *
147 * .. Scalar Arguments ..
148  LOGICAL tsterr
149  INTEGER nn, nout, nrhs
150  REAL thresh
151 * ..
152 * .. Array Arguments ..
153  LOGICAL dotype( * )
154  INTEGER iwork( * ), nval( * )
155  REAL a( * ), af( * ), b( * ), rwork( * ), work( * ),
156  $ x( * ), xact( * )
157 * ..
158 *
159 * =====================================================================
160 *
161 * .. Parameters ..
162  REAL one, zero
163  parameter ( one = 1.0e+0, zero = 0.0e+0 )
164  INTEGER ntypes
165  parameter ( ntypes = 12 )
166  INTEGER ntests
167  parameter ( ntests = 6 )
168 * ..
169 * .. Local Scalars ..
170  LOGICAL trfcon, zerot
171  CHARACTER dist, fact, trans, type
172  CHARACTER*3 path
173  INTEGER i, ifact, imat, in, info, itran, ix, izero, j,
174  $ k, k1, kl, koff, ku, lda, m, mode, n, nerrs,
175  $ nfail, nimat, nrun, nt
176  REAL ainvnm, anorm, anormi, anormo, cond, rcond,
177  $ rcondc, rcondi, rcondo
178 * ..
179 * .. Local Arrays ..
180  CHARACTER transs( 3 )
181  INTEGER iseed( 4 ), iseedy( 4 )
182  REAL result( ntests ), z( 3 )
183 * ..
184 * .. External Functions ..
185  REAL sasum, sget06, slangt
186  EXTERNAL sasum, sget06, slangt
187 * ..
188 * .. External Subroutines ..
189  EXTERNAL aladhd, alaerh, alasvm, scopy, serrvx, sget04,
192  $ slatms, sscal
193 * ..
194 * .. Intrinsic Functions ..
195  INTRINSIC max
196 * ..
197 * .. Scalars in Common ..
198  LOGICAL lerr, ok
199  CHARACTER*32 srnamt
200  INTEGER infot, nunit
201 * ..
202 * .. Common blocks ..
203  COMMON / infoc / infot, nunit, ok, lerr
204  COMMON / srnamc / srnamt
205 * ..
206 * .. Data statements ..
207  DATA iseedy / 0, 0, 0, 1 / , transs / 'N', 'T',
208  $ 'C' /
209 * ..
210 * .. Executable Statements ..
211 *
212  path( 1: 1 ) = 'Single precision'
213  path( 2: 3 ) = 'GT'
214  nrun = 0
215  nfail = 0
216  nerrs = 0
217  DO 10 i = 1, 4
218  iseed( i ) = iseedy( i )
219  10 CONTINUE
220 *
221 * Test the error exits
222 *
223  IF( tsterr )
224  $ CALL serrvx( path, nout )
225  infot = 0
226 *
227  DO 140 in = 1, nn
228 *
229 * Do for each value of N in NVAL.
230 *
231  n = nval( in )
232  m = max( n-1, 0 )
233  lda = max( 1, n )
234  nimat = ntypes
235  IF( n.LE.0 )
236  $ nimat = 1
237 *
238  DO 130 imat = 1, nimat
239 *
240 * Do the tests only if DOTYPE( IMAT ) is true.
241 *
242  IF( .NOT.dotype( imat ) )
243  $ GO TO 130
244 *
245 * Set up parameters with SLATB4.
246 *
247  CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
248  $ cond, dist )
249 *
250  zerot = imat.GE.8 .AND. imat.LE.10
251  IF( imat.LE.6 ) THEN
252 *
253 * Types 1-6: generate matrices of known condition number.
254 *
255  koff = max( 2-ku, 3-max( 1, n ) )
256  srnamt = 'SLATMS'
257  CALL slatms( n, n, dist, iseed, TYPE, rwork, mode, cond,
258  $ anorm, kl, ku, 'Z', af( koff ), 3, work,
259  $ info )
260 *
261 * Check the error code from SLATMS.
262 *
263  IF( info.NE.0 ) THEN
264  CALL alaerh( path, 'SLATMS', info, 0, ' ', n, n, kl,
265  $ ku, -1, imat, nfail, nerrs, nout )
266  GO TO 130
267  END IF
268  izero = 0
269 *
270  IF( n.GT.1 ) THEN
271  CALL scopy( n-1, af( 4 ), 3, a, 1 )
272  CALL scopy( n-1, af( 3 ), 3, a( n+m+1 ), 1 )
273  END IF
274  CALL scopy( n, af( 2 ), 3, a( m+1 ), 1 )
275  ELSE
276 *
277 * Types 7-12: generate tridiagonal matrices with
278 * unknown condition numbers.
279 *
280  IF( .NOT.zerot .OR. .NOT.dotype( 7 ) ) THEN
281 *
282 * Generate a matrix with elements from [-1,1].
283 *
284  CALL slarnv( 2, iseed, n+2*m, a )
285  IF( anorm.NE.one )
286  $ CALL sscal( n+2*m, anorm, a, 1 )
287  ELSE IF( izero.GT.0 ) THEN
288 *
289 * Reuse the last matrix by copying back the zeroed out
290 * elements.
291 *
292  IF( izero.EQ.1 ) THEN
293  a( n ) = z( 2 )
294  IF( n.GT.1 )
295  $ a( 1 ) = z( 3 )
296  ELSE IF( izero.EQ.n ) THEN
297  a( 3*n-2 ) = z( 1 )
298  a( 2*n-1 ) = z( 2 )
299  ELSE
300  a( 2*n-2+izero ) = z( 1 )
301  a( n-1+izero ) = z( 2 )
302  a( izero ) = z( 3 )
303  END IF
304  END IF
305 *
306 * If IMAT > 7, set one column of the matrix to 0.
307 *
308  IF( .NOT.zerot ) THEN
309  izero = 0
310  ELSE IF( imat.EQ.8 ) THEN
311  izero = 1
312  z( 2 ) = a( n )
313  a( n ) = zero
314  IF( n.GT.1 ) THEN
315  z( 3 ) = a( 1 )
316  a( 1 ) = zero
317  END IF
318  ELSE IF( imat.EQ.9 ) THEN
319  izero = n
320  z( 1 ) = a( 3*n-2 )
321  z( 2 ) = a( 2*n-1 )
322  a( 3*n-2 ) = zero
323  a( 2*n-1 ) = zero
324  ELSE
325  izero = ( n+1 ) / 2
326  DO 20 i = izero, n - 1
327  a( 2*n-2+i ) = zero
328  a( n-1+i ) = zero
329  a( i ) = zero
330  20 CONTINUE
331  a( 3*n-2 ) = zero
332  a( 2*n-1 ) = zero
333  END IF
334  END IF
335 *
336  DO 120 ifact = 1, 2
337  IF( ifact.EQ.1 ) THEN
338  fact = 'F'
339  ELSE
340  fact = 'N'
341  END IF
342 *
343 * Compute the condition number for comparison with
344 * the value returned by SGTSVX.
345 *
346  IF( zerot ) THEN
347  IF( ifact.EQ.1 )
348  $ GO TO 120
349  rcondo = zero
350  rcondi = zero
351 *
352  ELSE IF( ifact.EQ.1 ) THEN
353  CALL scopy( n+2*m, a, 1, af, 1 )
354 *
355 * Compute the 1-norm and infinity-norm of A.
356 *
357  anormo = slangt( '1', n, a, a( m+1 ), a( n+m+1 ) )
358  anormi = slangt( 'I', n, a, a( m+1 ), a( n+m+1 ) )
359 *
360 * Factor the matrix A.
361 *
362  CALL sgttrf( n, af, af( m+1 ), af( n+m+1 ),
363  $ af( n+2*m+1 ), iwork, info )
364 *
365 * Use SGTTRS to solve for one column at a time of
366 * inv(A), computing the maximum column sum as we go.
367 *
368  ainvnm = zero
369  DO 40 i = 1, n
370  DO 30 j = 1, n
371  x( j ) = zero
372  30 CONTINUE
373  x( i ) = one
374  CALL sgttrs( 'No transpose', n, 1, af, af( m+1 ),
375  $ af( n+m+1 ), af( n+2*m+1 ), iwork, x,
376  $ lda, info )
377  ainvnm = max( ainvnm, sasum( n, x, 1 ) )
378  40 CONTINUE
379 *
380 * Compute the 1-norm condition number of A.
381 *
382  IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
383  rcondo = one
384  ELSE
385  rcondo = ( one / anormo ) / ainvnm
386  END IF
387 *
388 * Use SGTTRS to solve for one column at a time of
389 * inv(A'), computing the maximum column sum as we go.
390 *
391  ainvnm = zero
392  DO 60 i = 1, n
393  DO 50 j = 1, n
394  x( j ) = zero
395  50 CONTINUE
396  x( i ) = one
397  CALL sgttrs( 'Transpose', n, 1, af, af( m+1 ),
398  $ af( n+m+1 ), af( n+2*m+1 ), iwork, x,
399  $ lda, info )
400  ainvnm = max( ainvnm, sasum( n, x, 1 ) )
401  60 CONTINUE
402 *
403 * Compute the infinity-norm condition number of A.
404 *
405  IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
406  rcondi = one
407  ELSE
408  rcondi = ( one / anormi ) / ainvnm
409  END IF
410  END IF
411 *
412  DO 110 itran = 1, 3
413  trans = transs( itran )
414  IF( itran.EQ.1 ) THEN
415  rcondc = rcondo
416  ELSE
417  rcondc = rcondi
418  END IF
419 *
420 * Generate NRHS random solution vectors.
421 *
422  ix = 1
423  DO 70 j = 1, nrhs
424  CALL slarnv( 2, iseed, n, xact( ix ) )
425  ix = ix + lda
426  70 CONTINUE
427 *
428 * Set the right hand side.
429 *
430  CALL slagtm( trans, n, nrhs, one, a, a( m+1 ),
431  $ a( n+m+1 ), xact, lda, zero, b, lda )
432 *
433  IF( ifact.EQ.2 .AND. itran.EQ.1 ) THEN
434 *
435 * --- Test SGTSV ---
436 *
437 * Solve the system using Gaussian elimination with
438 * partial pivoting.
439 *
440  CALL scopy( n+2*m, a, 1, af, 1 )
441  CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
442 *
443  srnamt = 'SGTSV '
444  CALL sgtsv( n, nrhs, af, af( m+1 ), af( n+m+1 ), x,
445  $ lda, info )
446 *
447 * Check error code from SGTSV .
448 *
449  IF( info.NE.izero )
450  $ CALL alaerh( path, 'SGTSV ', info, izero, ' ',
451  $ n, n, 1, 1, nrhs, imat, nfail,
452  $ nerrs, nout )
453  nt = 1
454  IF( izero.EQ.0 ) THEN
455 *
456 * Check residual of computed solution.
457 *
458  CALL slacpy( 'Full', n, nrhs, b, lda, work,
459  $ lda )
460  CALL sgtt02( trans, n, nrhs, a, a( m+1 ),
461  $ a( n+m+1 ), x, lda, work, lda,
462  $ result( 2 ) )
463 *
464 * Check solution from generated exact solution.
465 *
466  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
467  $ result( 3 ) )
468  nt = 3
469  END IF
470 *
471 * Print information about the tests that did not pass
472 * the threshold.
473 *
474  DO 80 k = 2, nt
475  IF( result( k ).GE.thresh ) THEN
476  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
477  $ CALL aladhd( nout, path )
478  WRITE( nout, fmt = 9999 )'SGTSV ', n, imat,
479  $ k, result( k )
480  nfail = nfail + 1
481  END IF
482  80 CONTINUE
483  nrun = nrun + nt - 1
484  END IF
485 *
486 * --- Test SGTSVX ---
487 *
488  IF( ifact.GT.1 ) THEN
489 *
490 * Initialize AF to zero.
491 *
492  DO 90 i = 1, 3*n - 2
493  af( i ) = zero
494  90 CONTINUE
495  END IF
496  CALL slaset( 'Full', n, nrhs, zero, zero, x, lda )
497 *
498 * Solve the system and compute the condition number and
499 * error bounds using SGTSVX.
500 *
501  srnamt = 'SGTSVX'
502  CALL sgtsvx( fact, trans, n, nrhs, a, a( m+1 ),
503  $ a( n+m+1 ), af, af( m+1 ), af( n+m+1 ),
504  $ af( n+2*m+1 ), iwork, b, lda, x, lda,
505  $ rcond, rwork, rwork( nrhs+1 ), work,
506  $ iwork( n+1 ), info )
507 *
508 * Check the error code from SGTSVX.
509 *
510  IF( info.NE.izero )
511  $ CALL alaerh( path, 'SGTSVX', info, izero,
512  $ fact // trans, n, n, 1, 1, nrhs, imat,
513  $ nfail, nerrs, nout )
514 *
515  IF( ifact.GE.2 ) THEN
516 *
517 * Reconstruct matrix from factors and compute
518 * residual.
519 *
520  CALL sgtt01( n, a, a( m+1 ), a( n+m+1 ), af,
521  $ af( m+1 ), af( n+m+1 ), af( n+2*m+1 ),
522  $ iwork, work, lda, rwork, result( 1 ) )
523  k1 = 1
524  ELSE
525  k1 = 2
526  END IF
527 *
528  IF( info.EQ.0 ) THEN
529  trfcon = .false.
530 *
531 * Check residual of computed solution.
532 *
533  CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
534  CALL sgtt02( trans, n, nrhs, a, a( m+1 ),
535  $ a( n+m+1 ), x, lda, work, lda,
536  $ result( 2 ) )
537 *
538 * Check solution from generated exact solution.
539 *
540  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
541  $ result( 3 ) )
542 *
543 * Check the error bounds from iterative refinement.
544 *
545  CALL sgtt05( trans, n, nrhs, a, a( m+1 ),
546  $ a( n+m+1 ), b, lda, x, lda, xact, lda,
547  $ rwork, rwork( nrhs+1 ), result( 4 ) )
548  nt = 5
549  END IF
550 *
551 * Print information about the tests that did not pass
552 * the threshold.
553 *
554  DO 100 k = k1, nt
555  IF( result( k ).GE.thresh ) THEN
556  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
557  $ CALL aladhd( nout, path )
558  WRITE( nout, fmt = 9998 )'SGTSVX', fact, trans,
559  $ n, imat, k, result( k )
560  nfail = nfail + 1
561  END IF
562  100 CONTINUE
563 *
564 * Check the reciprocal of the condition number.
565 *
566  result( 6 ) = sget06( rcond, rcondc )
567  IF( result( 6 ).GE.thresh ) THEN
568  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
569  $ CALL aladhd( nout, path )
570  WRITE( nout, fmt = 9998 )'SGTSVX', fact, trans, n,
571  $ imat, k, result( k )
572  nfail = nfail + 1
573  END IF
574  nrun = nrun + nt - k1 + 2
575 *
576  110 CONTINUE
577  120 CONTINUE
578  130 CONTINUE
579  140 CONTINUE
580 *
581 * Print a summary of the results.
582 *
583  CALL alasvm( path, nout, nfail, nrun, nerrs )
584 *
585  9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test ', i2,
586  $ ', ratio = ', g12.5 )
587  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N =',
588  $ i5, ', type ', i2, ', test ', i2, ', ratio = ', g12.5 )
589  RETURN
590 *
591 * End of SDRVGT
592 *
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
real function slangt(NORM, N, DL, D, DU)
SLANGT returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: slangt.f:108
subroutine sgtt02(TRANS, N, NRHS, DL, D, DU, X, LDX, B, LDB, RESID)
SGTT02
Definition: sgtt02.f:126
subroutine slarnv(IDIST, ISEED, N, X)
SLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: slarnv.f:99
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
subroutine sgttrf(N, DL, D, DU, DU2, IPIV, INFO)
SGTTRF
Definition: sgttrf.f:126
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 slagtm(TRANS, N, NRHS, ALPHA, DL, D, DU, X, LDX, BETA, B, LDB)
SLAGTM performs a matrix-matrix product of the form C = αAB+βC, where A is a tridiagonal matrix...
Definition: slagtm.f:147
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 sgtt01(N, DL, D, DU, DLF, DF, DUF, DU2, IPIV, WORK, LDWORK, RWORK, RESID)
SGTT01
Definition: sgtt01.f:136
subroutine sget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
SGET04
Definition: sget04.f:104
real function sasum(N, SX, INCX)
SASUM
Definition: sasum.f:54
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine sgtt05(TRANS, N, NRHS, DL, D, DU, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
SGTT05
Definition: sgtt05.f:167
subroutine sgttrs(TRANS, N, NRHS, DL, D, DU, DU2, IPIV, B, LDB, INFO)
SGTTRS
Definition: sgttrs.f:140
subroutine serrvx(PATH, NUNIT)
SERRVX
Definition: serrvx.f:57
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53
subroutine sgtsvx(FACT, TRANS, N, NRHS, DL, D, DU, DLF, DF, DUF, DU2, IPIV, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
SGTSVX computes the solution to system of linear equations A * X = B for GT matrices ...
Definition: sgtsvx.f:295
subroutine sgtsv(N, NRHS, DL, D, DU, B, LDB, INFO)
SGTSV computes the solution to system of linear equations A * X = B for GT matrices ...
Definition: sgtsv.f:129

Here is the call graph for this function:

Here is the caller graph for this function: