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

◆ zchkgt()

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

ZCHKGT

Purpose:
 ZCHKGT tests ZGTTRF, -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 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.
[out]A
          A is COMPLEX*16 array, dimension (NMAX*4)
[out]AF
          AF is COMPLEX*16 array, dimension (NMAX*4)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NSMAX)
          where NSMAX is the largest entry in NSVAL.
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NSMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NSMAX))
[out]RWORK
          RWORK is DOUBLE PRECISION 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 zchkgt.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 DOUBLE PRECISION THRESH
156* ..
157* .. Array Arguments ..
158 LOGICAL DOTYPE( * )
159 INTEGER IWORK( * ), NSVAL( * ), NVAL( * )
160 DOUBLE PRECISION RWORK( * )
161 COMPLEX*16 A( * ), AF( * ), B( * ), WORK( * ), X( * ),
162 $ XACT( * )
163* ..
164*
165* =====================================================================
166*
167* .. Parameters ..
168 DOUBLE PRECISION ONE, ZERO
169 parameter( one = 1.0d+0, zero = 0.0d+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 DOUBLE PRECISION AINVNM, ANORM, COND, RCOND, RCONDC, RCONDI,
183 $ RCONDO
184* ..
185* .. Local Arrays ..
186 CHARACTER TRANSS( 3 )
187 INTEGER ISEED( 4 ), ISEEDY( 4 )
188 DOUBLE PRECISION RESULT( NTESTS )
189 COMPLEX*16 Z( 3 )
190* ..
191* .. External Functions ..
192 DOUBLE PRECISION DGET06, DZASUM, ZLANGT
193 EXTERNAL dget06, dzasum, zlangt
194* ..
195* .. External Subroutines ..
196 EXTERNAL alaerh, alahd, alasum, zcopy, zdscal, zerrge,
199 $ zlatms
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 ) = 'Zomplex 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 zerrge( 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 ZLATB4.
253*
254 CALL zlatb4( 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 = 'ZLATMS'
264 CALL zlatms( 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 ZLATMS.
269*
270 IF( info.NE.0 ) THEN
271 CALL alaerh( path, 'ZLATMS', 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 zcopy( n-1, af( 4 ), 3, a, 1 )
279 CALL zcopy( n-1, af( 3 ), 3, a( n+m+1 ), 1 )
280 END IF
281 CALL zcopy( 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 zlarnv( 2, iseed, n+2*m, a )
293 IF( anorm.NE.one )
294 $ CALL zdscal( 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 zcopy( n+2*m, a, 1, af, 1 )
349 srnamt = 'ZGTTRF'
350 CALL zgttrf( n, af, af( m+1 ), af( n+m+1 ), af( n+2*m+1 ),
351 $ iwork, info )
352*
353* Check error code from ZGTTRF.
354*
355 IF( info.NE.izero )
356 $ CALL alaerh( path, 'ZGTTRF', info, izero, ' ', n, n, 1,
357 $ 1, -1, imat, nfail, nerrs, nout )
358 trfcon = info.NE.0
359*
360 CALL zgtt01( 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 = zlangt( norm, n, a, a( m+1 ), a( n+m+1 ) )
382*
383 IF( .NOT.trfcon ) THEN
384*
385* Use ZGTTRS 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 zgttrs( 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, dzasum( 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 = 'ZGTCON'
421 CALL zgtcon( 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 ZGTCON.
426*
427 IF( info.NE.0 )
428 $ CALL alaerh( path, 'ZGTCON', info, 0, norm, n, n, -1,
429 $ -1, -1, imat, nfail, nerrs, nout )
430*
431 result( 7 ) = dget06( 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 zlarnv( 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 zlagtm( trans, n, nrhs, one, a, a( m+1 ),
472 $ a( n+m+1 ), xact, lda, zero, b, lda )
473*
474*+ TEST 2
475* Solve op(A) * X = B and compute the residual.
476*
477 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
478 srnamt = 'ZGTTRS'
479 CALL zgttrs( trans, n, nrhs, af, af( m+1 ),
480 $ af( n+m+1 ), af( n+2*m+1 ), iwork, x,
481 $ lda, info )
482*
483* Check error code from ZGTTRS.
484*
485 IF( info.NE.0 )
486 $ CALL alaerh( path, 'ZGTTRS', info, 0, trans, n, n,
487 $ -1, -1, nrhs, imat, nfail, nerrs,
488 $ nout )
489*
490 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
491 CALL zgtt02( trans, n, nrhs, a, a( m+1 ), a( n+m+1 ),
492 $ x, lda, work, lda, result( 2 ) )
493*
494*+ TEST 3
495* Check solution from generated exact solution.
496*
497 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
498 $ result( 3 ) )
499*
500*+ TESTS 4, 5, and 6
501* Use iterative refinement to improve the solution.
502*
503 srnamt = 'ZGTRFS'
504 CALL zgtrfs( trans, n, nrhs, a, a( m+1 ), a( n+m+1 ),
505 $ af, af( m+1 ), af( n+m+1 ),
506 $ af( n+2*m+1 ), iwork, b, lda, x, lda,
507 $ rwork, rwork( nrhs+1 ), work,
508 $ rwork( 2*nrhs+1 ), info )
509*
510* Check error code from ZGTRFS.
511*
512 IF( info.NE.0 )
513 $ CALL alaerh( path, 'ZGTRFS', info, 0, trans, n, n,
514 $ -1, -1, nrhs, imat, nfail, nerrs,
515 $ nout )
516*
517 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
518 $ result( 4 ) )
519 CALL zgtt05( trans, n, nrhs, a, a( m+1 ), a( n+m+1 ),
520 $ b, lda, x, lda, xact, lda, rwork,
521 $ rwork( nrhs+1 ), result( 5 ) )
522*
523* Print information about the tests that did not pass the
524* threshold.
525*
526 DO 70 k = 2, 6
527 IF( result( k ).GE.thresh ) THEN
528 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
529 $ CALL alahd( nout, path )
530 WRITE( nout, fmt = 9998 )trans, n, nrhs, imat,
531 $ k, result( k )
532 nfail = nfail + 1
533 END IF
534 70 CONTINUE
535 nrun = nrun + 5
536 80 CONTINUE
537 90 CONTINUE
538 100 CONTINUE
539 110 CONTINUE
540*
541* Print a summary of the results.
542*
543 CALL alasum( path, nout, nfail, nrun, nerrs )
544*
545 9999 FORMAT( 12x, 'N =', i5, ',', 10x, ' type ', i2, ', test(', i2,
546 $ ') = ', g12.5 )
547 9998 FORMAT( ' TRANS=''', a1, ''', N =', i5, ', NRHS=', i3, ', type ',
548 $ i2, ', test(', i2, ') = ', g12.5 )
549 9997 FORMAT( ' NORM =''', a1, ''', N =', i5, ',', 10x, ' type ', i2,
550 $ ', test(', i2, ') = ', g12.5 )
551 RETURN
552*
553* End of ZCHKGT
554*
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
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
double precision function dzasum(n, zx, incx)
DZASUM
Definition dzasum.f:72
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgtcon(norm, n, dl, d, du, du2, ipiv, anorm, rcond, work, info)
ZGTCON
Definition zgtcon.f:141
subroutine zgtrfs(trans, n, nrhs, dl, d, du, dlf, df, duf, du2, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZGTRFS
Definition zgtrfs.f:210
subroutine zgttrf(n, dl, d, du, du2, ipiv, info)
ZGTTRF
Definition zgttrf.f:124
subroutine zgttrs(trans, n, nrhs, dl, d, du, du2, ipiv, b, ldb, info)
ZGTTRS
Definition zgttrs.f:138
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlagtm(trans, n, nrhs, alpha, dl, d, du, x, ldx, beta, b, ldb)
ZLAGTM performs a matrix-matrix product of the form C = αAB+βC, where A is a tridiagonal matrix,...
Definition zlagtm.f:145
double precision function zlangt(norm, n, dl, d, du)
ZLANGT returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlangt.f:106
subroutine zlarnv(idist, iseed, n, x)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition zlarnv.f:99
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zerrge(path, nunit)
ZERRGE
Definition zerrge.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zgtt01(n, dl, d, du, dlf, df, duf, du2, ipiv, work, ldwork, rwork, resid)
ZGTT01
Definition zgtt01.f:134
subroutine zgtt02(trans, n, nrhs, dl, d, du, x, ldx, b, ldb, resid)
ZGTT02
Definition zgtt02.f:124
subroutine zgtt05(trans, n, nrhs, dl, d, du, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
ZGTT05
Definition zgtt05.f:165
subroutine zlatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
ZLATB4
Definition zlatb4.f:121
subroutine zlatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
ZLATMS
Definition zlatms.f:332
Here is the call graph for this function:
Here is the caller graph for this function: