LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cdrvsp ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
real  THRESH,
logical  TSTERR,
integer  NMAX,
complex, dimension( * )  A,
complex, dimension( * )  AFAC,
complex, dimension( * )  AINV,
complex, dimension( * )  B,
complex, dimension( * )  X,
complex, dimension( * )  XACT,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

CDRVSP

Purpose:
 CDRVSP tests the driver routines CSPSV and -SVX.
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 COMPLEX array, dimension
                      (NMAX*(NMAX+1)/2)
[out]AFAC
          AFAC is COMPLEX array, dimension
                      (NMAX*(NMAX+1)/2)
[out]AINV
          AINV is COMPLEX array, dimension
                      (NMAX*(NMAX+1)/2)
[out]B
          B is COMPLEX array, dimension (NMAX*NRHS)
[out]X
          X is COMPLEX array, dimension (NMAX*NRHS)
[out]XACT
          XACT is COMPLEX array, dimension (NMAX*NRHS)
[out]WORK
          WORK is COMPLEX array, dimension
                      (NMAX*max(2,NRHS))
[out]RWORK
          RWORK is REAL 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 2011

Definition at line 159 of file cdrvsp.f.

159 *
160 * -- LAPACK test routine (version 3.4.0) --
161 * -- LAPACK is a software package provided by Univ. of Tennessee, --
162 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
163 * November 2011
164 *
165 * .. Scalar Arguments ..
166  LOGICAL tsterr
167  INTEGER nmax, nn, nout, nrhs
168  REAL thresh
169 * ..
170 * .. Array Arguments ..
171  LOGICAL dotype( * )
172  INTEGER iwork( * ), nval( * )
173  REAL rwork( * )
174  COMPLEX a( * ), afac( * ), ainv( * ), b( * ),
175  $ work( * ), x( * ), xact( * )
176 * ..
177 *
178 * =====================================================================
179 *
180 * .. Parameters ..
181  REAL one, zero
182  parameter ( one = 1.0e+0, zero = 0.0e+0 )
183  INTEGER ntypes, ntests
184  parameter ( ntypes = 11, ntests = 6 )
185  INTEGER nfact
186  parameter ( nfact = 2 )
187 * ..
188 * .. Local Scalars ..
189  LOGICAL zerot
190  CHARACTER dist, fact, packit, TYPE, uplo, xtype
191  CHARACTER*3 path
192  INTEGER i, i1, i2, ifact, imat, in, info, ioff, iuplo,
193  $ izero, j, k, k1, kl, ku, lda, mode, n, nb,
194  $ nbmin, nerrs, nfail, nimat, npp, nrun, nt
195  REAL ainvnm, anorm, cndnum, rcond, rcondc
196 * ..
197 * .. Local Arrays ..
198  CHARACTER facts( nfact )
199  INTEGER iseed( 4 ), iseedy( 4 )
200  REAL result( ntests )
201 * ..
202 * .. External Functions ..
203  REAL clansp, sget06
204  EXTERNAL clansp, sget06
205 * ..
206 * .. External Subroutines ..
207  EXTERNAL aladhd, alaerh, alasvm, ccopy, cerrvx, cget04,
210  $ csptri, xlaenv
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 facts / 'F', 'N' /
227 * ..
228 * .. Executable Statements ..
229 *
230 * Initialize constants and the random number seed.
231 *
232  path( 1: 1 ) = 'Complex precision'
233  path( 2: 3 ) = 'SP'
234  nrun = 0
235  nfail = 0
236  nerrs = 0
237  DO 10 i = 1, 4
238  iseed( i ) = iseedy( i )
239  10 CONTINUE
240 *
241 * Test the error exits
242 *
243  IF( tsterr )
244  $ CALL cerrvx( path, nout )
245  infot = 0
246 *
247 * Set the block size and minimum block size for testing.
248 *
249  nb = 1
250  nbmin = 2
251  CALL xlaenv( 1, nb )
252  CALL xlaenv( 2, nbmin )
253 *
254 * Do for each value of N in NVAL
255 *
256  DO 180 in = 1, nn
257  n = nval( in )
258  lda = max( n, 1 )
259  npp = n*( n+1 ) / 2
260  xtype = 'N'
261  nimat = ntypes
262  IF( n.LE.0 )
263  $ nimat = 1
264 *
265  DO 170 imat = 1, nimat
266 *
267 * Do the tests only if DOTYPE( IMAT ) is true.
268 *
269  IF( .NOT.dotype( imat ) )
270  $ GO TO 170
271 *
272 * Skip types 3, 4, 5, or 6 if the matrix size is too small.
273 *
274  zerot = imat.GE.3 .AND. imat.LE.6
275  IF( zerot .AND. n.LT.imat-2 )
276  $ GO TO 170
277 *
278 * Do first for UPLO = 'U', then for UPLO = 'L'
279 *
280  DO 160 iuplo = 1, 2
281  IF( iuplo.EQ.1 ) THEN
282  uplo = 'U'
283  packit = 'C'
284  ELSE
285  uplo = 'L'
286  packit = 'R'
287  END IF
288 *
289  IF( imat.NE.ntypes ) THEN
290 *
291 * Set up parameters with CLATB4 and generate a test
292 * matrix with CLATMS.
293 *
294  CALL clatb4( path, imat, n, n, TYPE, kl, ku, anorm,
295  $ mode, cndnum, dist )
296 *
297  srnamt = 'CLATMS'
298  CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
299  $ cndnum, anorm, kl, ku, packit, a, lda,
300  $ work, info )
301 *
302 * Check error code from CLATMS.
303 *
304  IF( info.NE.0 ) THEN
305  CALL alaerh( path, 'CLATMS', info, 0, uplo, n, n,
306  $ -1, -1, -1, imat, nfail, nerrs, nout )
307  GO TO 160
308  END IF
309 *
310 * For types 3-6, zero one or more rows and columns of
311 * the matrix to test that INFO is returned correctly.
312 *
313  IF( zerot ) THEN
314  IF( imat.EQ.3 ) THEN
315  izero = 1
316  ELSE IF( imat.EQ.4 ) THEN
317  izero = n
318  ELSE
319  izero = n / 2 + 1
320  END IF
321 *
322  IF( imat.LT.6 ) THEN
323 *
324 * Set row and column IZERO to zero.
325 *
326  IF( iuplo.EQ.1 ) THEN
327  ioff = ( izero-1 )*izero / 2
328  DO 20 i = 1, izero - 1
329  a( ioff+i ) = zero
330  20 CONTINUE
331  ioff = ioff + izero
332  DO 30 i = izero, n
333  a( ioff ) = zero
334  ioff = ioff + i
335  30 CONTINUE
336  ELSE
337  ioff = izero
338  DO 40 i = 1, izero - 1
339  a( ioff ) = zero
340  ioff = ioff + n - i
341  40 CONTINUE
342  ioff = ioff - izero
343  DO 50 i = izero, n
344  a( ioff+i ) = zero
345  50 CONTINUE
346  END IF
347  ELSE
348  IF( iuplo.EQ.1 ) THEN
349 *
350 * Set the first IZERO rows and columns to zero.
351 *
352  ioff = 0
353  DO 70 j = 1, n
354  i2 = min( j, izero )
355  DO 60 i = 1, i2
356  a( ioff+i ) = zero
357  60 CONTINUE
358  ioff = ioff + j
359  70 CONTINUE
360  ELSE
361 *
362 * Set the last IZERO rows and columns to zero.
363 *
364  ioff = 0
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 + n - j
371  90 CONTINUE
372  END IF
373  END IF
374  ELSE
375  izero = 0
376  END IF
377  ELSE
378 *
379 * Use a special block diagonal matrix to test alternate
380 * code for the 2-by-2 blocks.
381 *
382  CALL clatsp( uplo, n, a, iseed )
383  END IF
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 * Compute the condition number for comparison with
392 * the value returned by CSPSVX.
393 *
394  IF( zerot ) THEN
395  IF( ifact.EQ.1 )
396  $ GO TO 150
397  rcondc = zero
398 *
399  ELSE IF( ifact.EQ.1 ) THEN
400 *
401 * Compute the 1-norm of A.
402 *
403  anorm = clansp( '1', uplo, n, a, rwork )
404 *
405 * Factor the matrix A.
406 *
407  CALL ccopy( npp, a, 1, afac, 1 )
408  CALL csptrf( uplo, n, afac, iwork, info )
409 *
410 * Compute inv(A) and take its norm.
411 *
412  CALL ccopy( npp, afac, 1, ainv, 1 )
413  CALL csptri( uplo, n, ainv, iwork, work, info )
414  ainvnm = clansp( '1', uplo, n, ainv, rwork )
415 *
416 * Compute the 1-norm condition number of A.
417 *
418  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
419  rcondc = one
420  ELSE
421  rcondc = ( one / anorm ) / ainvnm
422  END IF
423  END IF
424 *
425 * Form an exact solution and set the right hand side.
426 *
427  srnamt = 'CLARHS'
428  CALL clarhs( path, xtype, uplo, ' ', n, n, kl, ku,
429  $ nrhs, a, lda, xact, lda, b, lda, iseed,
430  $ info )
431  xtype = 'C'
432 *
433 * --- Test CSPSV ---
434 *
435  IF( ifact.EQ.2 ) THEN
436  CALL ccopy( npp, a, 1, afac, 1 )
437  CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
438 *
439 * Factor the matrix and solve the system using CSPSV.
440 *
441  srnamt = 'CSPSV '
442  CALL cspsv( uplo, n, nrhs, afac, iwork, x, lda,
443  $ info )
444 *
445 * Adjust the expected value of INFO to account for
446 * pivoting.
447 *
448  k = izero
449  IF( k.GT.0 ) THEN
450  100 CONTINUE
451  IF( iwork( k ).LT.0 ) THEN
452  IF( iwork( k ).NE.-k ) THEN
453  k = -iwork( k )
454  GO TO 100
455  END IF
456  ELSE IF( iwork( k ).NE.k ) THEN
457  k = iwork( k )
458  GO TO 100
459  END IF
460  END IF
461 *
462 * Check error code from CSPSV .
463 *
464  IF( info.NE.k ) THEN
465  CALL alaerh( path, 'CSPSV ', info, k, uplo, n,
466  $ n, -1, -1, nrhs, imat, nfail,
467  $ nerrs, nout )
468  GO TO 120
469  ELSE IF( info.NE.0 ) THEN
470  GO TO 120
471  END IF
472 *
473 * Reconstruct matrix from factors and compute
474 * residual.
475 *
476  CALL cspt01( uplo, n, a, afac, iwork, ainv, lda,
477  $ rwork, result( 1 ) )
478 *
479 * Compute residual of the computed solution.
480 *
481  CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
482  CALL cspt02( uplo, n, nrhs, a, x, lda, work, lda,
483  $ rwork, result( 2 ) )
484 *
485 * Check solution from generated exact solution.
486 *
487  CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
488  $ result( 3 ) )
489  nt = 3
490 *
491 * Print information about the tests that did not pass
492 * the threshold.
493 *
494  DO 110 k = 1, nt
495  IF( result( k ).GE.thresh ) THEN
496  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
497  $ CALL aladhd( nout, path )
498  WRITE( nout, fmt = 9999 )'CSPSV ', uplo, n,
499  $ imat, k, result( k )
500  nfail = nfail + 1
501  END IF
502  110 CONTINUE
503  nrun = nrun + nt
504  120 CONTINUE
505  END IF
506 *
507 * --- Test CSPSVX ---
508 *
509  IF( ifact.EQ.2 .AND. npp.GT.0 )
510  $ CALL claset( 'Full', npp, 1, cmplx( zero ),
511  $ cmplx( zero ), afac, npp )
512  CALL claset( 'Full', n, nrhs, cmplx( zero ),
513  $ cmplx( zero ), x, lda )
514 *
515 * Solve the system and compute the condition number and
516 * error bounds using CSPSVX.
517 *
518  srnamt = 'CSPSVX'
519  CALL cspsvx( fact, uplo, n, nrhs, a, afac, iwork, b,
520  $ lda, x, lda, rcond, rwork,
521  $ rwork( nrhs+1 ), work, rwork( 2*nrhs+1 ),
522  $ info )
523 *
524 * Adjust the expected value of INFO to account for
525 * pivoting.
526 *
527  k = izero
528  IF( k.GT.0 ) THEN
529  130 CONTINUE
530  IF( iwork( k ).LT.0 ) THEN
531  IF( iwork( k ).NE.-k ) THEN
532  k = -iwork( k )
533  GO TO 130
534  END IF
535  ELSE IF( iwork( k ).NE.k ) THEN
536  k = iwork( k )
537  GO TO 130
538  END IF
539  END IF
540 *
541 * Check the error code from CSPSVX.
542 *
543  IF( info.NE.k ) THEN
544  CALL alaerh( path, 'CSPSVX', info, k, fact // uplo,
545  $ n, n, -1, -1, nrhs, imat, nfail,
546  $ nerrs, nout )
547  GO TO 150
548  END IF
549 *
550  IF( info.EQ.0 ) THEN
551  IF( ifact.GE.2 ) THEN
552 *
553 * Reconstruct matrix from factors and compute
554 * residual.
555 *
556  CALL cspt01( uplo, n, a, afac, iwork, ainv, lda,
557  $ rwork( 2*nrhs+1 ), result( 1 ) )
558  k1 = 1
559  ELSE
560  k1 = 2
561  END IF
562 *
563 * Compute residual of the computed solution.
564 *
565  CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
566  CALL cspt02( uplo, n, nrhs, a, x, lda, work, lda,
567  $ rwork( 2*nrhs+1 ), result( 2 ) )
568 *
569 * Check solution from generated exact solution.
570 *
571  CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
572  $ result( 3 ) )
573 *
574 * Check the error bounds from iterative refinement.
575 *
576  CALL cppt05( uplo, n, nrhs, a, b, lda, x, lda,
577  $ xact, lda, rwork, rwork( nrhs+1 ),
578  $ result( 4 ) )
579  ELSE
580  k1 = 6
581  END IF
582 *
583 * Compare RCOND from CSPSVX with the computed value
584 * in RCONDC.
585 *
586  result( 6 ) = sget06( rcond, rcondc )
587 *
588 * Print information about the tests that did not pass
589 * the threshold.
590 *
591  DO 140 k = k1, 6
592  IF( result( k ).GE.thresh ) THEN
593  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
594  $ CALL aladhd( nout, path )
595  WRITE( nout, fmt = 9998 )'CSPSVX', fact, uplo,
596  $ n, imat, k, result( k )
597  nfail = nfail + 1
598  END IF
599  140 CONTINUE
600  nrun = nrun + 7 - k1
601 *
602  150 CONTINUE
603 *
604  160 CONTINUE
605  170 CONTINUE
606  180 CONTINUE
607 *
608 * Print a summary of the results.
609 *
610  CALL alasvm( path, nout, nfail, nrun, nerrs )
611 *
612  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
613  $ ', test ', i2, ', ratio =', g12.5 )
614  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N =', i5,
615  $ ', type ', i2, ', test ', i2, ', ratio =', g12.5 )
616  RETURN
617 *
618 * End of CDRVSP
619 *
real function clansp(NORM, UPLO, N, AP, WORK)
CLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric matrix supplied in packed form.
Definition: clansp.f:117
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine clarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
CLARHS
Definition: clarhs.f:211
subroutine cspsv(UPLO, N, NRHS, AP, IPIV, B, LDB, INFO)
CSPSV computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: cspsv.f:164
subroutine cppt05(UPLO, N, NRHS, AP, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
CPPT05
Definition: cppt05.f:159
subroutine cspt02(UPLO, N, NRHS, A, X, LDX, B, LDB, RWORK, RESID)
CSPT02
Definition: cspt02.f:125
subroutine cerrvx(PATH, NUNIT)
CERRVX
Definition: cerrvx.f:57
subroutine clatsp(UPLO, N, X, ISEED)
CLATSP
Definition: clatsp.f:86
subroutine cspt01(UPLO, N, A, AFAC, IPIV, C, LDC, RWORK, RESID)
CSPT01
Definition: cspt01.f:114
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine csptrf(UPLO, N, AP, IPIV, INFO)
CSPTRF
Definition: csptrf.f:160
subroutine cspsvx(FACT, UPLO, N, NRHS, AP, AFP, IPIV, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO)
CSPSVX computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: cspsvx.f:279
subroutine claset(UPLO, M, N, ALPHA, BETA, A, LDA)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: claset.f:108
subroutine clatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
CLATMS
Definition: clatms.f:334
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:52
subroutine csptri(UPLO, N, AP, IPIV, WORK, INFO)
CSPTRI
Definition: csptri.f:111
subroutine clatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
CLATB4
Definition: clatb4.f:123
subroutine cget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
CGET04
Definition: cget04.f:104

Here is the call graph for this function:

Here is the caller graph for this function: