LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sdrvsy_aa()

subroutine sdrvsy_aa ( logical, dimension( * )  dotype,
integer  nn,
integer, dimension( * )  nval,
integer  nrhs,
real  thresh,
logical  tsterr,
integer  nmax,
real, dimension( * )  a,
real, dimension( * )  afac,
real, dimension( * )  ainv,
real, dimension( * )  b,
real, dimension( * )  x,
real, dimension( * )  xact,
real, dimension( * )  work,
real, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer  nout 
)

SDRVSY_AA

Purpose:
 SDRVSY_AA tests the driver routine SSYSV_AA.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix dimension N.
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[in]THRESH
          THRESH is REAL
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]AINV
          AINV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NRHS)
[out]X
          X is REAL array, dimension (NMAX*NRHS)
[out]XACT
          XACT is REAL array, dimension (NMAX*NRHS)
[out]WORK
          WORK is REAL array, dimension (NMAX*max(2,NRHS))
[out]RWORK
          RWORK is REAL array, dimension (NMAX+2*NRHS)
[out]IWORK
          IWORK is INTEGER array, dimension (2*NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 149 of file sdrvsy_aa.f.

152*
153* -- LAPACK test routine --
154* -- LAPACK is a software package provided by Univ. of Tennessee, --
155* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
156*
157* .. Scalar Arguments ..
158 LOGICAL TSTERR
159 INTEGER NMAX, NN, NOUT, NRHS
160 REAL THRESH
161* ..
162* .. Array Arguments ..
163 LOGICAL DOTYPE( * )
164 INTEGER IWORK( * ), NVAL( * )
165 REAL A( * ), AFAC( * ), AINV( * ), B( * ),
166 $ RWORK( * ), WORK( * ), X( * ), XACT( * )
167* ..
168*
169* =====================================================================
170*
171* .. Parameters ..
172 REAL ONE, ZERO
173 parameter( one = 1.0e+0, zero = 0.0e+0 )
174 INTEGER NTYPES, NTESTS
175 parameter( ntypes = 10, ntests = 3 )
176 INTEGER NFACT
177 parameter( nfact = 2 )
178* ..
179* .. Local Scalars ..
180 LOGICAL ZEROT
181 CHARACTER DIST, FACT, TYPE, UPLO, XTYPE
182 CHARACTER*3 MATPATH, PATH
183 INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
184 $ IZERO, J, K, KL, KU, LDA, LWORK, MODE, N,
185 $ NB, NBMIN, NERRS, NFAIL, NIMAT, NRUN, NT
186 REAL ANORM, CNDNUM
187* ..
188* .. Local Arrays ..
189 CHARACTER FACTS( NFACT ), UPLOS( 2 )
190 INTEGER ISEED( 4 ), ISEEDY( 4 )
191 REAL RESULT( NTESTS )
192* ..
193* .. External Functions ..
194 REAL DGET06, SLANSY
195 EXTERNAL dget06, slansy
196* ..
197* .. External Subroutines ..
198 EXTERNAL aladhd, alaerh, alasvm, serrvx, sget04, slacpy,
201* ..
202* .. Scalars in Common ..
203 LOGICAL LERR, OK
204 CHARACTER*32 SRNAMT
205 INTEGER INFOT, NUNIT
206* ..
207* .. Common blocks ..
208 COMMON / infoc / infot, nunit, ok, lerr
209 COMMON / srnamc / srnamt
210* ..
211* .. Intrinsic Functions ..
212 INTRINSIC max, min
213* ..
214* .. Data statements ..
215 DATA iseedy / 1988, 1989, 1990, 1991 /
216 DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
217* ..
218* .. Executable Statements ..
219*
220* Initialize constants and the random number seed.
221*
222* Test path
223*
224 path( 1: 1 ) = 'Single precision'
225 path( 2: 3 ) = 'SA'
226*
227* Path to generate matrices
228*
229 matpath( 1: 1 ) = 'Single precision'
230 matpath( 2: 3 ) = 'SY'
231*
232 nrun = 0
233 nfail = 0
234 nerrs = 0
235 DO 10 i = 1, 4
236 iseed( i ) = iseedy( i )
237 10 CONTINUE
238*
239* Test the error exits
240*
241 IF( tsterr )
242 $ CALL serrvx( path, nout )
243 infot = 0
244*
245* Set the block size and minimum block size for testing.
246*
247 nb = 1
248 nbmin = 2
249 CALL xlaenv( 1, nb )
250 CALL xlaenv( 2, nbmin )
251*
252* Do for each value of N in NVAL
253*
254 DO 180 in = 1, nn
255 n = nval( in )
256 lwork = max( 3*n-2, n*(1+nb) )
257 lwork = max( lwork, 1 )
258 lda = max( n, 1 )
259 xtype = 'N'
260 nimat = ntypes
261 IF( n.LE.0 )
262 $ nimat = 1
263*
264 DO 170 imat = 1, nimat
265*
266* Do the tests only if DOTYPE( IMAT ) is true.
267*
268 IF( .NOT.dotype( imat ) )
269 $ GO TO 170
270*
271* Skip types 3, 4, 5, or 6 if the matrix size is too small.
272*
273 zerot = imat.GE.3 .AND. imat.LE.6
274 IF( zerot .AND. n.LT.imat-2 )
275 $ GO TO 170
276*
277* Do first for UPLO = 'U', then for UPLO = 'L'
278*
279 DO 160 iuplo = 1, 2
280 uplo = uplos( iuplo )
281*
282* Set up parameters with SLATB4 and generate a test matrix
283* with SLATMS.
284*
285 CALL slatb4( matpath, imat, n, n, TYPE, KL, KU, ANORM,
286 $ MODE, CNDNUM, DIST )
287*
288 srnamt = 'SLATMS'
289 CALL slatms( n, n, dist, iseed, TYPE, RWORK, MODE,
290 $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA, WORK,
291 $ INFO )
292*
293* Check error code from SLATMS.
294*
295 IF( info.NE.0 ) THEN
296 CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
297 $ -1, -1, imat, nfail, nerrs, nout )
298 GO TO 160
299 END IF
300*
301* For types 3-6, zero one or more rows and columns of the
302* matrix to test that INFO is returned correctly.
303*
304 IF( zerot ) THEN
305 IF( imat.EQ.3 ) THEN
306 izero = 1
307 ELSE IF( imat.EQ.4 ) THEN
308 izero = n
309 ELSE
310 izero = n / 2 + 1
311 END IF
312*
313 IF( imat.LT.6 ) THEN
314*
315* Set row and column IZERO to zero.
316*
317 IF( iuplo.EQ.1 ) THEN
318 ioff = ( izero-1 )*lda
319 DO 20 i = 1, izero - 1
320 a( ioff+i ) = zero
321 20 CONTINUE
322 ioff = ioff + izero
323 DO 30 i = izero, n
324 a( ioff ) = zero
325 ioff = ioff + lda
326 30 CONTINUE
327 ELSE
328 ioff = izero
329 DO 40 i = 1, izero - 1
330 a( ioff ) = zero
331 ioff = ioff + lda
332 40 CONTINUE
333 ioff = ioff - izero
334 DO 50 i = izero, n
335 a( ioff+i ) = zero
336 50 CONTINUE
337 END IF
338 ELSE
339 ioff = 0
340 IF( iuplo.EQ.1 ) THEN
341*
342* Set the first IZERO rows and columns to zero.
343*
344 DO 70 j = 1, n
345 i2 = min( j, izero )
346 DO 60 i = 1, i2
347 a( ioff+i ) = zero
348 60 CONTINUE
349 ioff = ioff + lda
350 70 CONTINUE
351 izero = 1
352 ELSE
353*
354* Set the last IZERO rows and columns to zero.
355*
356 DO 90 j = 1, n
357 i1 = max( j, izero )
358 DO 80 i = i1, n
359 a( ioff+i ) = zero
360 80 CONTINUE
361 ioff = ioff + lda
362 90 CONTINUE
363 END IF
364 END IF
365 ELSE
366 izero = 0
367 END IF
368*
369 DO 150 ifact = 1, nfact
370*
371* Do first for FACT = 'F', then for other values.
372*
373 fact = facts( ifact )
374*
375* Form an exact solution and set the right hand side.
376*
377 srnamt = 'SLARHS'
378 CALL slarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
379 $ nrhs, a, lda, xact, lda, b, lda, iseed,
380 $ info )
381 xtype = 'C'
382*
383* --- Test SSYSV_AA ---
384*
385 IF( ifact.EQ.2 ) THEN
386 CALL slacpy( uplo, n, n, a, lda, afac, lda )
387 CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
388*
389* Factor the matrix and solve the system using SSYSV_AA.
390*
391 srnamt = 'SSYSV_AA'
392 CALL ssysv_aa( uplo, n, nrhs, afac, lda, iwork,
393 $ x, lda, work, lwork, info )
394*
395* Adjust the expected value of INFO to account for
396* pivoting.
397*
398 IF( izero.GT.0 ) THEN
399 j = 1
400 k = izero
401 100 CONTINUE
402 IF( j.EQ.k ) THEN
403 k = iwork( j )
404 ELSE IF( iwork( j ).EQ.k ) THEN
405 k = j
406 END IF
407 IF( j.LT.k ) THEN
408 j = j + 1
409 GO TO 100
410 END IF
411 ELSE
412 k = 0
413 END IF
414*
415* Check error code from SSYSV_AA .
416*
417 IF( info.NE.k ) THEN
418 CALL alaerh( path, 'SSYSV_AA ', info, k,
419 $ uplo, n, n, -1, -1, nrhs,
420 $ imat, nfail, nerrs, nout )
421 GO TO 120
422 ELSE IF( info.NE.0 ) THEN
423 GO TO 120
424 END IF
425*
426* Reconstruct matrix from factors and compute
427* residual.
428*
429 CALL ssyt01_aa( uplo, n, a, lda, afac, lda,
430 $ iwork, ainv, lda, rwork,
431 $ result( 1 ) )
432*
433* Compute residual of the computed solution.
434*
435 CALL slacpy( 'Full', n, nrhs, b, lda, work, lda )
436 CALL spot02( uplo, n, nrhs, a, lda, x, lda, work,
437 $ lda, rwork, result( 2 ) )
438 nt = 2
439*
440* Print information about the tests that did not pass
441* the threshold.
442*
443 DO 110 k = 1, nt
444 IF( result( k ).GE.thresh ) THEN
445 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
446 $ CALL aladhd( nout, path )
447 WRITE( nout, fmt = 9999 )'SSYSV_AA ',
448 $ uplo, n, imat, k, result( k )
449 nfail = nfail + 1
450 END IF
451 110 CONTINUE
452 nrun = nrun + nt
453 120 CONTINUE
454 END IF
455*
456 150 CONTINUE
457*
458 160 CONTINUE
459 170 CONTINUE
460 180 CONTINUE
461*
462* Print a summary of the results.
463*
464 CALL alasvm( path, nout, nfail, nrun, nerrs )
465*
466 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
467 $ ', test ', i2, ', ratio =', g12.5 )
468 RETURN
469*
470* End of SDRVSY_AA
471*
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
double precision function dget06(rcond, rcondc)
DGET06
Definition dget06.f:55
subroutine ssysv_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
SSYSV_AA computes the solution to system of linear equations A * X = B for SY matrices
Definition ssysv_aa.f:162
subroutine ssytrf_aa(uplo, n, a, lda, ipiv, work, lwork, info)
SSYTRF_AA
Definition ssytrf_aa.f:132
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
real function slansy(norm, uplo, n, a, lda, work)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansy.f:122
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 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_aa(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
SSYT01_AA
Definition ssyt01_aa.f:124
Here is the call graph for this function:
Here is the caller graph for this function: