LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ zdrvsy_aa_2stage()

subroutine zdrvsy_aa_2stage ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
double precision  THRESH,
logical  TSTERR,
integer  NMAX,
complex*16, dimension( * )  A,
complex*16, dimension( * )  AFAC,
complex*16, dimension( * )  AINV,
complex*16, dimension( * )  B,
complex*16, dimension( * )  X,
complex*16, dimension( * )  XACT,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

ZDRVSY_AA_2STAGE

Purpose:
 ZDRVSY_AA_2STAGE tests the driver routine ZSYSV_AA_2STAGE.
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 DOUBLE PRECISION
          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 COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]AINV
          AINV is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]X
          X is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]XACT
          XACT is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]WORK
          WORK is COMPLEX*16 array, dimension (NMAX*max(2,NRHS))
[out]RWORK
          RWORK is COMPLEX*16 array, dimension (NMAX+2*NRHS)
[out]IWORK
          IWORK is INTEGER array, dimension (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.
Date
November 2017

Definition at line 157 of file zdrvsy_aa_2stage.f.

157 *
158 * -- LAPACK test routine (version 3.8.0) --
159 * -- LAPACK is a software package provided by Univ. of Tennessee, --
160 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
161 * November 2017
162 *
163 * .. Scalar Arguments ..
164  LOGICAL TSTERR
165  INTEGER NMAX, NN, NOUT, NRHS
166  DOUBLE PRECISION THRESH
167 * ..
168 * .. Array Arguments ..
169  LOGICAL DOTYPE( * )
170  INTEGER IWORK( * ), NVAL( * )
171  DOUBLE PRECISION RWORK( * )
172  COMPLEX*16 A( * ), AFAC( * ), AINV( * ), B( * ),
173  $ WORK( * ), X( * ), XACT( * )
174 * ..
175 *
176 * =====================================================================
177 *
178 * .. Parameters ..
179  DOUBLE PRECISION ZERO
180  parameter( zero = 0.0d+0 )
181  COMPLEX*16 CZERO
182  parameter( czero = ( 0.0d+0, 0.0d+0 ) )
183  INTEGER NTYPES, NTESTS
184  parameter( ntypes = 10, ntests = 3 )
185  INTEGER NFACT
186  parameter( nfact = 2 )
187 * ..
188 * .. Local Scalars ..
189  LOGICAL ZEROT
190  CHARACTER DIST, FACT, TYPE, UPLO, XTYPE
191  CHARACTER*3 MATPATH, PATH
192  INTEGER I, I1, I2, IFACT, IMAT, IN, INFO, IOFF, IUPLO,
193  $ IZERO, J, K, KL, KU, LDA, LWORK, MODE, N,
194  $ NB, NBMIN, NERRS, NFAIL, NIMAT, NRUN, NT
195  DOUBLE PRECISION ANORM, CNDNUM
196 * ..
197 * .. Local Arrays ..
198  CHARACTER FACTS( NFACT ), UPLOS( 2 )
199  INTEGER ISEED( 4 ), ISEEDY( 4 )
200  DOUBLE PRECISION RESULT( NTESTS )
201 * ..
202 * .. External Functions ..
203  DOUBLE PRECISION DGET06, ZLANSY
204  EXTERNAL dget06, zlansy
205 * ..
206 * .. External Subroutines ..
207  EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx,
211 * ..
212 * .. Scalars in Common ..
213  LOGICAL LERR, OK
214  CHARACTER*32 SRNAMT
215  INTEGER INFOT, NUNIT
216 * ..
217 * .. Common blocks ..
218  COMMON / infoc / infot, nunit, ok, lerr
219  COMMON / srnamc / srnamt
220 * ..
221 * .. Intrinsic Functions ..
222  INTRINSIC cmplx, max, min
223 * ..
224 * .. Data statements ..
225  DATA iseedy / 1988, 1989, 1990, 1991 /
226  DATA uplos / 'U', 'L' / , facts / 'F', 'N' /
227 * ..
228 * .. Executable Statements ..
229 *
230 * Initialize constants and the random number seed.
231 *
232 * Test path
233 *
234  path( 1: 1 ) = 'Zomplex precision'
235  path( 2: 3 ) = 'H2'
236 *
237 * Path to generate matrices
238 *
239  matpath( 1: 1 ) = 'Zomplex precision'
240  matpath( 2: 3 ) = 'SY'
241 *
242  nrun = 0
243  nfail = 0
244  nerrs = 0
245  DO 10 i = 1, 4
246  iseed( i ) = iseedy( i )
247  10 CONTINUE
248 *
249 * Test the error exits
250 *
251  IF( tsterr )
252  $ CALL zerrvx( path, nout )
253  infot = 0
254 *
255 * Set the block size and minimum block size for testing.
256 *
257  nb = 1
258  nbmin = 2
259  CALL xlaenv( 1, nb )
260  CALL xlaenv( 2, nbmin )
261 *
262 * Do for each value of N in NVAL
263 *
264  DO 180 in = 1, nn
265  n = nval( in )
266  lda = max( n, 1 )
267  xtype = 'N'
268  nimat = ntypes
269  IF( n.LE.0 )
270  $ nimat = 1
271 *
272  DO 170 imat = 1, nimat
273 *
274 * Do the tests only if DOTYPE( IMAT ) is true.
275 *
276  IF( .NOT.dotype( imat ) )
277  $ GO TO 170
278 *
279 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
280 *
281  zerot = imat.GE.3 .AND. imat.LE.6
282  IF( zerot .AND. n.LT.imat-2 )
283  $ GO TO 170
284 *
285 * Do first for UPLO = 'U', then for UPLO = 'L'
286 *
287  DO 160 iuplo = 1, 2
288  uplo = uplos( iuplo )
289 *
290 * Begin generate the test matrix A.
291 *
292 * Set up parameters with ZLATB4 for the matrix generator
293 * based on the type of matrix to be generated.
294 *
295  CALL zlatb4( matpath, imat, n, n, TYPE, KL, KU, ANORM,
296  $ MODE, CNDNUM, DIST )
297 *
298 * Generate a matrix with ZLATMS.
299 *
300  srnamt = 'ZLATMS'
301  CALL zlatms( n, n, dist, iseed, TYPE, RWORK, MODE,
302  $ CNDNUM, ANORM, KL, KU, UPLO, A, LDA,
303  $ WORK, INFO )
304 *
305 * Check error code from ZLATMS and handle error.
306 *
307  IF( info.NE.0 ) THEN
308  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n,
309  $ -1, -1, -1, imat, nfail, nerrs, nout )
310  GO TO 160
311  END IF
312 *
313 * For types 3-6, zero one or more rows and columns of
314 * the matrix to test that INFO is returned correctly.
315 *
316  IF( zerot ) THEN
317  IF( imat.EQ.3 ) THEN
318  izero = 1
319  ELSE IF( imat.EQ.4 ) THEN
320  izero = n
321  ELSE
322  izero = n / 2 + 1
323  END IF
324 *
325  IF( imat.LT.6 ) THEN
326 *
327 * Set row and column IZERO to zero.
328 *
329  IF( iuplo.EQ.1 ) THEN
330  ioff = ( izero-1 )*lda
331  DO 20 i = 1, izero - 1
332  a( ioff+i ) = czero
333  20 CONTINUE
334  ioff = ioff + izero
335  DO 30 i = izero, n
336  a( ioff ) = czero
337  ioff = ioff + lda
338  30 CONTINUE
339  ELSE
340  ioff = izero
341  DO 40 i = 1, izero - 1
342  a( ioff ) = czero
343  ioff = ioff + lda
344  40 CONTINUE
345  ioff = ioff - izero
346  DO 50 i = izero, n
347  a( ioff+i ) = czero
348  50 CONTINUE
349  END IF
350  ELSE
351  ioff = 0
352  IF( iuplo.EQ.1 ) THEN
353 *
354 * Set the first IZERO rows and columns to zero.
355 *
356  DO 70 j = 1, n
357  i2 = min( j, izero )
358  DO 60 i = 1, i2
359  a( ioff+i ) = czero
360  60 CONTINUE
361  ioff = ioff + lda
362  70 CONTINUE
363  izero = 1
364  ELSE
365 *
366 * Set the first IZERO rows and columns to zero.
367 *
368  ioff = 0
369  DO 90 j = 1, n
370  i1 = max( j, izero )
371  DO 80 i = i1, n
372  a( ioff+i ) = czero
373  80 CONTINUE
374  ioff = ioff + lda
375  90 CONTINUE
376  END IF
377  END IF
378  ELSE
379  izero = 0
380  END IF
381 *
382 * End generate the test matrix A.
383 *
384 *
385  DO 150 ifact = 1, nfact
386 *
387 * Do first for FACT = 'F', then for other values.
388 *
389  fact = facts( ifact )
390 *
391 * Form an exact solution and set the right hand side.
392 *
393  srnamt = 'ZLARHS'
394  CALL zlarhs( matpath, xtype, uplo, ' ', n, n, kl, ku,
395  $ nrhs, a, lda, xact, lda, b, lda, iseed,
396  $ info )
397  xtype = 'C'
398 *
399 * --- Test ZSYSV_AA_2STAGE ---
400 *
401  IF( ifact.EQ.2 ) THEN
402  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
403  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
404 *
405 * Factor the matrix and solve the system using ZSYSV_AA.
406 *
407  srnamt = 'ZSYSV_AA_2STAGE '
408  lwork = min(n*nb, 3*nmax*nmax)
409  CALL zsysv_aa_2stage( uplo, n, nrhs, afac, lda,
410  $ ainv, (3*nb+1)*n,
411  $ iwork, iwork( 1+n ),
412  $ x, lda, work, lwork, info )
413 *
414 * Adjust the expected value of INFO to account for
415 * pivoting.
416 *
417  IF( izero.GT.0 ) THEN
418  j = 1
419  k = izero
420  100 CONTINUE
421  IF( j.EQ.k ) THEN
422  k = iwork( j )
423  ELSE IF( iwork( j ).EQ.k ) THEN
424  k = j
425  END IF
426  IF( j.LT.k ) THEN
427  j = j + 1
428  GO TO 100
429  END IF
430  ELSE
431  k = 0
432  END IF
433 *
434 * Check error code from ZSYSV_AA_2STAGE .
435 *
436  IF( info.NE.k ) THEN
437  CALL alaerh( path, 'ZSYSV_AA_2STAGE', info, k,
438  $ uplo, n, n, -1, -1, nrhs,
439  $ imat, nfail, nerrs, nout )
440  GO TO 120
441  ELSE IF( info.NE.0 ) THEN
442  GO TO 120
443  END IF
444 *
445 * Compute residual of the computed solution.
446 *
447  CALL zlacpy( 'Full', n, nrhs, b, lda, work, lda )
448  CALL zsyt02( uplo, n, nrhs, a, lda, x, lda, work,
449  $ lda, rwork, result( 1 ) )
450 *
451 * Reconstruct matrix from factors and compute
452 * residual.
453 *
454 c CALL ZSY01_AA( UPLO, N, A, LDA, AFAC, LDA,
455 c $ IWORK, AINV, LDA, RWORK,
456 c $ RESULT( 2 ) )
457 c NT = 2
458  nt = 1
459 *
460 * Print information about the tests that did not pass
461 * the threshold.
462 *
463  DO 110 k = 1, nt
464  IF( result( k ).GE.thresh ) THEN
465  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
466  $ CALL aladhd( nout, path )
467  WRITE( nout, fmt = 9999 )'ZSYSV_AA_2STAGE ',
468  $ uplo, n, imat, k, result( k )
469  nfail = nfail + 1
470  END IF
471  110 CONTINUE
472  nrun = nrun + nt
473  120 CONTINUE
474  END IF
475 *
476  150 CONTINUE
477 *
478  160 CONTINUE
479  170 CONTINUE
480  180 CONTINUE
481 *
482 * Print a summary of the results.
483 *
484  CALL alasvm( path, nout, nfail, nrun, nerrs )
485 *
486  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
487  $ ', test ', i2, ', ratio =', g12.5 )
488  RETURN
489 *
490 * End of ZDRVSY_AA_2STAGE
491 *
Here is the call graph for this function:
Here is the caller graph for this function:
dget06
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
zerrvx
subroutine zerrvx(PATH, NUNIT)
ZERRVX
Definition: zerrvx.f:57
zlansy
double precision function zlansy(NORM, UPLO, N, A, LDA, WORK)
ZLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition: zlansy.f:125
alasvm
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
zlarhs
subroutine zlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
ZLARHS
Definition: zlarhs.f:211
zsyt01_aa
subroutine zsyt01_aa(UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C, LDC, RWORK, RESID)
ZSYT01
Definition: zsyt01_aa.f:128
zget04
subroutine zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
zlacpy
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
zsyt02
subroutine zsyt02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZSYT02
Definition: zsyt02.f:129
aladhd
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:92
alaerh
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
zsytrf_aa_2stage
subroutine zsytrf_aa_2stage(UPLO, N, A, LDA, TB, LTB, IPIV, IPIV2, WORK, LWORK, INFO)
ZSYTRF_AA_2STAGE
Definition: zsytrf_aa_2stage.f:162
zsysv_aa_2stage
subroutine zsysv_aa_2stage(UPLO, N, NRHS, A, LDA, TB, LTB, IPIV, IPIV2, B, LDB, WORK, LWORK, INFO)
ZSYSV_AA_2STAGE computes the solution to system of linear equations A * X = B for SY matrices
Definition: zsysv_aa_2stage.f:188
zlatb4
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
xlaenv
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
zlatms
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334