LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sdrvsy_rook.f
Go to the documentation of this file.
1*> \brief \b SDRVSY_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 SDRVSY_ROOK( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
12* $ NMAX, A, AFAC, AINV, B, X, XACT, WORK,
13* $ RWORK, 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 A( * ), AFAC( * ), AINV( * ), B( * ),
24* $ RWORK( * ), WORK( * ), X( * ), XACT( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*>
33*> SDRVSY_ROOK tests the driver routines SSYSV_ROOK.
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 side vectors to be generated for
63*> each linear system.
64*> \endverbatim
65*>
66*> \param[in] THRESH
67*> \verbatim
68*> THRESH is DOUBLE PRECISION
69*> The threshold value for the test ratios. A result is
70*> included in the output file if RESULT >= THRESH. To have
71*> every test ratio printed, use THRESH = 0.
72*> \endverbatim
73*>
74*> \param[in] TSTERR
75*> \verbatim
76*> TSTERR is LOGICAL
77*> Flag that indicates whether error exits are to be tested.
78*> \endverbatim
79*>
80*> \param[in] NMAX
81*> \verbatim
82*> NMAX is INTEGER
83*> The maximum value permitted for N, used in dimensioning the
84*> work arrays.
85*> \endverbatim
86*>
87*> \param[out] A
88*> \verbatim
89*> A is REAL array, dimension (NMAX*NMAX)
90*> \endverbatim
91*>
92*> \param[out] AFAC
93*> \verbatim
94*> AFAC is REAL array, dimension (NMAX*NMAX)
95*> \endverbatim
96*>
97*> \param[out] AINV
98*> \verbatim
99*> AINV is REAL array, dimension (NMAX*NMAX)
100*> \endverbatim
101*>
102*> \param[out] B
103*> \verbatim
104*> B is REAL array, dimension (NMAX*NRHS)
105*> \endverbatim
106*>
107*> \param[out] X
108*> \verbatim
109*> X is REAL array, dimension (NMAX*NRHS)
110*> \endverbatim
111*>
112*> \param[out] XACT
113*> \verbatim
114*> XACT is REAL array, dimension (NMAX*NRHS)
115*> \endverbatim
116*>
117*> \param[out] WORK
118*> \verbatim
119*> WORK is REAL array, dimension
120*> (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 (2*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 double_lin
148*
149* =====================================================================
150 SUBROUTINE sdrvsy_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 A( * ), AFAC( * ), AINV( * ), B( * ),
167 $ rwork( * ), work( * ), x( * ), xact( * )
168* ..
169*
170* =====================================================================
171*
172* .. Parameters ..
173 REAL ONE, ZERO
174 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
175 INTEGER NTYPES, NTESTS
176 parameter( ntypes = 10, ntests = 3 )
177 INTEGER NFACT
178 parameter( nfact = 2 )
179* ..
180* .. Local Scalars ..
181 LOGICAL ZEROT
182 CHARACTER DIST, FACT, TYPE, UPLO, XTYPE
183 CHARACTER*3 PATH, MATPATH
184 INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
185 $ izero, j, k, kl, ku, lda, lwork, mode, n,
186 $ nb, nbmin, nerrs, nfail, nimat, nrun, nt
187 REAL AINVNM, ANORM, CNDNUM, RCONDC
188* ..
189* .. Local Arrays ..
190 CHARACTER FACTS( NFACT ), UPLOS( 2 )
191 INTEGER ISEED( 4 ), ISEEDY( 4 )
192 REAL RESULT( NTESTS )
193* ..
194* .. External Functions ..
195 REAL SLANSY
196 EXTERNAL SLANSY
197* ..
198* .. External Subroutines ..
199 EXTERNAL aladhd, alaerh, alasvm, serrvx, sget04, slacpy,
202 $ ssytri_rook,
203 $ xlaenv
204* ..
205* .. Scalars in Common ..
206 LOGICAL LERR, OK
207 CHARACTER*32 SRNAMT
208 INTEGER INFOT, NUNIT
209* ..
210* .. Common blocks ..
211 COMMON / infoc / infot, nunit, ok, lerr
212 COMMON / srnamc / srnamt
213* ..
214* .. Intrinsic Functions ..
215 INTRINSIC max, min
216* ..
217* .. Data statements ..
218 DATA iseedy / 1988, 1989, 1990, 1991 /
219 DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
220* ..
221* .. Executable Statements ..
222*
223* Initialize constants and the random number seed.
224*
225* Test path
226*
227 path( 1: 1 ) = 'Single precision'
228 path( 2: 3 ) = 'SR'
229*
230* Path to generate matrices
231*
232 matpath( 1: 1 ) = 'Single precision'
233 matpath( 2: 3 ) = 'SY'
234*
235 nrun = 0
236 nfail = 0
237 nerrs = 0
238 DO 10 i = 1, 4
239 iseed( i ) = iseedy( i )
240 10 CONTINUE
241 lwork = max( 2*nmax, nmax*nrhs )
242*
243* Test the error exits
244*
245 IF( tsterr )
246 $ CALL serrvx( path, nout )
247 infot = 0
248*
249* Set the block size and minimum block size for which the block
250* routine should be used, which will be later returned by ILAENV.
251*
252 nb = 1
253 nbmin = 2
254 CALL xlaenv( 1, nb )
255 CALL xlaenv( 2, nbmin )
256*
257* Do for each value of N in NVAL
258*
259 DO 180 in = 1, nn
260 n = nval( in )
261 lda = max( n, 1 )
262 xtype = 'N'
263 nimat = ntypes
264 IF( n.LE.0 )
265 $ nimat = 1
266*
267 DO 170 imat = 1, nimat
268*
269* Do the tests only if DOTYPE( IMAT ) is true.
270*
271 IF( .NOT.dotype( imat ) )
272 $ GO TO 170
273*
274* Skip types 3, 4, 5, or 6 if the matrix size is too small.
275*
276 zerot = imat.GE.3 .AND. imat.LE.6
277 IF( zerot .AND. n.LT.imat-2 )
278 $ GO TO 170
279*
280* Do first for UPLO = 'U', then for UPLO = 'L'
281*
282 DO 160 iuplo = 1, 2
283 uplo = uplos( iuplo )
284*
285* Begin generate the test matrix A.
286*
287* Set up parameters with SLATB4 for the matrix generator
288* based on the type of matrix to be generated.
289*
290 CALL slatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
291 $ mode, cndnum, dist )
292*
293* Generate a matrix with SLATMS.
294*
295 srnamt = 'SLATMS'
296 CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
297 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
298 $ info )
299*
300* Check error code from SLATMS and handle error.
301*
302 IF( info.NE.0 ) THEN
303 CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
304 $ -1, -1, imat, nfail, nerrs, nout )
305*
306* Skip all tests for this generated matrix
307*
308 GO TO 160
309 END IF
310*
311* For types 3-6, zero one or more rows and columns of the
312* matrix to test that INFO is returned correctly.
313*
314 IF( zerot ) THEN
315 IF( imat.EQ.3 ) THEN
316 izero = 1
317 ELSE IF( imat.EQ.4 ) THEN
318 izero = n
319 ELSE
320 izero = n / 2 + 1
321 END IF
322*
323 IF( imat.LT.6 ) THEN
324*
325* Set row and column IZERO to zero.
326*
327 IF( iuplo.EQ.1 ) THEN
328 ioff = ( izero-1 )*lda
329 DO 20 i = 1, izero - 1
330 a( ioff+i ) = zero
331 20 CONTINUE
332 ioff = ioff + izero
333 DO 30 i = izero, n
334 a( ioff ) = zero
335 ioff = ioff + lda
336 30 CONTINUE
337 ELSE
338 ioff = izero
339 DO 40 i = 1, izero - 1
340 a( ioff ) = zero
341 ioff = ioff + lda
342 40 CONTINUE
343 ioff = ioff - izero
344 DO 50 i = izero, n
345 a( ioff+i ) = zero
346 50 CONTINUE
347 END IF
348 ELSE
349 ioff = 0
350 IF( iuplo.EQ.1 ) THEN
351*
352* Set the first IZERO rows and columns to zero.
353*
354 DO 70 j = 1, n
355 i2 = min( j, izero )
356 DO 60 i = 1, i2
357 a( ioff+i ) = zero
358 60 CONTINUE
359 ioff = ioff + lda
360 70 CONTINUE
361 ELSE
362*
363* Set the last IZERO rows and columns to zero.
364*
365 DO 90 j = 1, n
366 i1 = max( j, izero )
367 DO 80 i = i1, n
368 a( ioff+i ) = zero
369 80 CONTINUE
370 ioff = ioff + lda
371 90 CONTINUE
372 END IF
373 END IF
374 ELSE
375 izero = 0
376 END IF
377*
378* End generate the test matrix A.
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 DSYSVX_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 = slansy( '1', uplo, n, a, lda, rwork )
399*
400* Factor the matrix A.
401*
402 CALL slacpy( uplo, n, n, a, lda, afac, lda )
403 CALL ssytrf_rook( uplo, n, afac, lda, iwork, work,
404 $ lwork, info )
405*
406* Compute inv(A) and take its norm.
407*
408 CALL slacpy( uplo, n, n, afac, lda, ainv, lda )
409 lwork = (n+nb+1)*(nb+3)
410 CALL ssytri_rook( uplo, n, ainv, lda, iwork,
411 $ work, info )
412 ainvnm = slansy( '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 = 'SLARHS'
426 CALL slarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
427 $ nrhs, a, lda, xact, lda, b, lda, iseed,
428 $ info )
429 xtype = 'C'
430*
431* --- Test SSYSV_ROOK ---
432*
433 IF( ifact.EQ.2 ) THEN
434 CALL slacpy( uplo, n, n, a, lda, afac, lda )
435 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
436*
437* Factor the matrix and solve the system using
438* SSYSV_ROOK.
439*
440 srnamt = 'SSYSV_ROOK'
441 CALL ssysv_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 SSYSV_ROOK and handle error.
462*
463 IF( info.NE.k ) THEN
464 CALL alaerh( path, 'SSYSV_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 ssyt01_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 slacpy( 'Full', n, nrhs, b, lda, work, lda )
482 CALL spot02( 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 sget04( 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 )'SSYSV_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 SDRVSY_ROOK
523*
524 END
subroutine alasvm(type, nout, nfail, nrun, nerrs)
ALASVM
Definition alasvm.f:73
subroutine slarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
SLARHS
Definition slarhs.f:205
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 ssysv_rook(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
SSYSV_ROOK computes the solution to system of linear equations A * X = B for SY matrices
Definition ssysv_rook.f:204
subroutine ssytrf_rook(uplo, n, a, lda, ipiv, work, lwork, info)
SSYTRF_ROOK
subroutine ssytri_rook(uplo, n, a, lda, ipiv, work, info)
SSYTRI_ROOK
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
subroutine sdrvsy_rook(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
SDRVSY_ROOK
subroutine serrvx(path, nunit)
SERRVX
Definition serrvx.f:55
subroutine sget04(n, nrhs, x, ldx, xact, ldxact, rcond, resid)
SGET04
Definition sget04.f:102
subroutine slatb4(path, imat, m, n, type, kl, ku, anorm, mode, cndnum, dist)
SLATB4
Definition slatb4.f:120
subroutine slatms(m, n, dist, iseed, sym, d, mode, cond, dmax, kl, ku, pack, a, lda, work, info)
SLATMS
Definition slatms.f:321
subroutine spot02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
SPOT02
Definition spot02.f:127
subroutine spot05(uplo, n, nrhs, a, lda, b, ldb, x, ldx, xact, ldxact, ferr, berr, reslts)
SPOT05
Definition spot05.f:164
subroutine ssyt01_rook(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
SSYT01_ROOK