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

◆ cchkgt()

subroutine cchkgt ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nns,
integer, dimension( * )  nsval,
real  thresh,
logical  tsterr,
complex, dimension( * )  a,
complex, dimension( * )  af,
complex, dimension( * )  b,
complex, dimension( * )  x,
complex, dimension( * )  xact,
complex, dimension( * )  work,
real, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

CCHKGT

Purpose:
 CCHKGT tests CGTTRF, -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 COMPLEX array, dimension (NMAX*4)
[out]AF
          AF is COMPLEX array, dimension (NMAX*4)
[out]B
          B is COMPLEX array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is REAL array, dimension
                      (max(NMAX)+2*NSMAX)
[out]IWORK
          IWORK is INTEGER array, dimension (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 145 of file cchkgt.f.

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