LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sdrvsy_rk.f
Go to the documentation of this file.
1*> \brief \b SDRVSY_RK
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_RK( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
12* $ NMAX, A, AFAC, E, 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( * ), E( * ), AINV( * ), B( * ),
24* $ RWORK( * ), WORK( * ), X( * ), XACT( * )
25* ..
26*
27*
28*> \par Purpose:
29* =============
30*>
31*> \verbatim
32*> SDRVSY_RK tests the driver routines SSYSV_RK.
33*> \endverbatim
34*
35* Arguments:
36* ==========
37*
38*> \param[in] DOTYPE
39*> \verbatim
40*> DOTYPE is LOGICAL array, dimension (NTYPES)
41*> The matrix types to be used for testing. Matrices of type j
42*> (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
43*> .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
44*> \endverbatim
45*>
46*> \param[in] NN
47*> \verbatim
48*> NN is INTEGER
49*> The number of values of N contained in the vector NVAL.
50*> \endverbatim
51*>
52*> \param[in] NVAL
53*> \verbatim
54*> NVAL is INTEGER array, dimension (NN)
55*> The values of the matrix dimension N.
56*> \endverbatim
57*>
58*> \param[in] NRHS
59*> \verbatim
60*> NRHS is INTEGER
61*> The number of right hand side vectors to be generated for
62*> each linear system.
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[in] NMAX
80*> \verbatim
81*> NMAX is INTEGER
82*> The maximum value permitted for N, used in dimensioning the
83*> work arrays.
84*> \endverbatim
85*>
86*> \param[out] A
87*> \verbatim
88*> A is REAL array, dimension (NMAX*NMAX)
89*> \endverbatim
90*>
91*> \param[out] AFAC
92*> \verbatim
93*> AFAC is REAL array, dimension (NMAX*NMAX)
94*> \endverbatim
95*>
96*> \param[out] E
97*> \verbatim
98*> E is REAL array, dimension (NMAX)
99*> \endverbatim
100*>
101*> \param[out] AINV
102*> \verbatim
103*> AINV is REAL array, dimension (NMAX*NMAX)
104*> \endverbatim
105*>
106*> \param[out] B
107*> \verbatim
108*> B is REAL array, dimension (NMAX*NRHS)
109*> \endverbatim
110*>
111*> \param[out] X
112*> \verbatim
113*> X is REAL array, dimension (NMAX*NRHS)
114*> \endverbatim
115*>
116*> \param[out] XACT
117*> \verbatim
118*> XACT is REAL array, dimension (NMAX*NRHS)
119*> \endverbatim
120*>
121*> \param[out] WORK
122*> \verbatim
123*> WORK is REAL array, dimension (NMAX*max(2,NRHS))
124*> \endverbatim
125*>
126*> \param[out] RWORK
127*> \verbatim
128*> RWORK is REAL array, dimension (NMAX+2*NRHS)
129*> \endverbatim
130*>
131*> \param[out] IWORK
132*> \verbatim
133*> IWORK is INTEGER array, dimension (2*NMAX)
134*> \endverbatim
135*>
136*> \param[in] NOUT
137*> \verbatim
138*> NOUT is INTEGER
139*> The unit number for output.
140*> \endverbatim
141*
142* Authors:
143* ========
144*
145*> \author Univ. of Tennessee
146*> \author Univ. of California Berkeley
147*> \author Univ. of Colorado Denver
148*> \author NAG Ltd.
149*
150*> \ingroup single_lin
151*
152* =====================================================================
153 SUBROUTINE sdrvsy_rk( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR,
154 $ NMAX, A, AFAC, E, AINV, B, X, XACT, WORK,
155 $ RWORK, IWORK, NOUT )
156*
157* -- LAPACK test routine --
158* -- LAPACK is a software package provided by Univ. of Tennessee, --
159* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
160*
161* .. Scalar Arguments ..
162 LOGICAL TSTERR
163 INTEGER NMAX, NN, NOUT, NRHS
164 REAL THRESH
165* ..
166* .. Array Arguments ..
167 LOGICAL DOTYPE( * )
168 INTEGER IWORK( * ), NVAL( * )
169 REAL A( * ), AFAC( * ), AINV( * ), B( * ), E( * ),
170 $ rwork( * ), work( * ), x( * ), xact( * )
171* ..
172*
173* =====================================================================
174*
175* .. Parameters ..
176 REAL ONE, ZERO
177 PARAMETER ( ONE = 1.0e+0, zero = 0.0e+0 )
178 INTEGER NTYPES, NTESTS
179 parameter( ntypes = 10, ntests = 3 )
180 INTEGER NFACT
181 parameter( nfact = 2 )
182* ..
183* .. Local Scalars ..
184 LOGICAL ZEROT
185 CHARACTER DIST, FACT, TYPE, UPLO, XTYPE
186 CHARACTER*3 PATH, MATPATH
187 INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
188 $ izero, j, k, kl, ku, lda, lwork, mode, n,
189 $ nb, nbmin, nerrs, nfail, nimat, nrun, nt
190 REAL AINVNM, ANORM, CNDNUM, RCONDC
191* ..
192* .. Local Arrays ..
193 CHARACTER FACTS( NFACT ), UPLOS( 2 )
194 INTEGER ISEED( 4 ), ISEEDY( 4 )
195 REAL RESULT( NTESTS )
196* ..
197* .. External Functions ..
198 REAL SLANSY
199 EXTERNAL SLANSY
200* ..
201* .. External Subroutines ..
202 EXTERNAL aladhd, alaerh, alasvm, serrvx, sget04, slacpy,
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 ) = 'Single precision'
229 path( 2: 3 ) = 'SK'
230*
231* Path to generate matrices
232*
233 matpath( 1: 1 ) = 'Single precision'
234 matpath( 2: 3 ) = 'SY'
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 serrvx( 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 SLATB4 for the matrix generator
289* based on the type of matrix to be generated.
290*
291 CALL slatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
292 $ mode, cndnum, dist )
293*
294* Generate a matrix with SLATMS.
295*
296 srnamt = 'SLATMS'
297 CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
298 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
299 $ info )
300*
301* Check error code from SLATMS and handle error.
302*
303 IF( info.NE.0 ) THEN
304 CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
305 $ -1, -1, imat, nfail, nerrs, nout )
306*
307* Skip all tests for this generated matrix
308*
309 GO TO 160
310 END IF
311*
312* For types 3-6, zero one or more rows and columns of the
313* matrix to test that INFO is returned correctly.
314*
315 IF( zerot ) THEN
316 IF( imat.EQ.3 ) THEN
317 izero = 1
318 ELSE IF( imat.EQ.4 ) THEN
319 izero = n
320 ELSE
321 izero = n / 2 + 1
322 END IF
323*
324 IF( imat.LT.6 ) THEN
325*
326* Set row and column IZERO to zero.
327*
328 IF( iuplo.EQ.1 ) THEN
329 ioff = ( izero-1 )*lda
330 DO 20 i = 1, izero - 1
331 a( ioff+i ) = zero
332 20 CONTINUE
333 ioff = ioff + izero
334 DO 30 i = izero, n
335 a( ioff ) = zero
336 ioff = ioff + lda
337 30 CONTINUE
338 ELSE
339 ioff = izero
340 DO 40 i = 1, izero - 1
341 a( ioff ) = zero
342 ioff = ioff + lda
343 40 CONTINUE
344 ioff = ioff - izero
345 DO 50 i = izero, n
346 a( ioff+i ) = zero
347 50 CONTINUE
348 END IF
349 ELSE
350 ioff = 0
351 IF( iuplo.EQ.1 ) THEN
352*
353* Set the first IZERO rows and columns to zero.
354*
355 DO 70 j = 1, n
356 i2 = min( j, izero )
357 DO 60 i = 1, i2
358 a( ioff+i ) = zero
359 60 CONTINUE
360 ioff = ioff + lda
361 70 CONTINUE
362 ELSE
363*
364* Set the last IZERO rows and columns to zero.
365*
366 DO 90 j = 1, n
367 i1 = max( j, izero )
368 DO 80 i = i1, n
369 a( ioff+i ) = zero
370 80 CONTINUE
371 ioff = ioff + lda
372 90 CONTINUE
373 END IF
374 END IF
375 ELSE
376 izero = 0
377 END IF
378*
379* End generate the test matrix A.
380*
381 DO 150 ifact = 1, nfact
382*
383* Do first for FACT = 'F', then for other values.
384*
385 fact = facts( ifact )
386*
387* Compute the condition number
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_rk( uplo, n, afac, lda, e, 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*
411* We need to compute the inverse to compute
412* RCONDC that is used later in TEST3.
413*
414 CALL ssytri_3( uplo, n, ainv, lda, e, iwork,
415 $ work, lwork, info )
416 ainvnm = slansy( '1', uplo, n, ainv, lda, rwork )
417*
418* Compute the 1-norm condition number of A.
419*
420 IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
421 rcondc = one
422 ELSE
423 rcondc = ( one / anorm ) / ainvnm
424 END IF
425 END IF
426*
427* Form an exact solution and set the right hand side.
428*
429 srnamt = 'SLARHS'
430 CALL slarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
431 $ nrhs, a, lda, xact, lda, b, lda, iseed,
432 $ info )
433 xtype = 'C'
434*
435* --- Test SSYSV_RK ---
436*
437 IF( ifact.EQ.2 ) THEN
438 CALL slacpy( uplo, n, n, a, lda, afac, lda )
439 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
440*
441* Factor the matrix and solve the system using
442* SSYSV_RK.
443*
444 srnamt = 'SSYSV_RK'
445 CALL ssysv_rk( uplo, n, nrhs, afac, lda, e, iwork,
446 $ x, lda, work, lwork, info )
447*
448* Adjust the expected value of INFO to account for
449* pivoting.
450*
451 k = izero
452 IF( k.GT.0 ) THEN
453 100 CONTINUE
454 IF( iwork( k ).LT.0 ) THEN
455 IF( iwork( k ).NE.-k ) THEN
456 k = -iwork( k )
457 GO TO 100
458 END IF
459 ELSE IF( iwork( k ).NE.k ) THEN
460 k = iwork( k )
461 GO TO 100
462 END IF
463 END IF
464*
465* Check error code from SSYSV_RK and handle error.
466*
467 IF( info.NE.k ) THEN
468 CALL alaerh( path, 'SSYSV_RK', info, k, uplo,
469 $ n, n, -1, -1, nrhs, imat, nfail,
470 $ nerrs, nout )
471 GO TO 120
472 ELSE IF( info.NE.0 ) THEN
473 GO TO 120
474 END IF
475*
476*+ TEST 1 Reconstruct matrix from factors and compute
477* residual.
478*
479 CALL ssyt01_3( uplo, n, a, lda, afac, lda, e,
480 $ iwork, ainv, lda, rwork,
481 $ result( 1 ) )
482*
483*+ TEST 2 Compute residual of the computed solution.
484*
485 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
486 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
487 $ lda, rwork, result( 2 ) )
488*
489*+ TEST 3
490* Check solution from generated exact solution.
491*
492 CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
493 $ result( 3 ) )
494 nt = 3
495*
496* Print information about the tests that did not pass
497* the threshold.
498*
499 DO 110 k = 1, nt
500 IF( result( k ).GE.thresh ) THEN
501 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
502 $ CALL aladhd( nout, path )
503 WRITE( nout, fmt = 9999 )'SSYSV_RK', uplo,
504 $ n, imat, k, result( k )
505 nfail = nfail + 1
506 END IF
507 110 CONTINUE
508 nrun = nrun + nt
509 120 CONTINUE
510 END IF
511*
512 150 CONTINUE
513*
514 160 CONTINUE
515 170 CONTINUE
516 180 CONTINUE
517*
518* Print a summary of the results.
519*
520 CALL alasvm( path, nout, nfail, nrun, nerrs )
521*
522 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
523 $ ', test ', i2, ', ratio =', g12.5 )
524 RETURN
525*
526* End of SDRVSY_RK
527*
528 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_rk(uplo, n, nrhs, a, lda, e, ipiv, b, ldb, work, lwork, info)
SSYSV_RK computes the solution to system of linear equations A * X = B for SY matrices
Definition ssysv_rk.f:228
subroutine ssytrf_rk(uplo, n, a, lda, e, ipiv, work, lwork, info)
SSYTRF_RK computes the factorization of a real symmetric indefinite matrix using the bounded Bunch-Ka...
Definition ssytrf_rk.f:259
subroutine ssytri_3(uplo, n, a, lda, e, ipiv, work, lwork, info)
SSYTRI_3
Definition ssytri_3.f:170
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 sdrvsy_rk(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, e, ainv, b, x, xact, work, rwork, iwork, nout)
SDRVSY_RK
Definition sdrvsy_rk.f:156
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 ssyt01_3(uplo, n, a, lda, afac, ldafac, e, ipiv, c, ldc, rwork, resid)
SSYT01_3
Definition ssyt01_3.f:140