LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvhe_rook.f
Go to the documentation of this file.
1*> \brief \b ZDRVHE_ROOK
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 ZDRVHE_ROOK( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
12* NMAX, A, AFAC, AINV, B, X, XACT, WORK, RWORK,
13* IWORK, NOUT )
14*
15* .. Scalar Arguments ..
16* LOGICAL TSTERR
17* INTEGER NMAX, NN, NOUT, NRHS
18* DOUBLE PRECISION THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NVAL( * )
23* DOUBLE PRECISION RWORK( * )
24* COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
25* $ WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> ZDRVHE_ROOK tests the driver routines ZHESV_ROOK.
35*> \endverbatim
36*
37* Arguments:
38* ==========
39*
40*> \param[in] DOTYPE
41*> \verbatim
42*> DOTYPE is LOGICAL array, dimension (NTYPES)
43*> The matrix types to be used for testing. Matrices of type j
44*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
45*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
46*> \endverbatim
47*>
48*> \param[in] NN
49*> \verbatim
50*> NN is INTEGER
51*> The number of values of N contained in the vector NVAL.
52*> \endverbatim
53*>
54*> \param[in] NVAL
55*> \verbatim
56*> NVAL is INTEGER array, dimension (NN)
57*> The values of the matrix dimension N.
58*> \endverbatim
59*>
60*> \param[in] NRHS
61*> \verbatim
62*> NRHS is INTEGER
63*> The number of right hand side vectors to be generated for
64*> each linear system.
65*> \endverbatim
66*>
67*> \param[in] THRESH
68*> \verbatim
69*> THRESH is DOUBLE PRECISION
70*> The threshold value for the test ratios. A result is
71*> included in the output file if RESULT >= THRESH. To have
72*> every test ratio printed, use THRESH = 0.
73*> \endverbatim
74*>
75*> \param[in] TSTERR
76*> \verbatim
77*> TSTERR is LOGICAL
78*> Flag that indicates whether error exits are to be tested.
79*> \endverbatim
80*>
81*> \param[in] NMAX
82*> \verbatim
83*> NMAX is INTEGER
84*> The maximum value permitted for N, used in dimensioning the
85*> work arrays.
86*> \endverbatim
87*>
88*> \param[out] A
89*> \verbatim
90*> A is COMPLEX*16 array, dimension (NMAX*NMAX)
91*> \endverbatim
92*>
93*> \param[out] AFAC
94*> \verbatim
95*> AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
96*> \endverbatim
97*>
98*> \param[out] AINV
99*> \verbatim
100*> AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
101*> \endverbatim
102*>
103*> \param[out] B
104*> \verbatim
105*> B is COMPLEX*16 array, dimension (NMAX*NRHS)
106*> \endverbatim
107*>
108*> \param[out] X
109*> \verbatim
110*> X is COMPLEX*16 array, dimension (NMAX*NRHS)
111*> \endverbatim
112*>
113*> \param[out] XACT
114*> \verbatim
115*> XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] WORK
119*> \verbatim
120*> WORK is COMPLEX*16 array, dimension (NMAX*max(2,NRHS))
121*> \endverbatim
122*>
123*> \param[out] RWORK
124*> \verbatim
125*> RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
126*> \endverbatim
127*>
128*> \param[out] IWORK
129*> \verbatim
130*> IWORK is INTEGER array, dimension (NMAX)
131*> \endverbatim
132*>
133*> \param[in] NOUT
134*> \verbatim
135*> NOUT is INTEGER
136*> The unit number for output.
137*> \endverbatim
138*
139* Authors:
140* ========
141*
142*> \author Univ. of Tennessee
143*> \author Univ. of California Berkeley
144*> \author Univ. of Colorado Denver
145*> \author NAG Ltd.
146*
147*> \ingroup complex16_lin
148*
149* =====================================================================
150 SUBROUTINE zdrvhe_rook( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
151 $ NMAX, A, AFAC, AINV, B, X, XACT, WORK,
152 $ RWORK, IWORK, NOUT )
153*
154* -- LAPACK test routine --
155* -- LAPACK is a software package provided by Univ. of Tennessee, --
156* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
157*
158* .. Scalar Arguments ..
159 LOGICAL TSTERR
160 INTEGER NMAX, NN, NOUT, NRHS
161 DOUBLE PRECISION THRESH
162* ..
163* .. Array Arguments ..
164 LOGICAL DOTYPE( * )
165 INTEGER IWORK( * ), NVAL( * )
166 DOUBLE PRECISION RWORK( * )
167 COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
168 $ work( * ), x( * ), xact( * )
169* ..
170*
171* =====================================================================
172*
173* .. Parameters ..
174 DOUBLE PRECISION ONE, ZERO
175 PARAMETER ( ONE = 1.0d+0, zero = 0.0d+0 )
176 INTEGER NTYPES, NTESTS
177 parameter( ntypes = 10, ntests = 3 )
178 INTEGER NFACT
179 parameter( nfact = 2 )
180* ..
181* .. Local Scalars ..
182 LOGICAL ZEROT
183 CHARACTER DIST, FACT, TYPE, UPLO, XTYPE
184 CHARACTER*3 MATPATH, PATH
185 INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
186 $ izero, j, k, kl, ku, lda, lwork, mode, n,
187 $ nb, nbmin, nerrs, nfail, nimat, nrun, nt
188 DOUBLE PRECISION AINVNM, ANORM, CNDNUM, RCONDC
189* ..
190* .. Local Arrays ..
191 CHARACTER FACTS( NFACT ), UPLOS( 2 )
192 INTEGER ISEED( 4 ), ISEEDY( 4 )
193 DOUBLE PRECISION RESULT( NTESTS )
194
195* ..
196* .. External Functions ..
197 DOUBLE PRECISION ZLANHE
198 EXTERNAL ZLANHE
199* ..
200* .. External Subroutines ..
201 EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx,
205* ..
206* .. Scalars in Common ..
207 LOGICAL LERR, OK
208 CHARACTER*32 SRNAMT
209 INTEGER INFOT, NUNIT
210* ..
211* .. Common blocks ..
212 COMMON / infoc / infot, nunit, ok, lerr
213 COMMON / srnamc / srnamt
214* ..
215* .. Intrinsic Functions ..
216 INTRINSIC max, min
217* ..
218* .. Data statements ..
219 DATA iseedy / 1988, 1989, 1990, 1991 /
220 DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
221* ..
222* .. Executable Statements ..
223*
224* Initialize constants and the random number seed.
225*
226* Test path
227*
228 path( 1: 1 ) = 'Zomplex precision'
229 path( 2: 3 ) = 'HR'
230*
231* Path to generate matrices
232*
233 matpath( 1: 1 ) = 'Zomplex precision'
234 matpath( 2: 3 ) = 'HE'
235*
236 nrun = 0
237 nfail = 0
238 nerrs = 0
239 DO 10 i = 1, 4
240 iseed( i ) = iseedy( i )
241 10 CONTINUE
242 lwork = max( 2*nmax, nmax*nrhs )
243*
244* Test the error exits
245*
246 IF( tsterr )
247 $ CALL zerrvx( path, nout )
248 infot = 0
249*
250* Set the block size and minimum block size for which the block
251* routine should be used, which will be later returned by ILAENV.
252*
253 nb = 1
254 nbmin = 2
255 CALL xlaenv( 1, nb )
256 CALL xlaenv( 2, nbmin )
257*
258* Do for each value of N in NVAL
259*
260 DO 180 in = 1, nn
261 n = nval( in )
262 lda = max( n, 1 )
263 xtype = 'N'
264 nimat = ntypes
265 IF( n.LE.0 )
266 $ nimat = 1
267*
268 DO 170 imat = 1, nimat
269*
270* Do the tests only if DOTYPE( IMAT ) is true.
271*
272 IF( .NOT.dotype( imat ) )
273 $ GO TO 170
274*
275* Skip types 3, 4, 5, or 6 if the matrix size is too small.
276*
277 zerot = imat.GE.3 .AND. imat.LE.6
278 IF( zerot .AND. n.LT.imat-2 )
279 $ GO TO 170
280*
281* Do first for UPLO = 'U', then for UPLO = 'L'
282*
283 DO 160 iuplo = 1, 2
284 uplo = uplos( iuplo )
285*
286* Begin generate the test matrix A.
287*
288* Set up parameters with ZLATB4 for the matrix generator
289* based on the type of matrix to be generated.
290*
291 CALL zlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
292 $ mode, cndnum, dist )
293*
294* Generate a matrix with ZLATMS.
295*
296 srnamt = 'ZLATMS'
297 CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
298 $ cndnum, anorm, kl, ku, uplo, a, lda,
299 $ work, info )
300*
301* Check error code from ZLATMS and handle error.
302*
303 IF( info.NE.0 ) THEN
304 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
305 $ -1, -1, -1, imat, nfail, nerrs, nout )
306 GO TO 160
307 END IF
308*
309* For types 3-6, zero one or more rows and columns of
310* the matrix to test that INFO is returned correctly.
311*
312 IF( zerot ) THEN
313 IF( imat.EQ.3 ) THEN
314 izero = 1
315 ELSE IF( imat.EQ.4 ) THEN
316 izero = n
317 ELSE
318 izero = n / 2 + 1
319 END IF
320*
321 IF( imat.LT.6 ) THEN
322*
323* Set row and column IZERO to zero.
324*
325 IF( iuplo.EQ.1 ) THEN
326 ioff = ( izero-1 )*lda
327 DO 20 i = 1, izero - 1
328 a( ioff+i ) = zero
329 20 CONTINUE
330 ioff = ioff + izero
331 DO 30 i = izero, n
332 a( ioff ) = zero
333 ioff = ioff + lda
334 30 CONTINUE
335 ELSE
336 ioff = izero
337 DO 40 i = 1, izero - 1
338 a( ioff ) = zero
339 ioff = ioff + lda
340 40 CONTINUE
341 ioff = ioff - izero
342 DO 50 i = izero, n
343 a( ioff+i ) = zero
344 50 CONTINUE
345 END IF
346 ELSE
347 IF( iuplo.EQ.1 ) THEN
348*
349* Set the first IZERO rows and columns to zero.
350*
351 ioff = 0
352 DO 70 j = 1, n
353 i2 = min( j, izero )
354 DO 60 i = 1, i2
355 a( ioff+i ) = zero
356 60 CONTINUE
357 ioff = ioff + lda
358 70 CONTINUE
359 ELSE
360*
361* Set the first IZERO rows and columns to zero.
362*
363 ioff = 0
364 DO 90 j = 1, n
365 i1 = max( j, izero )
366 DO 80 i = i1, n
367 a( ioff+i ) = zero
368 80 CONTINUE
369 ioff = ioff + lda
370 90 CONTINUE
371 END IF
372 END IF
373 ELSE
374 izero = 0
375 END IF
376*
377* End generate the test matrix A.
378*
379*
380 DO 150 ifact = 1, nfact
381*
382* Do first for FACT = 'F', then for other values.
383*
384 fact = facts( ifact )
385*
386* Compute the condition number for comparison with
387* the value returned by ZHESVX_ROOK.
388*
389 IF( zerot ) THEN
390 IF( ifact.EQ.1 )
391 $ GO TO 150
392 rcondc = zero
393*
394 ELSE IF( ifact.EQ.1 ) THEN
395*
396* Compute the 1-norm of A.
397*
398 anorm = zlanhe( '1', uplo, n, a, lda, rwork )
399*
400* Factor the matrix A.
401*
402
403 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
404 CALL zhetrf_rook( uplo, n, afac, lda, iwork, work,
405 $ lwork, info )
406*
407* Compute inv(A) and take its norm.
408*
409 CALL zlacpy( uplo, n, n, afac, lda, ainv, lda )
410 lwork = (n+nb+1)*(nb+3)
411 CALL zhetri_rook( uplo, n, ainv, lda, iwork,
412 $ work, info )
413 ainvnm = zlanhe( '1', uplo, n, ainv, lda, rwork )
414*
415* Compute the 1-norm condition number of A.
416*
417 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
418 rcondc = one
419 ELSE
420 rcondc = ( one / anorm ) / ainvnm
421 END IF
422 END IF
423*
424* Form an exact solution and set the right hand side.
425*
426 srnamt = 'ZLARHS'
427 CALL zlarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
428 $ nrhs, a, lda, xact, lda, b, lda, iseed,
429 $ info )
430 xtype = 'C'
431*
432* --- Test ZHESV_ROOK ---
433*
434 IF( ifact.EQ.2 ) THEN
435 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
436 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
437*
438* Factor the matrix and solve the system using
439* ZHESV_ROOK.
440*
441 srnamt = 'ZHESV_ROOK'
442 CALL zhesv_rook( uplo, n, nrhs, afac, lda, iwork,
443 $ x, lda, work, lwork, info )
444*
445* Adjust the expected value of INFO to account for
446* pivoting.
447*
448 k = izero
449 IF( k.GT.0 ) THEN
450 100 CONTINUE
451 IF( iwork( k ).LT.0 ) THEN
452 IF( iwork( k ).NE.-k ) THEN
453 k = -iwork( k )
454 GO TO 100
455 END IF
456 ELSE IF( iwork( k ).NE.k ) THEN
457 k = iwork( k )
458 GO TO 100
459 END IF
460 END IF
461*
462* Check error code from ZHESV_ROOK and handle error.
463*
464 IF( info.NE.k ) THEN
465 CALL alaerh( path, 'ZHESV_ROOK', info, k, uplo,
466 $ n, n, -1, -1, nrhs, imat, nfail,
467 $ nerrs, nout )
468 GO TO 120
469 ELSE IF( info.NE.0 ) THEN
470 GO TO 120
471 END IF
472*
473*+ TEST 1 Reconstruct matrix from factors and compute
474* residual.
475*
476 CALL zhet01_rook( uplo, n, a, lda, afac, lda,
477 $ iwork, ainv, lda, rwork,
478 $ result( 1 ) )
479*
480*+ TEST 2 Compute residual of the computed solution.
481*
482 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
483 CALL zpot02( uplo, n, nrhs, a, lda, x, lda, work,
484 $ lda, rwork, result( 2 ) )
485*
486*+ TEST 3
487* Check solution from generated exact solution.
488*
489 CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
490 $ result( 3 ) )
491 nt = 3
492*
493* Print information about the tests that did not pass
494* the threshold.
495*
496 DO 110 k = 1, nt
497 IF( result( k ).GE.thresh ) THEN
498 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
499 $ CALL aladhd( nout, path )
500 WRITE( nout, fmt = 9999 )'ZHESV_ROOK', uplo,
501 $ n, imat, k, result( k )
502 nfail = nfail + 1
503 END IF
504 110 CONTINUE
505 nrun = nrun + nt
506 120 CONTINUE
507 END IF
508*
509 150 CONTINUE
510*
511 160 CONTINUE
512 170 CONTINUE
513 180 CONTINUE
514*
515* Print a summary of the results.
516*
517 CALL alasvm( path, nout, nfail, nrun, nerrs )
518*
519 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
520 $ ', test ', i2, ', ratio =', g12.5 )
521 RETURN
522*
523* End of ZDRVHE_ROOK
524*
525 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.f:208
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 zhesv_rook(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
ZHESV_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using the ...
Definition zhesv_rook.f:205
subroutine zhetrf_rook(uplo, n, a, lda, ipiv, work, lwork, info)
ZHETRF_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bun...
subroutine zhetri_rook(uplo, n, a, lda, ipiv, work, info)
ZHETRI_ROOK computes the inverse of HE matrix using the factorization obtained with the bounded Bunch...
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 zdrvhe_rook(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
ZDRVHE_ROOK
subroutine zerrvx(path, nunit)
ZERRVX
Definition zerrvx.f:55
subroutine zget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
ZGET04
Definition zget04.f:102
subroutine zhet01_rook(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
ZHET01_ROOK
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
subroutine zpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZPOT02
Definition zpot02.f:127