LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cdrvhe_rook.f
Go to the documentation of this file.
1*> \brief \b CDRVHE_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 CDRVHE_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* REAL THRESH
19* ..
20* .. Array Arguments ..
21* LOGICAL DOTYPE( * )
22* INTEGER IWORK( * ), NVAL( * )
23* REAL RWORK( * )
24* COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ),
25* $ WORK( * ), X( * ), XACT( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> CDRVHE_ROOK tests the driver routines CHESV_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 REAL
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 array, dimension (NMAX*NMAX)
91*> \endverbatim
92*>
93*> \param[out] AFAC
94*> \verbatim
95*> AFAC is COMPLEX array, dimension (NMAX*NMAX)
96*> \endverbatim
97*>
98*> \param[out] AINV
99*> \verbatim
100*> AINV is COMPLEX array, dimension (NMAX*NMAX)
101*> \endverbatim
102*>
103*> \param[out] B
104*> \verbatim
105*> B is COMPLEX array, dimension (NMAX*NRHS)
106*> \endverbatim
107*>
108*> \param[out] X
109*> \verbatim
110*> X is COMPLEX array, dimension (NMAX*NRHS)
111*> \endverbatim
112*>
113*> \param[out] XACT
114*> \verbatim
115*> XACT is COMPLEX array, dimension (NMAX*NRHS)
116*> \endverbatim
117*>
118*> \param[out] WORK
119*> \verbatim
120*> WORK is COMPLEX array, dimension (NMAX*max(2,NRHS))
121*> \endverbatim
122*>
123*> \param[out] RWORK
124*> \verbatim
125*> RWORK is REAL 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 complex_lin
148*
149* =====================================================================
150 SUBROUTINE cdrvhe_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 REAL THRESH
162* ..
163* .. Array Arguments ..
164 LOGICAL DOTYPE( * )
165 INTEGER IWORK( * ), NVAL( * )
166 REAL RWORK( * )
167 COMPLEX A( * ), AFAC( * ), AINV( * ), B( * ),
168 $ work( * ), x( * ), xact( * )
169* ..
170*
171* =====================================================================
172*
173* .. Parameters ..
174 REAL ONE, ZERO
175 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+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 REAL AINVNM, ANORM, CNDNUM, RCONDC
189* ..
190* .. Local Arrays ..
191 CHARACTER FACTS( NFACT ), UPLOS( 2 )
192 INTEGER ISEED( 4 ), ISEEDY( 4 )
193 REAL RESULT( NTESTS )
194
195* ..
196* .. External Functions ..
197 REAL CLANHE
198 EXTERNAL CLANHE
199* ..
200* .. External Subroutines ..
201 EXTERNAL aladhd, alaerh, alasvm, xlaenv, cerrvx,
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 ) = 'Complex precision'
229 path( 2: 3 ) = 'HR'
230*
231* Path to generate matrices
232*
233 matpath( 1: 1 ) = 'Complex 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 cerrvx( 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 CLATB4 for the matrix generator
289* based on the type of matrix to be generated.
290*
291 CALL clatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
292 $ mode, cndnum, dist )
293*
294* Generate a matrix with CLATMS.
295*
296 srnamt = 'CLATMS'
297 CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
298 $ cndnum, anorm, kl, ku, uplo, a, lda,
299 $ work, info )
300*
301* Check error code from CLATMS and handle error.
302*
303 IF( info.NE.0 ) THEN
304 CALL alaerh( path, 'CLATMS', 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 CHESVX_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 = clanhe( '1', uplo, n, a, lda, rwork )
399*
400* Factor the matrix A.
401*
402 CALL clacpy( uplo, n, n, a, lda, afac, lda )
403 CALL chetrf_rook( uplo, n, afac, lda, iwork, work,
404 $ lwork, info )
405*
406* Compute inv(A) and take its norm.
407*
408 CALL clacpy( uplo, n, n, afac, lda, ainv, lda )
409 lwork = (n+nb+1)*(nb+3)
410 CALL chetri_rook( uplo, n, ainv, lda, iwork,
411 $ work, info )
412 ainvnm = clanhe( '1', uplo, n, ainv, lda, rwork )
413*
414* Compute the 1-norm condition number of A.
415*
416 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
417 rcondc = one
418 ELSE
419 rcondc = ( one / anorm ) / ainvnm
420 END IF
421 END IF
422*
423* Form an exact solution and set the right hand side.
424*
425 srnamt = 'CLARHS'
426 CALL clarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
427 $ nrhs, a, lda, xact, lda, b, lda, iseed,
428 $ info )
429 xtype = 'C'
430*
431* --- Test CHESV_ROOK ---
432*
433 IF( ifact.EQ.2 ) THEN
434 CALL clacpy( uplo, n, n, a, lda, afac, lda )
435 CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
436*
437* Factor the matrix and solve the system using
438* CHESV_ROOK.
439*
440 srnamt = 'CHESV_ROOK'
441 CALL chesv_rook( uplo, n, nrhs, afac, lda, iwork,
442 $ x, lda, work, lwork, info )
443*
444* Adjust the expected value of INFO to account for
445* pivoting.
446*
447 k = izero
448 IF( k.GT.0 ) THEN
449 100 CONTINUE
450 IF( iwork( k ).LT.0 ) THEN
451 IF( iwork( k ).NE.-k ) THEN
452 k = -iwork( k )
453 GO TO 100
454 END IF
455 ELSE IF( iwork( k ).NE.k ) THEN
456 k = iwork( k )
457 GO TO 100
458 END IF
459 END IF
460*
461* Check error code from CHESV_ROOK and handle error.
462*
463 IF( info.NE.k ) THEN
464 CALL alaerh( path, 'CHESV_ROOK', info, k, uplo,
465 $ n, n, -1, -1, nrhs, imat, nfail,
466 $ nerrs, nout )
467 GO TO 120
468 ELSE IF( info.NE.0 ) THEN
469 GO TO 120
470 END IF
471*
472*+ TEST 1 Reconstruct matrix from factors and compute
473* residual.
474*
475 CALL chet01_rook( uplo, n, a, lda, afac, lda,
476 $ iwork, ainv, lda, rwork,
477 $ result( 1 ) )
478*
479*+ TEST 2 Compute residual of the computed solution.
480*
481 CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
482 CALL cpot02( uplo, n, nrhs, a, lda, x, lda, work,
483 $ lda, rwork, result( 2 ) )
484*
485*+ TEST 3
486* Check solution from generated exact solution.
487*
488 CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
489 $ result( 3 ) )
490 nt = 3
491*
492* Print information about the tests that did not pass
493* the threshold.
494*
495 DO 110 k = 1, nt
496 IF( result( k ).GE.thresh ) THEN
497 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
498 $ CALL aladhd( nout, path )
499 WRITE( nout, fmt = 9999 )'CHESV_ROOK', uplo,
500 $ n, imat, k, result( k )
501 nfail = nfail + 1
502 END IF
503 110 CONTINUE
504 nrun = nrun + nt
505 120 CONTINUE
506 END IF
507*
508 150 CONTINUE
509*
510 160 CONTINUE
511 170 CONTINUE
512 180 CONTINUE
513*
514* Print a summary of the results.
515*
516 CALL alasvm( path, nout, nfail, nrun, nerrs )
517*
518 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
519 $ ', test ', i2, ', ratio =', g12.5 )
520 RETURN
521*
522* End of CDRVHE_ROOK
523*
524 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine clarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
CLARHS
Definition clarhs.f:208
subroutine xlaenv(ispec, nvalue)
XLAENV
Definition xlaenv.f:81
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 cdrvhe_rook(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
CDRVHE_ROOK
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 chet01_rook(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
CHET01_ROOK
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 cpot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
CPOT02
Definition cpot02.f:127
subroutine chesv_rook(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
CHESV_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using the ...
Definition chesv_rook.f:205
subroutine chetrf_rook(uplo, n, a, lda, ipiv, work, lwork, info)
CHETRF_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bun...
subroutine chetri_rook(uplo, n, a, lda, ipiv, work, info)
CHETRI_ROOK computes the inverse of HE matrix using the factorization obtained with the bounded Bunch...
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