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