LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zdrvsy_aa.f
Go to the documentation of this file.
1*> \brief \b ZDRVSY_AA
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 ZDRVSY_AA( DOTYPE, NN, NVAL, NRHS, THRESH, TSTERR, NMAX,
12* A, AFAC, AINV, B, X, XACT, WORK, RWORK, IWORK,
13* 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*> ZDRVSY_AA tests the driver routine ZSYSV_AA.
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 COMPLEX*16
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 COMPLEX*16 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 complex16_lin
148*
149* =====================================================================
150 SUBROUTINE zdrvsy_aa( 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 ZERO
175 PARAMETER ( ZERO = 0.0d+0 )
176 COMPLEX*16 CZERO
177 parameter( czero = 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 MATPATH, PATH
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 DOUBLE PRECISION ANORM, CNDNUM
191* ..
192* .. Local Arrays ..
193 CHARACTER FACTS( NFACT ), UPLOS( 2 )
194 INTEGER ISEED( 4 ), ISEEDY( 4 )
195 DOUBLE PRECISION RESULT( NTESTS )
196* ..
197* .. External Functions ..
198 DOUBLE PRECISION DGET06, ZLANSY
199 EXTERNAL DGET06, ZLANSY
200* ..
201* .. External Subroutines ..
202 EXTERNAL aladhd, alaerh, alasvm, zerrvx, zget04, zlacpy,
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 ) = 'SA'
230*
231* Path to generate matrices
232*
233 matpath( 1: 1 ) = 'Zomplex 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*
243* Test the error exits
244*
245 IF( tsterr )
246 $ CALL zerrvx( path, nout )
247 infot = 0
248*
249* Set the block size and minimum block size for testing.
250*
251 nb = 1
252 nbmin = 2
253 CALL xlaenv( 1, nb )
254 CALL xlaenv( 2, nbmin )
255*
256* Do for each value of N in NVAL
257*
258 DO 180 in = 1, nn
259 n = nval( in )
260 lwork = max( 3*n-2, n*(1+nb) )
261 lwork = max( lwork, 1 )
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* Set up parameters with ZLATB4 and generate a test matrix
287* with ZLATMS.
288*
289 CALL zlatb4( matpath, imat, n, n, TYPE, kl, ku, anorm,
290 $ mode, cndnum, dist )
291*
292 srnamt = 'ZLATMS'
293 CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
294 $ cndnum, anorm, kl, ku, uplo, a, lda, work,
295 $ info )
296*
297* Check error code from ZLATMS.
298*
299 IF( info.NE.0 ) THEN
300 CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
301 $ -1, -1, imat, nfail, nerrs, nout )
302 GO TO 160
303 END IF
304*
305* For types 3-6, zero one or more rows and columns of the
306* matrix to test that INFO is returned correctly.
307*
308 IF( zerot ) THEN
309 IF( imat.EQ.3 ) THEN
310 izero = 1
311 ELSE IF( imat.EQ.4 ) THEN
312 izero = n
313 ELSE
314 izero = n / 2 + 1
315 END IF
316*
317 IF( imat.LT.6 ) THEN
318*
319* Set row and column IZERO to zero.
320*
321 IF( iuplo.EQ.1 ) THEN
322 ioff = ( izero-1 )*lda
323 DO 20 i = 1, izero - 1
324 a( ioff+i ) = czero
325 20 CONTINUE
326 ioff = ioff + izero
327 DO 30 i = izero, n
328 a( ioff ) = czero
329 ioff = ioff + lda
330 30 CONTINUE
331 ELSE
332 ioff = izero
333 DO 40 i = 1, izero - 1
334 a( ioff ) = czero
335 ioff = ioff + lda
336 40 CONTINUE
337 ioff = ioff - izero
338 DO 50 i = izero, n
339 a( ioff+i ) = czero
340 50 CONTINUE
341 END IF
342 ELSE
343 ioff = 0
344 IF( iuplo.EQ.1 ) THEN
345*
346* Set the first IZERO rows and columns to zero.
347*
348 DO 70 j = 1, n
349 i2 = min( j, izero )
350 DO 60 i = 1, i2
351 a( ioff+i ) = czero
352 60 CONTINUE
353 ioff = ioff + lda
354 70 CONTINUE
355 izero = 1
356 ELSE
357*
358* Set the last IZERO rows and columns to zero.
359*
360 DO 90 j = 1, n
361 i1 = max( j, izero )
362 DO 80 i = i1, n
363 a( ioff+i ) = czero
364 80 CONTINUE
365 ioff = ioff + lda
366 90 CONTINUE
367 END IF
368 END IF
369 ELSE
370 izero = 0
371 END IF
372*
373 DO 150 ifact = 1, nfact
374*
375* Do first for FACT = 'F', then for other values.
376*
377 fact = facts( ifact )
378*
379* Form an exact solution and set the right hand side.
380*
381 srnamt = 'ZLARHS'
382 CALL zlarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
383 $ nrhs, a, lda, xact, lda, b, lda, iseed,
384 $ info )
385 xtype = 'C'
386*
387* --- Test ZSYSV_AA ---
388*
389 IF( ifact.EQ.2 ) THEN
390 CALL zlacpy( uplo, n, n, a, lda, afac, lda )
391 CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
392*
393* Factor the matrix and solve the system using ZSYSV_AA.
394*
395 srnamt = 'ZSYSV_AA'
396 CALL zsysv_aa( uplo, n, nrhs, afac, lda, iwork,
397 $ x, lda, work, lwork, info )
398*
399* Adjust the expected value of INFO to account for
400* pivoting.
401*
402 IF( izero.GT.0 ) THEN
403 j = 1
404 k = izero
405 100 CONTINUE
406 IF( j.EQ.k ) THEN
407 k = iwork( j )
408 ELSE IF( iwork( j ).EQ.k ) THEN
409 k = j
410 END IF
411 IF( j.LT.k ) THEN
412 j = j + 1
413 GO TO 100
414 END IF
415 ELSE
416 k = 0
417 END IF
418*
419* Check error code from ZSYSV_AA .
420*
421 IF( info.NE.k ) THEN
422 CALL alaerh( path, 'ZSYSV_AA ', info, k,
423 $ uplo, n, n, -1, -1, nrhs,
424 $ imat, nfail, nerrs, nout )
425 GO TO 120
426 ELSE IF( info.NE.0 ) THEN
427 GO TO 120
428 END IF
429*
430* Reconstruct matrix from factors and compute
431* residual.
432*
433 CALL zsyt01_aa( uplo, n, a, lda, afac, lda,
434 $ iwork, ainv, lda, rwork,
435 $ result( 1 ) )
436*
437* Compute residual of the computed solution.
438*
439 CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
440 CALL zsyt02( uplo, n, nrhs, a, lda, x, lda, work,
441 $ lda, rwork, result( 2 ) )
442 nt = 2
443*
444* Print information about the tests that did not pass
445* the threshold.
446*
447 DO 110 k = 1, nt
448 IF( result( k ).GE.thresh ) THEN
449 IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
450 $ CALL aladhd( nout, path )
451 WRITE( nout, fmt = 9999 )'ZSYSV_AA ',
452 $ uplo, n, imat, k, result( k )
453 nfail = nfail + 1
454 END IF
455 110 CONTINUE
456 nrun = nrun + nt
457 120 CONTINUE
458 END IF
459*
460 150 CONTINUE
461*
462 160 CONTINUE
463 170 CONTINUE
464 180 CONTINUE
465*
466* Print a summary of the results.
467*
468 CALL alasvm( path, nout, nfail, nrun, nerrs )
469*
470 9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
471 $ ', test ', i2, ', ratio =', g12.5 )
472 RETURN
473*
474* End of ZDRVSY_AA
475*
476 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 zsysv_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
ZSYSV_AA computes the solution to system of linear equations A * X = B for SY matrices
Definition zsysv_aa.f:162
subroutine zsytrf_aa(uplo, n, a, lda, ipiv, work, lwork, info)
ZSYTRF_AA
Definition zsytrf_aa.f:132
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 zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
subroutine zdrvsy_aa(dotype, nn, nval, nrhs, thresh, tsterr, nmax, a, afac, ainv, b, x, xact, work, rwork, iwork, nout)
ZDRVSY_AA
Definition zdrvsy_aa.f:153
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 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 zsyt01_aa(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
ZSYT01
Definition zsyt01_aa.f:124
subroutine zsyt02(uplo, n, nrhs, a, lda, x, ldx, b, ldb, rwork, resid)
ZSYT02
Definition zsyt02.f:127