LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sdrvgt.f
Go to the documentation of this file.
1*> \brief \b SDRVGT
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE SDRVGT( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, AF,
12* B, X, XACT, WORK, RWORK, IWORK, NOUT )
13*
14* .. Scalar Arguments ..
15* LOGICAL TSTERR
16* INTEGER NN, NOUT, NRHS
17* REAL THRESH
18* ..
19* .. Array Arguments ..
20* LOGICAL DOTYPE( * )
21* INTEGER IWORK( * ), NVAL( * )
22* REAL A( * ), AF( * ), B( * ), RWORK( * ), WORK( * ),
23* $ X( * ), XACT( * )
24* ..
25*
26*
27*> \par Purpose:
28* =============
29*>
30*> \verbatim
31*>
32*> SDRVGT tests SGTSV and -SVX.
33*> \endverbatim
34*
35* Arguments:
36* ==========
37*
38*> \param[in] DOTYPE
39*> \verbatim
40*> DOTYPE is LOGICAL array, dimension (NTYPES)
41*> The matrix types to be used for testing. Matrices of type j
42*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
43*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
44*> \endverbatim
45*>
46*> \param[in] NN
47*> \verbatim
48*> NN is INTEGER
49*> The number of values of N contained in the vector NVAL.
50*> \endverbatim
51*>
52*> \param[in] NVAL
53*> \verbatim
54*> NVAL is INTEGER array, dimension (NN)
55*> The values of the matrix dimension N.
56*> \endverbatim
57*>
58*> \param[in] NRHS
59*> \verbatim
60*> NRHS is INTEGER
61*> The number of right hand sides, NRHS >= 0.
62*> \endverbatim
63*>
64*> \param[in] THRESH
65*> \verbatim
66*> THRESH is REAL
67*> The threshold value for the test ratios. A result is
68*> included in the output file if RESULT >= THRESH. To have
69*> every test ratio printed, use THRESH = 0.
70*> \endverbatim
71*>
72*> \param[in] TSTERR
73*> \verbatim
74*> TSTERR is LOGICAL
75*> Flag that indicates whether error exits are to be tested.
76*> \endverbatim
77*>
78*> \param[out] A
79*> \verbatim
80*> A is REAL array, dimension (NMAX*4)
81*> \endverbatim
82*>
83*> \param[out] AF
84*> \verbatim
85*> AF is REAL array, dimension (NMAX*4)
86*> \endverbatim
87*>
88*> \param[out] B
89*> \verbatim
90*> B is REAL array, dimension (NMAX*NRHS)
91*> \endverbatim
92*>
93*> \param[out] X
94*> \verbatim
95*> X is REAL array, dimension (NMAX*NRHS)
96*> \endverbatim
97*>
98*> \param[out] XACT
99*> \verbatim
100*> XACT is REAL array, dimension (NMAX*NRHS)
101*> \endverbatim
102*>
103*> \param[out] WORK
104*> \verbatim
105*> WORK is REAL array, dimension
106*> (NMAX*max(3,NRHS))
107*> \endverbatim
108*>
109*> \param[out] RWORK
110*> \verbatim
111*> RWORK is REAL array, dimension
112*> (max(NMAX,2*NRHS))
113*> \endverbatim
114*>
115*> \param[out] IWORK
116*> \verbatim
117*> IWORK is INTEGER array, dimension (2*NMAX)
118*> \endverbatim
119*>
120*> \param[in] NOUT
121*> \verbatim
122*> NOUT is INTEGER
123*> The unit number for output.
124*> \endverbatim
125*
126* Authors:
127* ========
128*
129*> \author Univ. of Tennessee
130*> \author Univ. of California Berkeley
131*> \author Univ. of Colorado Denver
132*> \author NAG Ltd.
133*
134*> \ingroup single_lin
135*
136* =====================================================================
137 SUBROUTINE sdrvgt( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, A, AF,
138 $ B, X, XACT, WORK, RWORK, IWORK, NOUT )
139*
140* -- LAPACK test routine --
141* -- LAPACK is a software package provided by Univ. of Tennessee, --
142* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
143*
144* .. Scalar Arguments ..
145 LOGICAL TSTERR
146 INTEGER NN, NOUT, NRHS
147 REAL THRESH
148* ..
149* .. Array Arguments ..
150 LOGICAL DOTYPE( * )
151 INTEGER IWORK( * ), NVAL( * )
152 REAL A( * ), AF( * ), B( * ), RWORK( * ), WORK( * ),
153 $ x( * ), xact( * )
154* ..
155*
156* =====================================================================
157*
158* .. Parameters ..
159 REAL ONE, ZERO
160 parameter( one = 1.0e+0, zero = 0.0e+0 )
161 INTEGER NTYPES
162 parameter( ntypes = 12 )
163 INTEGER NTESTS
164 parameter( ntests = 6 )
165* ..
166* .. Local Scalars ..
167 LOGICAL TRFCON, ZEROT
168 CHARACTER DIST, FACT, TRANS, TYPE
169 CHARACTER*3 PATH
170 INTEGER I, IFACT, IMAT, IN, INFO, ITRAN, IX, IZERO, J,
171 $ k, k1, kl, koff, ku, lda, m, mode, n, nerrs,
172 $ nfail, nimat, nrun, nt
173 REAL AINVNM, ANORM, ANORMI, ANORMO, COND, RCOND,
174 $ rcondc, rcondi, rcondo
175* ..
176* .. Local Arrays ..
177 CHARACTER TRANSS( 3 )
178 INTEGER ISEED( 4 ), ISEEDY( 4 )
179 REAL RESULT( NTESTS ), Z( 3 )
180* ..
181* .. External Functions ..
182 REAL SASUM, SGET06, SLANGT
183 EXTERNAL sasum, sget06, slangt
184* ..
185* .. External Subroutines ..
186 EXTERNAL aladhd, alaerh, alasvm, scopy, serrvx, sget04,
189 $ slatms, sscal
190* ..
191* .. Intrinsic Functions ..
192 INTRINSIC max
193* ..
194* .. Scalars in Common ..
195 LOGICAL LERR, OK
196 CHARACTER*32 SRNAMT
197 INTEGER INFOT, NUNIT
198* ..
199* .. Common blocks ..
200 COMMON / infoc / infot, nunit, ok, lerr
201 COMMON / srnamc / srnamt
202* ..
203* .. Data statements ..
204 DATA iseedy / 0, 0, 0, 1 / , transs / 'N', 'T',
205 $ 'C' /
206* ..
207* .. Executable Statements ..
208*
209 path( 1: 1 ) = 'Single precision'
210 path( 2: 3 ) = 'GT'
211 nrun = 0
212 nfail = 0
213 nerrs = 0
214 DO 10 i = 1, 4
215 iseed( i ) = iseedy( i )
216 10 CONTINUE
217*
218* Test the error exits
219*
220 IF( tsterr )
221 $ CALL serrvx( path, nout )
222 infot = 0
223*
224 DO 140 in = 1, nn
225*
226* Do for each value of N in NVAL.
227*
228 n = nval( in )
229 m = max( n-1, 0 )
230 lda = max( 1, n )
231 nimat = ntypes
232 IF( n.LE.0 )
233 $ nimat = 1
234*
235 DO 130 imat = 1, nimat
236*
237* Do the tests only if DOTYPE( IMAT ) is true.
238*
239 IF( .NOT.dotype( imat ) )
240 $ GO TO 130
241*
242* Set up parameters with SLATB4.
243*
244 CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
245 $ cond, dist )
246*
247 zerot = imat.GE.8 .AND. imat.LE.10
248 IF( imat.LE.6 ) THEN
249*
250* Types 1-6: generate matrices of known condition number.
251*
252 koff = max( 2-ku, 3-max( 1, n ) )
253 srnamt = 'SLATMS'
254 CALL slatms( n, n, dist, iseed, TYPE, rwork, mode, cond,
255 $ anorm, kl, ku, 'Z', af( koff ), 3, work,
256 $ info )
257*
258* Check the error code from SLATMS.
259*
260 IF( info.NE.0 ) THEN
261 CALL alaerh( path, 'SLATMS', info, 0, ' ', n, n, kl,
262 $ ku, -1, imat, nfail, nerrs, nout )
263 GO TO 130
264 END IF
265 izero = 0
266*
267 IF( n.GT.1 ) THEN
268 CALL scopy( n-1, af( 4 ), 3, a, 1 )
269 CALL scopy( n-1, af( 3 ), 3, a( n+m+1 ), 1 )
270 END IF
271 CALL scopy( n, af( 2 ), 3, a( m+1 ), 1 )
272 ELSE
273*
274* Types 7-12: generate tridiagonal matrices with
275* unknown condition numbers.
276*
277 IF( .NOT.zerot .OR. .NOT.dotype( 7 ) ) THEN
278*
279* Generate a matrix with elements from [-1,1].
280*
281 CALL slarnv( 2, iseed, n+2*m, a )
282 IF( anorm.NE.one )
283 $ CALL sscal( n+2*m, anorm, a, 1 )
284 ELSE IF( izero.GT.0 ) THEN
285*
286* Reuse the last matrix by copying back the zeroed out
287* elements.
288*
289 IF( izero.EQ.1 ) THEN
290 a( n ) = z( 2 )
291 IF( n.GT.1 )
292 $ a( 1 ) = z( 3 )
293 ELSE IF( izero.EQ.n ) THEN
294 a( 3*n-2 ) = z( 1 )
295 a( 2*n-1 ) = z( 2 )
296 ELSE
297 a( 2*n-2+izero ) = z( 1 )
298 a( n-1+izero ) = z( 2 )
299 a( izero ) = z( 3 )
300 END IF
301 END IF
302*
303* If IMAT > 7, set one column of the matrix to 0.
304*
305 IF( .NOT.zerot ) THEN
306 izero = 0
307 ELSE IF( imat.EQ.8 ) THEN
308 izero = 1
309 z( 2 ) = a( n )
310 a( n ) = zero
311 IF( n.GT.1 ) THEN
312 z( 3 ) = a( 1 )
313 a( 1 ) = zero
314 END IF
315 ELSE IF( imat.EQ.9 ) THEN
316 izero = n
317 z( 1 ) = a( 3*n-2 )
318 z( 2 ) = a( 2*n-1 )
319 a( 3*n-2 ) = zero
320 a( 2*n-1 ) = zero
321 ELSE
322 izero = ( n+1 ) / 2
323 DO 20 i = izero, n - 1
324 a( 2*n-2+i ) = zero
325 a( n-1+i ) = zero
326 a( i ) = zero
327 20 CONTINUE
328 a( 3*n-2 ) = zero
329 a( 2*n-1 ) = zero
330 END IF
331 END IF
332*
333 DO 120 ifact = 1, 2
334 IF( ifact.EQ.1 ) THEN
335 fact = 'F'
336 ELSE
337 fact = 'N'
338 END IF
339*
340* Compute the condition number for comparison with
341* the value returned by SGTSVX.
342*
343 IF( zerot ) THEN
344 IF( ifact.EQ.1 )
345 $ GO TO 120
346 rcondo = zero
347 rcondi = zero
348*
349 ELSE IF( ifact.EQ.1 ) THEN
350 CALL scopy( n+2*m, a, 1, af, 1 )
351*
352* Compute the 1-norm and infinity-norm of A.
353*
354 anormo = slangt( '1', n, a, a( m+1 ), a( n+m+1 ) )
355 anormi = slangt( 'I', n, a, a( m+1 ), a( n+m+1 ) )
356*
357* Factor the matrix A.
358*
359 CALL sgttrf( n, af, af( m+1 ), af( n+m+1 ),
360 $ af( n+2*m+1 ), iwork, info )
361*
362* Use SGTTRS to solve for one column at a time of
363* inv(A), computing the maximum column sum as we go.
364*
365 ainvnm = zero
366 DO 40 i = 1, n
367 DO 30 j = 1, n
368 x( j ) = zero
369 30 CONTINUE
370 x( i ) = one
371 CALL sgttrs( 'No transpose', n, 1, af, af( m+1 ),
372 $ af( n+m+1 ), af( n+2*m+1 ), iwork, x,
373 $ lda, info )
374 ainvnm = max( ainvnm, sasum( n, x, 1 ) )
375 40 CONTINUE
376*
377* Compute the 1-norm condition number of A.
378*
379 IF( anormo.LE.zero .OR. ainvnm.LE.zero ) THEN
380 rcondo = one
381 ELSE
382 rcondo = ( one / anormo ) / ainvnm
383 END IF
384*
385* Use SGTTRS to solve for one column at a time of
386* inv(A'), computing the maximum column sum as we go.
387*
388 ainvnm = zero
389 DO 60 i = 1, n
390 DO 50 j = 1, n
391 x( j ) = zero
392 50 CONTINUE
393 x( i ) = one
394 CALL sgttrs( 'Transpose', n, 1, af, af( m+1 ),
395 $ af( n+m+1 ), af( n+2*m+1 ), iwork, x,
396 $ lda, info )
397 ainvnm = max( ainvnm, sasum( n, x, 1 ) )
398 60 CONTINUE
399*
400* Compute the infinity-norm condition number of A.
401*
402 IF( anormi.LE.zero .OR. ainvnm.LE.zero ) THEN
403 rcondi = one
404 ELSE
405 rcondi = ( one / anormi ) / ainvnm
406 END IF
407 END IF
408*
409 DO 110 itran = 1, 3
410 trans = transs( itran )
411 IF( itran.EQ.1 ) THEN
412 rcondc = rcondo
413 ELSE
414 rcondc = rcondi
415 END IF
416*
417* Generate NRHS random solution vectors.
418*
419 ix = 1
420 DO 70 j = 1, nrhs
421 CALL slarnv( 2, iseed, n, xact( ix ) )
422 ix = ix + lda
423 70 CONTINUE
424*
425* Set the right hand side.
426*
427 CALL slagtm( trans, n, nrhs, one, a, a( m+1 ),
428 $ a( n+m+1 ), xact, lda, zero, b, lda )
429*
430 IF( ifact.EQ.2 .AND. itran.EQ.1 ) THEN
431*
432* --- Test SGTSV ---
433*
434* Solve the system using Gaussian elimination with
435* partial pivoting.
436*
437 CALL scopy( n+2*m, a, 1, af, 1 )
438 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
439*
440 srnamt = 'SGTSV '
441 CALL sgtsv( n, nrhs, af, af( m+1 ), af( n+m+1 ), x,
442 $ lda, info )
443*
444* Check error code from SGTSV .
445*
446 IF( info.NE.izero )
447 $ CALL alaerh( path, 'SGTSV ', info, izero, ' ',
448 $ n, n, 1, 1, nrhs, imat, nfail,
449 $ nerrs, nout )
450 nt = 1
451 IF( izero.EQ.0 ) THEN
452*
453* Check residual of computed solution.
454*
455 CALL slacpy( 'Full', n, nrhs, b, lda, work,
456 $ lda )
457 CALL sgtt02( trans, n, nrhs, a, a( m+1 ),
458 $ a( n+m+1 ), x, lda, work, lda,
459 $ result( 2 ) )
460*
461* Check solution from generated exact solution.
462*
463 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
464 $ result( 3 ) )
465 nt = 3
466 END IF
467*
468* Print information about the tests that did not pass
469* the threshold.
470*
471 DO 80 k = 2, nt
472 IF( result( k ).GE.thresh ) THEN
473 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
474 $ CALL aladhd( nout, path )
475 WRITE( nout, fmt = 9999 )'SGTSV ', n, imat,
476 $ k, result( k )
477 nfail = nfail + 1
478 END IF
479 80 CONTINUE
480 nrun = nrun + nt - 1
481 END IF
482*
483* --- Test SGTSVX ---
484*
485 IF( ifact.GT.1 ) THEN
486*
487* Initialize AF to zero.
488*
489 DO 90 i = 1, 3*n - 2
490 af( i ) = zero
491 90 CONTINUE
492 END IF
493 CALL slaset( 'Full', n, nrhs, zero, zero, x, lda )
494*
495* Solve the system and compute the condition number and
496* error bounds using SGTSVX.
497*
498 srnamt = 'SGTSVX'
499 CALL sgtsvx( fact, trans, n, nrhs, a, a( m+1 ),
500 $ a( n+m+1 ), af, af( m+1 ), af( n+m+1 ),
501 $ af( n+2*m+1 ), iwork, b, lda, x, lda,
502 $ rcond, rwork, rwork( nrhs+1 ), work,
503 $ iwork( n+1 ), info )
504*
505* Check the error code from SGTSVX.
506*
507 IF( info.NE.izero )
508 $ CALL alaerh( path, 'SGTSVX', info, izero,
509 $ fact // trans, n, n, 1, 1, nrhs, imat,
510 $ nfail, nerrs, nout )
511*
512 IF( ifact.GE.2 ) THEN
513*
514* Reconstruct matrix from factors and compute
515* residual.
516*
517 CALL sgtt01( n, a, a( m+1 ), a( n+m+1 ), af,
518 $ af( m+1 ), af( n+m+1 ), af( n+2*m+1 ),
519 $ iwork, work, lda, rwork, result( 1 ) )
520 k1 = 1
521 ELSE
522 k1 = 2
523 END IF
524*
525 IF( info.EQ.0 ) THEN
526 trfcon = .false.
527*
528* Check residual of computed solution.
529*
530 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
531 CALL sgtt02( trans, n, nrhs, a, a( m+1 ),
532 $ a( n+m+1 ), x, lda, work, lda,
533 $ result( 2 ) )
534*
535* Check solution from generated exact solution.
536*
537 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
538 $ result( 3 ) )
539*
540* Check the error bounds from iterative refinement.
541*
542 CALL sgtt05( trans, n, nrhs, a, a( m+1 ),
543 $ a( n+m+1 ), b, lda, x, lda, xact, lda,
544 $ rwork, rwork( nrhs+1 ), result( 4 ) )
545 nt = 5
546 END IF
547*
548* Print information about the tests that did not pass
549* the threshold.
550*
551 DO 100 k = k1, nt
552 IF( result( k ).GE.thresh ) THEN
553 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
554 $ CALL aladhd( nout, path )
555 WRITE( nout, fmt = 9998 )'SGTSVX', fact, trans,
556 $ n, imat, k, result( k )
557 nfail = nfail + 1
558 END IF
559 100 CONTINUE
560*
561* Check the reciprocal of the condition number.
562*
563 result( 6 ) = sget06( rcond, rcondc )
564 IF( result( 6 ).GE.thresh ) THEN
565 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
566 $ CALL aladhd( nout, path )
567 WRITE( nout, fmt = 9998 )'SGTSVX', fact, trans, n,
568 $ imat, k, result( k )
569 nfail = nfail + 1
570 END IF
571 nrun = nrun + nt - k1 + 2
572*
573 110 CONTINUE
574 120 CONTINUE
575 130 CONTINUE
576 140 CONTINUE
577*
578* Print a summary of the results.
579*
580 CALL alasvm( path, nout, nfail, nrun, nerrs )
581*
582 9999 FORMAT( 1x, a, ', N =', i5, ', type ', i2, ', test ', i2,
583 $ ', ratio = ', g12.5 )
584 9998 FORMAT( 1x, a, ', FACT=''', a1, ''', TRANS=''', a1, ''', N =',
585 $ i5, ', type ', i2, ', test ', i2, ', ratio = ', g12.5 )
586 RETURN
587*
588* End of SDRVGT
589*
590 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine aladhd(iounit, path)
ALADHD
Definition aladhd.f:90
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
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:127
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:293
subroutine sgttrf(n, dl, d, du, du2, ipiv, info)
SGTTRF
Definition sgttrf.f:124
subroutine sgttrs(trans, n, nrhs, dl, d, du, du2, ipiv, b, ldb, info)
SGTTRS
Definition sgttrs.f:138
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
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:145
subroutine slarnv(idist, iseed, n, x)
SLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition slarnv.f:97
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:110
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sdrvgt(dotype, nn, nval, nrhs, thresh, tsterr, a, af, b, x, xact, work, rwork, iwork, nout)
SDRVGT
Definition sdrvgt.f:139
subroutine serrvx(path, nunit)
SERRVX
Definition serrvx.f:55
subroutine sget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
SGET04
Definition sget04.f:102
subroutine sgtt01(n, dl, d, du, dlf, df, duf, du2, ipiv, work, ldwork, rwork, resid)
SGTT01
Definition sgtt01.f:134
subroutine sgtt02(trans, n, nrhs, dl, d, du, x, ldx, b, ldb, resid)
SGTT02
Definition sgtt02.f:125
subroutine sgtt05(trans, n, nrhs, dl, d, du, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
SGTT05
Definition sgtt05.f:165
subroutine slatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
SLATB4
Definition slatb4.f:120
subroutine slatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
SLATMS
Definition slatms.f:321