LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cdrvhp ( 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 
)

CDRVHP

Purpose:
 CDRVHP tests the driver routines CHPSV 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 cdrvhp.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 = 10, 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 clanhp, sget06
204  EXTERNAL clanhp, sget06
205 * ..
206 * .. External Subroutines ..
207  EXTERNAL aladhd, alaerh, alasvm, ccopy, cerrvx, cget04,
210  $ cppt05, 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 ) = 'C'
233  path( 2: 3 ) = 'HP'
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 * Set up parameters with CLATB4 and generate a test matrix
290 * with CLATMS.
291 *
292  CALL clatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
293  $ cndnum, dist )
294 *
295  srnamt = 'CLATMS'
296  CALL clatms( n, n, dist, iseed, TYPE, rwork, mode,
297  $ cndnum, anorm, kl, ku, packit, a, lda, work,
298  $ info )
299 *
300 * Check error code from CLATMS.
301 *
302  IF( info.NE.0 ) THEN
303  CALL alaerh( path, 'CLATMS', info, 0, uplo, n, n, -1,
304  $ -1, -1, imat, nfail, nerrs, nout )
305  GO TO 160
306  END IF
307 *
308 * For types 3-6, zero one or more rows and columns of the
309 * matrix to test that INFO is returned correctly.
310 *
311  IF( zerot ) THEN
312  IF( imat.EQ.3 ) THEN
313  izero = 1
314  ELSE IF( imat.EQ.4 ) THEN
315  izero = n
316  ELSE
317  izero = n / 2 + 1
318  END IF
319 *
320  IF( imat.LT.6 ) THEN
321 *
322 * Set row and column IZERO to zero.
323 *
324  IF( iuplo.EQ.1 ) THEN
325  ioff = ( izero-1 )*izero / 2
326  DO 20 i = 1, izero - 1
327  a( ioff+i ) = zero
328  20 CONTINUE
329  ioff = ioff + izero
330  DO 30 i = izero, n
331  a( ioff ) = zero
332  ioff = ioff + i
333  30 CONTINUE
334  ELSE
335  ioff = izero
336  DO 40 i = 1, izero - 1
337  a( ioff ) = zero
338  ioff = ioff + n - i
339  40 CONTINUE
340  ioff = ioff - izero
341  DO 50 i = izero, n
342  a( ioff+i ) = zero
343  50 CONTINUE
344  END IF
345  ELSE
346  ioff = 0
347  IF( iuplo.EQ.1 ) THEN
348 *
349 * Set the first IZERO rows and columns to zero.
350 *
351  DO 70 j = 1, n
352  i2 = min( j, izero )
353  DO 60 i = 1, i2
354  a( ioff+i ) = zero
355  60 CONTINUE
356  ioff = ioff + j
357  70 CONTINUE
358  ELSE
359 *
360 * Set the last IZERO rows and columns to zero.
361 *
362  DO 90 j = 1, n
363  i1 = max( j, izero )
364  DO 80 i = i1, n
365  a( ioff+i ) = zero
366  80 CONTINUE
367  ioff = ioff + n - j
368  90 CONTINUE
369  END IF
370  END IF
371  ELSE
372  izero = 0
373  END IF
374 *
375 * Set the imaginary part of the diagonals.
376 *
377  IF( iuplo.EQ.1 ) THEN
378  CALL claipd( n, a, 2, 1 )
379  ELSE
380  CALL claipd( n, a, n, -1 )
381  END IF
382 *
383  DO 150 ifact = 1, nfact
384 *
385 * Do first for FACT = 'F', then for other values.
386 *
387  fact = facts( ifact )
388 *
389 * Compute the condition number for comparison with
390 * the value returned by CHPSVX.
391 *
392  IF( zerot ) THEN
393  IF( ifact.EQ.1 )
394  $ GO TO 150
395  rcondc = zero
396 *
397  ELSE IF( ifact.EQ.1 ) THEN
398 *
399 * Compute the 1-norm of A.
400 *
401  anorm = clanhp( '1', uplo, n, a, rwork )
402 *
403 * Factor the matrix A.
404 *
405  CALL ccopy( npp, a, 1, afac, 1 )
406  CALL chptrf( uplo, n, afac, iwork, info )
407 *
408 * Compute inv(A) and take its norm.
409 *
410  CALL ccopy( npp, afac, 1, ainv, 1 )
411  CALL chptri( uplo, n, ainv, iwork, work, info )
412  ainvnm = clanhp( '1', uplo, n, ainv, rwork )
413 *
414 * Compute the 1-norm condition number of A.
415 *
416  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
417  rcondc = one
418  ELSE
419  rcondc = ( one / anorm ) / ainvnm
420  END IF
421  END IF
422 *
423 * Form an exact solution and set the right hand side.
424 *
425  srnamt = 'CLARHS'
426  CALL clarhs( path, xtype, uplo, ' ', n, n, kl, ku,
427  $ nrhs, a, lda, xact, lda, b, lda, iseed,
428  $ info )
429  xtype = 'C'
430 *
431 * --- Test CHPSV ---
432 *
433  IF( ifact.EQ.2 ) THEN
434  CALL ccopy( npp, a, 1, afac, 1 )
435  CALL clacpy( 'Full', n, nrhs, b, lda, x, lda )
436 *
437 * Factor the matrix and solve the system using CHPSV.
438 *
439  srnamt = 'CHPSV '
440  CALL chpsv( uplo, n, nrhs, afac, iwork, x, lda,
441  $ info )
442 *
443 * Adjust the expected value of INFO to account for
444 * pivoting.
445 *
446  k = izero
447  IF( k.GT.0 ) THEN
448  100 CONTINUE
449  IF( iwork( k ).LT.0 ) THEN
450  IF( iwork( k ).NE.-k ) THEN
451  k = -iwork( k )
452  GO TO 100
453  END IF
454  ELSE IF( iwork( k ).NE.k ) THEN
455  k = iwork( k )
456  GO TO 100
457  END IF
458  END IF
459 *
460 * Check error code from CHPSV .
461 *
462  IF( info.NE.k ) THEN
463  CALL alaerh( path, 'CHPSV ', info, k, uplo, n,
464  $ n, -1, -1, nrhs, imat, nfail,
465  $ nerrs, nout )
466  GO TO 120
467  ELSE IF( info.NE.0 ) THEN
468  GO TO 120
469  END IF
470 *
471 * Reconstruct matrix from factors and compute
472 * residual.
473 *
474  CALL chpt01( uplo, n, a, afac, iwork, ainv, lda,
475  $ rwork, result( 1 ) )
476 *
477 * Compute residual of the computed solution.
478 *
479  CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
480  CALL cppt02( uplo, n, nrhs, a, x, lda, work, lda,
481  $ rwork, result( 2 ) )
482 *
483 * Check solution from generated exact solution.
484 *
485  CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
486  $ result( 3 ) )
487  nt = 3
488 *
489 * Print information about the tests that did not pass
490 * the threshold.
491 *
492  DO 110 k = 1, nt
493  IF( result( k ).GE.thresh ) THEN
494  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
495  $ CALL aladhd( nout, path )
496  WRITE( nout, fmt = 9999 )'CHPSV ', uplo, n,
497  $ imat, k, result( k )
498  nfail = nfail + 1
499  END IF
500  110 CONTINUE
501  nrun = nrun + nt
502  120 CONTINUE
503  END IF
504 *
505 * --- Test CHPSVX ---
506 *
507  IF( ifact.EQ.2 .AND. npp.GT.0 )
508  $ CALL claset( 'Full', npp, 1, cmplx( zero ),
509  $ cmplx( zero ), afac, npp )
510  CALL claset( 'Full', n, nrhs, cmplx( zero ),
511  $ cmplx( zero ), x, lda )
512 *
513 * Solve the system and compute the condition number and
514 * error bounds using CHPSVX.
515 *
516  srnamt = 'CHPSVX'
517  CALL chpsvx( fact, uplo, n, nrhs, a, afac, iwork, b,
518  $ lda, x, lda, rcond, rwork,
519  $ rwork( nrhs+1 ), work, rwork( 2*nrhs+1 ),
520  $ info )
521 *
522 * Adjust the expected value of INFO to account for
523 * pivoting.
524 *
525  k = izero
526  IF( k.GT.0 ) THEN
527  130 CONTINUE
528  IF( iwork( k ).LT.0 ) THEN
529  IF( iwork( k ).NE.-k ) THEN
530  k = -iwork( k )
531  GO TO 130
532  END IF
533  ELSE IF( iwork( k ).NE.k ) THEN
534  k = iwork( k )
535  GO TO 130
536  END IF
537  END IF
538 *
539 * Check the error code from CHPSVX.
540 *
541  IF( info.NE.k ) THEN
542  CALL alaerh( path, 'CHPSVX', info, k, fact // uplo,
543  $ n, n, -1, -1, nrhs, imat, nfail,
544  $ nerrs, nout )
545  GO TO 150
546  END IF
547 *
548  IF( info.EQ.0 ) THEN
549  IF( ifact.GE.2 ) THEN
550 *
551 * Reconstruct matrix from factors and compute
552 * residual.
553 *
554  CALL chpt01( uplo, n, a, afac, iwork, ainv, lda,
555  $ rwork( 2*nrhs+1 ), result( 1 ) )
556  k1 = 1
557  ELSE
558  k1 = 2
559  END IF
560 *
561 * Compute residual of the computed solution.
562 *
563  CALL clacpy( 'Full', n, nrhs, b, lda, work, lda )
564  CALL cppt02( uplo, n, nrhs, a, x, lda, work, lda,
565  $ rwork( 2*nrhs+1 ), result( 2 ) )
566 *
567 * Check solution from generated exact solution.
568 *
569  CALL cget04( n, nrhs, x, lda, xact, lda, rcondc,
570  $ result( 3 ) )
571 *
572 * Check the error bounds from iterative refinement.
573 *
574  CALL cppt05( uplo, n, nrhs, a, b, lda, x, lda,
575  $ xact, lda, rwork, rwork( nrhs+1 ),
576  $ result( 4 ) )
577  ELSE
578  k1 = 6
579  END IF
580 *
581 * Compare RCOND from CHPSVX with the computed value
582 * in RCONDC.
583 *
584  result( 6 ) = sget06( rcond, rcondc )
585 *
586 * Print information about the tests that did not pass
587 * the threshold.
588 *
589  DO 140 k = k1, 6
590  IF( result( k ).GE.thresh ) THEN
591  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
592  $ CALL aladhd( nout, path )
593  WRITE( nout, fmt = 9998 )'CHPSVX', fact, uplo,
594  $ n, imat, k, result( k )
595  nfail = nfail + 1
596  END IF
597  140 CONTINUE
598  nrun = nrun + 7 - k1
599 *
600  150 CONTINUE
601 *
602  160 CONTINUE
603  170 CONTINUE
604  180 CONTINUE
605 *
606 * Print a summary of the results.
607 *
608  CALL alasvm( path, nout, nfail, nrun, nerrs )
609 *
610  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i2,
611  $ ', test ', i2, ', ratio =', g12.5 )
612  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N =', i5,
613  $ ', type ', i2, ', test ', i2, ', ratio =', g12.5 )
614  RETURN
615 *
616 * End of CDRVHP
617 *
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 chpsv(UPLO, N, NRHS, AP, IPIV, B, LDB, INFO)
CHPSV computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: chpsv.f:164
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 claipd(N, A, INDA, VINDA)
CLAIPD
Definition: claipd.f:85
subroutine cppt05(UPLO, N, NRHS, AP, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
CPPT05
Definition: cppt05.f:159
subroutine cerrvx(PATH, NUNIT)
CERRVX
Definition: cerrvx.f:57
subroutine chptri(UPLO, N, AP, IPIV, WORK, INFO)
CHPTRI
Definition: chptri.f:111
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine chptrf(UPLO, N, AP, IPIV, INFO)
CHPTRF
Definition: chptrf.f:161
subroutine cppt02(UPLO, N, NRHS, A, X, LDX, B, LDB, RWORK, RESID)
CPPT02
Definition: cppt02.f:125
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
real function clanhp(NORM, UPLO, N, AP, WORK)
CLANHP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix supplied in packed form.
Definition: clanhp.f:119
subroutine clatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
CLATB4
Definition: clatb4.f:123
subroutine chpt01(UPLO, N, A, AFAC, IPIV, C, LDC, RWORK, RESID)
CHPT01
Definition: chpt01.f:115
subroutine cget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
CGET04
Definition: cget04.f:104
subroutine chpsvx(FACT, UPLO, N, NRHS, AP, AFP, IPIV, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO)
CHPSVX computes the solution to system of linear equations A * X = B for OTHER matrices ...
Definition: chpsvx.f:279

Here is the call graph for this function:

Here is the caller graph for this function: