LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ schkgt()

subroutine schkgt ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nns,
integer, dimension( * )  nsval,
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 
)

SCHKGT

Purpose:
 SCHKGT tests SGTTRF, -TRS, -RFS, and -CON
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]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[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*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is REAL array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is REAL array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is REAL array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL array, dimension
                      (max(NMAX,2*NSMAX))
[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.

Definition at line 144 of file schkgt.f.

146*
147* -- LAPACK test routine --
148* -- LAPACK is a software package provided by Univ. of Tennessee, --
149* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
150*
151* .. Scalar Arguments ..
152 LOGICAL TSTERR
153 INTEGER NN, NNS, NOUT
154 REAL THRESH
155* ..
156* .. Array Arguments ..
157 LOGICAL DOTYPE( * )
158 INTEGER IWORK( * ), NSVAL( * ), NVAL( * )
159 REAL A( * ), AF( * ), B( * ), RWORK( * ), WORK( * ),
160 $ X( * ), XACT( * )
161* ..
162*
163* =====================================================================
164*
165* .. Parameters ..
166 REAL ONE, ZERO
167 parameter( one = 1.0e+0, zero = 0.0e+0 )
168 INTEGER NTYPES
169 parameter( ntypes = 12 )
170 INTEGER NTESTS
171 parameter( ntests = 7 )
172* ..
173* .. Local Scalars ..
174 LOGICAL TRFCON, ZEROT
175 CHARACTER DIST, NORM, TRANS, TYPE
176 CHARACTER*3 PATH
177 INTEGER I, IMAT, IN, INFO, IRHS, ITRAN, IX, IZERO, J,
178 $ K, KL, KOFF, KU, LDA, M, MODE, N, NERRS, NFAIL,
179 $ NIMAT, NRHS, NRUN
180 REAL AINVNM, ANORM, COND, RCOND, RCONDC, RCONDI,
181 $ RCONDO
182* ..
183* .. Local Arrays ..
184 CHARACTER TRANSS( 3 )
185 INTEGER ISEED( 4 ), ISEEDY( 4 )
186 REAL RESULT( NTESTS ), Z( 3 )
187* ..
188* .. External Functions ..
189 REAL SASUM, SGET06, SLANGT
190 EXTERNAL sasum, sget06, slangt
191* ..
192* .. External Subroutines ..
193 EXTERNAL alaerh, alahd, alasum, scopy, serrge, sget04,
196 $ sscal
197* ..
198* .. Intrinsic Functions ..
199 INTRINSIC max
200* ..
201* .. Scalars in Common ..
202 LOGICAL LERR, OK
203 CHARACTER*32 SRNAMT
204 INTEGER INFOT, NUNIT
205* ..
206* .. Common blocks ..
207 COMMON / infoc / infot, nunit, ok, lerr
208 COMMON / srnamc / srnamt
209* ..
210* .. Data statements ..
211 DATA iseedy / 0, 0, 0, 1 / , transs / 'N', 'T',
212 $ 'C' /
213* ..
214* .. Executable Statements ..
215*
216 path( 1: 1 ) = 'Single precision'
217 path( 2: 3 ) = 'GT'
218 nrun = 0
219 nfail = 0
220 nerrs = 0
221 DO 10 i = 1, 4
222 iseed( i ) = iseedy( i )
223 10 CONTINUE
224*
225* Test the error exits
226*
227 IF( tsterr )
228 $ CALL serrge( path, nout )
229 infot = 0
230*
231 DO 110 in = 1, nn
232*
233* Do for each value of N in NVAL.
234*
235 n = nval( in )
236 m = max( n-1, 0 )
237 lda = max( 1, n )
238 nimat = ntypes
239 IF( n.LE.0 )
240 $ nimat = 1
241*
242 DO 100 imat = 1, nimat
243*
244* Do the tests only if DOTYPE( IMAT ) is true.
245*
246 IF( .NOT.dotype( imat ) )
247 $ GO TO 100
248*
249* Set up parameters with SLATB4.
250*
251 CALL slatb4( path, imat, n, n, TYPE, KL, KU, ANORM, MODE,
252 $ COND, DIST )
253*
254 zerot = imat.GE.8 .AND. imat.LE.10
255 IF( imat.LE.6 ) THEN
256*
257* Types 1-6: generate matrices of known condition number.
258*
259 koff = max( 2-ku, 3-max( 1, n ) )
260 srnamt = 'SLATMS'
261 CALL slatms( n, n, dist, iseed, TYPE, RWORK, MODE, COND,
262 $ ANORM, KL, KU, 'Z', AF( KOFF ), 3, WORK,
263 $ INFO )
264*
265* Check the error code from SLATMS.
266*
267 IF( info.NE.0 ) THEN
268 CALL alaerh( path, 'SLATMS', info, 0, ' ', n, n, kl,
269 $ ku, -1, imat, nfail, nerrs, nout )
270 GO TO 100
271 END IF
272 izero = 0
273*
274 IF( n.GT.1 ) THEN
275 CALL scopy( n-1, af( 4 ), 3, a, 1 )
276 CALL scopy( n-1, af( 3 ), 3, a( n+m+1 ), 1 )
277 END IF
278 CALL scopy( n, af( 2 ), 3, a( m+1 ), 1 )
279 ELSE
280*
281* Types 7-12: generate tridiagonal matrices with
282* unknown condition numbers.
283*
284 IF( .NOT.zerot .OR. .NOT.dotype( 7 ) ) THEN
285*
286* Generate a matrix with elements from [-1,1].
287*
288 CALL slarnv( 2, iseed, n+2*m, a )
289 IF( anorm.NE.one )
290 $ CALL sscal( n+2*m, anorm, a, 1 )
291 ELSE IF( izero.GT.0 ) THEN
292*
293* Reuse the last matrix by copying back the zeroed out
294* elements.
295*
296 IF( izero.EQ.1 ) THEN
297 a( n ) = z( 2 )
298 IF( n.GT.1 )
299 $ a( 1 ) = z( 3 )
300 ELSE IF( izero.EQ.n ) THEN
301 a( 3*n-2 ) = z( 1 )
302 a( 2*n-1 ) = z( 2 )
303 ELSE
304 a( 2*n-2+izero ) = z( 1 )
305 a( n-1+izero ) = z( 2 )
306 a( izero ) = z( 3 )
307 END IF
308 END IF
309*
310* If IMAT > 7, set one column of the matrix to 0.
311*
312 IF( .NOT.zerot ) THEN
313 izero = 0
314 ELSE IF( imat.EQ.8 ) THEN
315 izero = 1
316 z( 2 ) = a( n )
317 a( n ) = zero
318 IF( n.GT.1 ) THEN
319 z( 3 ) = a( 1 )
320 a( 1 ) = zero
321 END IF
322 ELSE IF( imat.EQ.9 ) THEN
323 izero = n
324 z( 1 ) = a( 3*n-2 )
325 z( 2 ) = a( 2*n-1 )
326 a( 3*n-2 ) = zero
327 a( 2*n-1 ) = zero
328 ELSE
329 izero = ( n+1 ) / 2
330 DO 20 i = izero, n - 1
331 a( 2*n-2+i ) = zero
332 a( n-1+i ) = zero
333 a( i ) = zero
334 20 CONTINUE
335 a( 3*n-2 ) = zero
336 a( 2*n-1 ) = zero
337 END IF
338 END IF
339*
340*+ TEST 1
341* Factor A as L*U and compute the ratio
342* norm(L*U - A) / (n * norm(A) * EPS )
343*
344 CALL scopy( n+2*m, a, 1, af, 1 )
345 srnamt = 'SGTTRF'
346 CALL sgttrf( n, af, af( m+1 ), af( n+m+1 ), af( n+2*m+1 ),
347 $ iwork, info )
348*
349* Check error code from SGTTRF.
350*
351 IF( info.NE.izero )
352 $ CALL alaerh( path, 'SGTTRF', info, izero, ' ', n, n, 1,
353 $ 1, -1, imat, nfail, nerrs, nout )
354 trfcon = info.NE.0
355*
356 CALL sgtt01( n, a, a( m+1 ), a( n+m+1 ), af, af( m+1 ),
357 $ af( n+m+1 ), af( n+2*m+1 ), iwork, work, lda,
358 $ rwork, result( 1 ) )
359*
360* Print the test ratio if it is .GE. THRESH.
361*
362 IF( result( 1 ).GE.thresh ) THEN
363 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
364 $ CALL alahd( nout, path )
365 WRITE( nout, fmt = 9999 )n, imat, 1, result( 1 )
366 nfail = nfail + 1
367 END IF
368 nrun = nrun + 1
369*
370 DO 50 itran = 1, 2
371 trans = transs( itran )
372 IF( itran.EQ.1 ) THEN
373 norm = 'O'
374 ELSE
375 norm = 'I'
376 END IF
377 anorm = slangt( norm, n, a, a( m+1 ), a( n+m+1 ) )
378*
379 IF( .NOT.trfcon ) THEN
380*
381* Use SGTTRS to solve for one column at a time of inv(A)
382* or inv(A^T), computing the maximum column sum as we
383* go.
384*
385 ainvnm = zero
386 DO 40 i = 1, n
387 DO 30 j = 1, n
388 x( j ) = zero
389 30 CONTINUE
390 x( i ) = one
391 CALL sgttrs( trans, n, 1, af, af( m+1 ),
392 $ af( n+m+1 ), af( n+2*m+1 ), iwork, x,
393 $ lda, info )
394 ainvnm = max( ainvnm, sasum( n, x, 1 ) )
395 40 CONTINUE
396*
397* Compute RCONDC = 1 / (norm(A) * norm(inv(A))
398*
399 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
400 rcondc = one
401 ELSE
402 rcondc = ( one / anorm ) / ainvnm
403 END IF
404 IF( itran.EQ.1 ) THEN
405 rcondo = rcondc
406 ELSE
407 rcondi = rcondc
408 END IF
409 ELSE
410 rcondc = zero
411 END IF
412*
413*+ TEST 7
414* Estimate the reciprocal of the condition number of the
415* matrix.
416*
417 srnamt = 'SGTCON'
418 CALL sgtcon( norm, n, af, af( m+1 ), af( n+m+1 ),
419 $ af( n+2*m+1 ), iwork, anorm, rcond, work,
420 $ iwork( n+1 ), info )
421*
422* Check error code from SGTCON.
423*
424 IF( info.NE.0 )
425 $ CALL alaerh( path, 'SGTCON', info, 0, norm, n, n, -1,
426 $ -1, -1, imat, nfail, nerrs, nout )
427*
428 result( 7 ) = sget06( rcond, rcondc )
429*
430* Print the test ratio if it is .GE. THRESH.
431*
432 IF( result( 7 ).GE.thresh ) THEN
433 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
434 $ CALL alahd( nout, path )
435 WRITE( nout, fmt = 9997 )norm, n, imat, 7,
436 $ result( 7 )
437 nfail = nfail + 1
438 END IF
439 nrun = nrun + 1
440 50 CONTINUE
441*
442* Skip the remaining tests if the matrix is singular.
443*
444 IF( trfcon )
445 $ GO TO 100
446*
447 DO 90 irhs = 1, nns
448 nrhs = nsval( irhs )
449*
450* Generate NRHS random solution vectors.
451*
452 ix = 1
453 DO 60 j = 1, nrhs
454 CALL slarnv( 2, iseed, n, xact( ix ) )
455 ix = ix + lda
456 60 CONTINUE
457*
458 DO 80 itran = 1, 3
459 trans = transs( itran )
460 IF( itran.EQ.1 ) THEN
461 rcondc = rcondo
462 ELSE
463 rcondc = rcondi
464 END IF
465*
466* Set the right hand side.
467*
468 CALL slagtm( trans, n, nrhs, one, a, a( m+1 ),
469 $ a( n+m+1 ), xact, lda, zero, b, lda )
470*
471*+ TEST 2
472* Solve op(A) * X = B and compute the residual.
473*
474 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
475 srnamt = 'SGTTRS'
476 CALL sgttrs( trans, n, nrhs, af, af( m+1 ),
477 $ af( n+m+1 ), af( n+2*m+1 ), iwork, x,
478 $ lda, info )
479*
480* Check error code from SGTTRS.
481*
482 IF( info.NE.0 )
483 $ CALL alaerh( path, 'SGTTRS', info, 0, trans, n, n,
484 $ -1, -1, nrhs, imat, nfail, nerrs,
485 $ nout )
486*
487 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
488 CALL sgtt02( trans, n, nrhs, a, a( m+1 ), a( n+m+1 ),
489 $ x, lda, work, lda, result( 2 ) )
490*
491*+ TEST 3
492* Check solution from generated exact solution.
493*
494 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
495 $ result( 3 ) )
496*
497*+ TESTS 4, 5, and 6
498* Use iterative refinement to improve the solution.
499*
500 srnamt = 'SGTRFS'
501 CALL sgtrfs( trans, n, nrhs, a, a( m+1 ), a( n+m+1 ),
502 $ af, af( m+1 ), af( n+m+1 ),
503 $ af( n+2*m+1 ), iwork, b, lda, x, lda,
504 $ rwork, rwork( nrhs+1 ), work,
505 $ iwork( n+1 ), info )
506*
507* Check error code from SGTRFS.
508*
509 IF( info.NE.0 )
510 $ CALL alaerh( path, 'SGTRFS', info, 0, trans, n, n,
511 $ -1, -1, nrhs, imat, nfail, nerrs,
512 $ nout )
513*
514 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
515 $ result( 4 ) )
516 CALL sgtt05( trans, n, nrhs, a, a( m+1 ), a( n+m+1 ),
517 $ b, lda, x, lda, xact, lda, rwork,
518 $ rwork( nrhs+1 ), result( 5 ) )
519*
520* Print information about the tests that did not pass
521* the threshold.
522*
523 DO 70 k = 2, 6
524 IF( result( k ).GE.thresh ) THEN
525 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
526 $ CALL alahd( nout, path )
527 WRITE( nout, fmt = 9998 )trans, n, nrhs, imat,
528 $ k, result( k )
529 nfail = nfail + 1
530 END IF
531 70 CONTINUE
532 nrun = nrun + 5
533 80 CONTINUE
534 90 CONTINUE
535*
536 100 CONTINUE
537 110 CONTINUE
538*
539* Print a summary of the results.
540*
541 CALL alasum( path, nout, nfail, nrun, nerrs )
542*
543 9999 FORMAT( 12x, 'N =', i5, ',', 10x, ' type ', i2, ', test(', i2,
544 $ ') = ', g12.5 )
545 9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
546 $ i2, ', test(', i2, ') = ', g12.5 )
547 9997 FORMAT( ' NORM =''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
548 $ ', test(', i2, ') = ', g12.5 )
549 RETURN
550*
551* End of SCHKGT
552*
subroutine alasum(type, nout, nfail, nrun, nerrs)
ALASUM
Definition alasum.f:73
subroutine alaerh(path, subnam, info, infoe, opts, m, n, kl, ku, n5, imat, nfail, nerrs, nout)
ALAERH
Definition alaerh.f:147
subroutine alahd(iounit, path)
ALAHD
Definition alahd.f:107
real function sasum(n, sx, incx)
SASUM
Definition sasum.f:72
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgtcon(norm, n, dl, d, du, du2, ipiv, anorm, rcond, work, iwork, info)
SGTCON
Definition sgtcon.f:146
subroutine sgtrfs(trans, n, nrhs, dl, d, du, dlf, df, duf, du2, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
SGTRFS
Definition sgtrfs.f:209
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
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:106
subroutine slarnv(idist, iseed, n, x)
SLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition slarnv.f:97
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine serrge(path, nunit)
SERRGE
Definition serrge.f:55
subroutine sget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
SGET04
Definition sget04.f:102
real function sget06(rcond, rcondc)
SGET06
Definition sget06.f:55
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
Here is the call graph for this function:
Here is the caller graph for this function: