LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zdrvpo ( 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( * )  ASAV,
complex*16, dimension( * )  B,
complex*16, dimension( * )  BSAV,
complex*16, dimension( * )  X,
complex*16, dimension( * )  XACT,
double precision, dimension( * )  S,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer  NOUT 
)

ZDRVPO

ZDRVPOX

Purpose:
 ZDRVPO tests the driver routines ZPOSV 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 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]ASAV
          ASAV is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV 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]S
          S is DOUBLE PRECISION array, dimension (NMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[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
Purpose:
 ZDRVPO tests the driver routines ZPOSV, -SVX, and -SVXX.

 Note that this file is used only when the XBLAS are available,
 otherwise zdrvpo.f defines this subroutine.
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]ASAV
          ASAV is COMPLEX*16 array, dimension (NMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV 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]S
          S is DOUBLE PRECISION array, dimension (NMAX)
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (NMAX+2*NRHS)
[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 2013

Definition at line 161 of file zdrvpo.f.

161 *
162 * -- LAPACK test routine (version 3.4.0) --
163 * -- LAPACK is a software package provided by Univ. of Tennessee, --
164 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
165 * November 2011
166 *
167 * .. Scalar Arguments ..
168  LOGICAL tsterr
169  INTEGER nmax, nn, nout, nrhs
170  DOUBLE PRECISION thresh
171 * ..
172 * .. Array Arguments ..
173  LOGICAL dotype( * )
174  INTEGER nval( * )
175  DOUBLE PRECISION rwork( * ), s( * )
176  COMPLEX*16 a( * ), afac( * ), asav( * ), b( * ),
177  $ bsav( * ), work( * ), x( * ), xact( * )
178 * ..
179 *
180 * =====================================================================
181 *
182 * .. Parameters ..
183  DOUBLE PRECISION one, zero
184  parameter ( one = 1.0d+0, zero = 0.0d+0 )
185  INTEGER ntypes
186  parameter ( ntypes = 9 )
187  INTEGER ntests
188  parameter ( ntests = 6 )
189 * ..
190 * .. Local Scalars ..
191  LOGICAL equil, nofact, prefac, zerot
192  CHARACTER dist, equed, fact, TYPE, uplo, xtype
193  CHARACTER*3 path
194  INTEGER i, iequed, ifact, imat, in, info, ioff, iuplo,
195  $ izero, k, k1, kl, ku, lda, mode, n, nb, nbmin,
196  $ nerrs, nfact, nfail, nimat, nrun, nt
197  DOUBLE PRECISION ainvnm, amax, anorm, cndnum, rcond, rcondc,
198  $ roldc, scond
199 * ..
200 * .. Local Arrays ..
201  CHARACTER equeds( 2 ), facts( 3 ), uplos( 2 )
202  INTEGER iseed( 4 ), iseedy( 4 )
203  DOUBLE PRECISION result( ntests )
204 * ..
205 * .. External Functions ..
206  LOGICAL lsame
207  DOUBLE PRECISION dget06, zlanhe
208  EXTERNAL lsame, dget06, zlanhe
209 * ..
210 * .. External Subroutines ..
211  EXTERNAL aladhd, alaerh, alasvm, xlaenv, zerrvx, zget04,
214  $ zpot05, zpotrf, zpotri
215 * ..
216 * .. Scalars in Common ..
217  LOGICAL lerr, ok
218  CHARACTER*32 srnamt
219  INTEGER infot, nunit
220 * ..
221 * .. Common blocks ..
222  COMMON / infoc / infot, nunit, ok, lerr
223  COMMON / srnamc / srnamt
224 * ..
225 * .. Intrinsic Functions ..
226  INTRINSIC dcmplx, max
227 * ..
228 * .. Data statements ..
229  DATA iseedy / 1988, 1989, 1990, 1991 /
230  DATA uplos / 'U', 'L' /
231  DATA facts / 'F', 'N', 'E' /
232  DATA equeds / 'N', 'Y' /
233 * ..
234 * .. Executable Statements ..
235 *
236 * Initialize constants and the random number seed.
237 *
238  path( 1: 1 ) = 'Zomplex precision'
239  path( 2: 3 ) = 'PO'
240  nrun = 0
241  nfail = 0
242  nerrs = 0
243  DO 10 i = 1, 4
244  iseed( i ) = iseedy( i )
245  10 CONTINUE
246 *
247 * Test the error exits
248 *
249  IF( tsterr )
250  $ CALL zerrvx( path, nout )
251  infot = 0
252 *
253 * Set the block size and minimum block size for testing.
254 *
255  nb = 1
256  nbmin = 2
257  CALL xlaenv( 1, nb )
258  CALL xlaenv( 2, nbmin )
259 *
260 * Do for each value of N in NVAL
261 *
262  DO 130 in = 1, nn
263  n = nval( in )
264  lda = max( n, 1 )
265  xtype = 'N'
266  nimat = ntypes
267  IF( n.LE.0 )
268  $ nimat = 1
269 *
270  DO 120 imat = 1, nimat
271 *
272 * Do the tests only if DOTYPE( IMAT ) is true.
273 *
274  IF( .NOT.dotype( imat ) )
275  $ GO TO 120
276 *
277 * Skip types 3, 4, or 5 if the matrix size is too small.
278 *
279  zerot = imat.GE.3 .AND. imat.LE.5
280  IF( zerot .AND. n.LT.imat-2 )
281  $ GO TO 120
282 *
283 * Do first for UPLO = 'U', then for UPLO = 'L'
284 *
285  DO 110 iuplo = 1, 2
286  uplo = uplos( iuplo )
287 *
288 * Set up parameters with ZLATB4 and generate a test matrix
289 * with ZLATMS.
290 *
291  CALL zlatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
292  $ cndnum, dist )
293 *
294  srnamt = 'ZLATMS'
295  CALL zlatms( n, n, dist, iseed, TYPE, rwork, mode,
296  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
297  $ info )
298 *
299 * Check error code from ZLATMS.
300 *
301  IF( info.NE.0 ) THEN
302  CALL alaerh( path, 'ZLATMS', info, 0, uplo, n, n, -1,
303  $ -1, -1, imat, nfail, nerrs, nout )
304  GO TO 110
305  END IF
306 *
307 * For types 3-5, zero one row and column of the matrix to
308 * test that INFO is returned correctly.
309 *
310  IF( zerot ) THEN
311  IF( imat.EQ.3 ) THEN
312  izero = 1
313  ELSE IF( imat.EQ.4 ) THEN
314  izero = n
315  ELSE
316  izero = n / 2 + 1
317  END IF
318  ioff = ( izero-1 )*lda
319 *
320 * Set row and column IZERO of A to 0.
321 *
322  IF( iuplo.EQ.1 ) THEN
323  DO 20 i = 1, izero - 1
324  a( ioff+i ) = zero
325  20 CONTINUE
326  ioff = ioff + izero
327  DO 30 i = izero, n
328  a( ioff ) = zero
329  ioff = ioff + lda
330  30 CONTINUE
331  ELSE
332  ioff = izero
333  DO 40 i = 1, izero - 1
334  a( ioff ) = zero
335  ioff = ioff + lda
336  40 CONTINUE
337  ioff = ioff - izero
338  DO 50 i = izero, n
339  a( ioff+i ) = zero
340  50 CONTINUE
341  END IF
342  ELSE
343  izero = 0
344  END IF
345 *
346 * Set the imaginary part of the diagonals.
347 *
348  CALL zlaipd( n, a, lda+1, 0 )
349 *
350 * Save a copy of the matrix A in ASAV.
351 *
352  CALL zlacpy( uplo, n, n, a, lda, asav, lda )
353 *
354  DO 100 iequed = 1, 2
355  equed = equeds( iequed )
356  IF( iequed.EQ.1 ) THEN
357  nfact = 3
358  ELSE
359  nfact = 1
360  END IF
361 *
362  DO 90 ifact = 1, nfact
363  fact = facts( ifact )
364  prefac = lsame( fact, 'F' )
365  nofact = lsame( fact, 'N' )
366  equil = lsame( fact, 'E' )
367 *
368  IF( zerot ) THEN
369  IF( prefac )
370  $ GO TO 90
371  rcondc = zero
372 *
373  ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
374 *
375 * Compute the condition number for comparison with
376 * the value returned by ZPOSVX (FACT = 'N' reuses
377 * the condition number from the previous iteration
378 * with FACT = 'F').
379 *
380  CALL zlacpy( uplo, n, n, asav, lda, afac, lda )
381  IF( equil .OR. iequed.GT.1 ) THEN
382 *
383 * Compute row and column scale factors to
384 * equilibrate the matrix A.
385 *
386  CALL zpoequ( n, afac, lda, s, scond, amax,
387  $ info )
388  IF( info.EQ.0 .AND. n.GT.0 ) THEN
389  IF( iequed.GT.1 )
390  $ scond = zero
391 *
392 * Equilibrate the matrix.
393 *
394  CALL zlaqhe( uplo, n, afac, lda, s, scond,
395  $ amax, equed )
396  END IF
397  END IF
398 *
399 * Save the condition number of the
400 * non-equilibrated system for use in ZGET04.
401 *
402  IF( equil )
403  $ roldc = rcondc
404 *
405 * Compute the 1-norm of A.
406 *
407  anorm = zlanhe( '1', uplo, n, afac, lda, rwork )
408 *
409 * Factor the matrix A.
410 *
411  CALL zpotrf( uplo, n, afac, lda, info )
412 *
413 * Form the inverse of A.
414 *
415  CALL zlacpy( uplo, n, n, afac, lda, a, lda )
416  CALL zpotri( uplo, n, a, lda, info )
417 *
418 * Compute the 1-norm condition number of A.
419 *
420  ainvnm = zlanhe( '1', uplo, n, a, lda, rwork )
421  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
422  rcondc = one
423  ELSE
424  rcondc = ( one / anorm ) / ainvnm
425  END IF
426  END IF
427 *
428 * Restore the matrix A.
429 *
430  CALL zlacpy( uplo, n, n, asav, lda, a, lda )
431 *
432 * Form an exact solution and set the right hand side.
433 *
434  srnamt = 'ZLARHS'
435  CALL zlarhs( path, xtype, uplo, ' ', n, n, kl, ku,
436  $ nrhs, a, lda, xact, lda, b, lda,
437  $ iseed, info )
438  xtype = 'C'
439  CALL zlacpy( 'Full', n, nrhs, b, lda, bsav, lda )
440 *
441  IF( nofact ) THEN
442 *
443 * --- Test ZPOSV ---
444 *
445 * Compute the L*L' or U'*U factorization of the
446 * matrix and solve the system.
447 *
448  CALL zlacpy( uplo, n, n, a, lda, afac, lda )
449  CALL zlacpy( 'Full', n, nrhs, b, lda, x, lda )
450 *
451  srnamt = 'ZPOSV '
452  CALL zposv( uplo, n, nrhs, afac, lda, x, lda,
453  $ info )
454 *
455 * Check error code from ZPOSV .
456 *
457  IF( info.NE.izero ) THEN
458  CALL alaerh( path, 'ZPOSV ', info, izero,
459  $ uplo, n, n, -1, -1, nrhs, imat,
460  $ nfail, nerrs, nout )
461  GO TO 70
462  ELSE IF( info.NE.0 ) THEN
463  GO TO 70
464  END IF
465 *
466 * Reconstruct matrix from factors and compute
467 * residual.
468 *
469  CALL zpot01( uplo, n, a, lda, afac, lda, rwork,
470  $ result( 1 ) )
471 *
472 * Compute residual of the computed solution.
473 *
474  CALL zlacpy( 'Full', n, nrhs, b, lda, work,
475  $ lda )
476  CALL zpot02( uplo, n, nrhs, a, lda, x, lda,
477  $ work, lda, rwork, result( 2 ) )
478 *
479 * Check solution from generated exact solution.
480 *
481  CALL zget04( n, nrhs, x, lda, xact, lda, rcondc,
482  $ result( 3 ) )
483  nt = 3
484 *
485 * Print information about the tests that did not
486 * pass the threshold.
487 *
488  DO 60 k = 1, nt
489  IF( result( k ).GE.thresh ) THEN
490  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
491  $ CALL aladhd( nout, path )
492  WRITE( nout, fmt = 9999 )'ZPOSV ', uplo,
493  $ n, imat, k, result( k )
494  nfail = nfail + 1
495  END IF
496  60 CONTINUE
497  nrun = nrun + nt
498  70 CONTINUE
499  END IF
500 *
501 * --- Test ZPOSVX ---
502 *
503  IF( .NOT.prefac )
504  $ CALL zlaset( uplo, n, n, dcmplx( zero ),
505  $ dcmplx( zero ), afac, lda )
506  CALL zlaset( 'Full', n, nrhs, dcmplx( zero ),
507  $ dcmplx( zero ), x, lda )
508  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
509 *
510 * Equilibrate the matrix if FACT='F' and
511 * EQUED='Y'.
512 *
513  CALL zlaqhe( uplo, n, a, lda, s, scond, amax,
514  $ equed )
515  END IF
516 *
517 * Solve the system and compute the condition number
518 * and error bounds using ZPOSVX.
519 *
520  srnamt = 'ZPOSVX'
521  CALL zposvx( fact, uplo, n, nrhs, a, lda, afac,
522  $ lda, equed, s, b, lda, x, lda, rcond,
523  $ rwork, rwork( nrhs+1 ), work,
524  $ rwork( 2*nrhs+1 ), info )
525 *
526 * Check the error code from ZPOSVX.
527 *
528  IF( info.NE.izero ) THEN
529  CALL alaerh( path, 'ZPOSVX', info, izero,
530  $ fact // uplo, n, n, -1, -1, nrhs,
531  $ imat, nfail, nerrs, nout )
532  GO TO 90
533  END IF
534 *
535  IF( info.EQ.0 ) THEN
536  IF( .NOT.prefac ) THEN
537 *
538 * Reconstruct matrix from factors and compute
539 * residual.
540 *
541  CALL zpot01( uplo, n, a, lda, afac, lda,
542  $ rwork( 2*nrhs+1 ), result( 1 ) )
543  k1 = 1
544  ELSE
545  k1 = 2
546  END IF
547 *
548 * Compute residual of the computed solution.
549 *
550  CALL zlacpy( 'Full', n, nrhs, bsav, lda, work,
551  $ lda )
552  CALL zpot02( uplo, n, nrhs, asav, lda, x, lda,
553  $ work, lda, rwork( 2*nrhs+1 ),
554  $ result( 2 ) )
555 *
556 * Check solution from generated exact solution.
557 *
558  IF( nofact .OR. ( prefac .AND. lsame( equed,
559  $ 'N' ) ) ) THEN
560  CALL zget04( n, nrhs, x, lda, xact, lda,
561  $ rcondc, result( 3 ) )
562  ELSE
563  CALL zget04( n, nrhs, x, lda, xact, lda,
564  $ roldc, result( 3 ) )
565  END IF
566 *
567 * Check the error bounds from iterative
568 * refinement.
569 *
570  CALL zpot05( uplo, n, nrhs, asav, lda, b, lda,
571  $ x, lda, xact, lda, rwork,
572  $ rwork( nrhs+1 ), result( 4 ) )
573  ELSE
574  k1 = 6
575  END IF
576 *
577 * Compare RCOND from ZPOSVX with the computed value
578 * in RCONDC.
579 *
580  result( 6 ) = dget06( rcond, rcondc )
581 *
582 * Print information about the tests that did not pass
583 * the threshold.
584 *
585  DO 80 k = k1, 6
586  IF( result( k ).GE.thresh ) THEN
587  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
588  $ CALL aladhd( nout, path )
589  IF( prefac ) THEN
590  WRITE( nout, fmt = 9997 )'ZPOSVX', fact,
591  $ uplo, n, equed, imat, k, result( k )
592  ELSE
593  WRITE( nout, fmt = 9998 )'ZPOSVX', fact,
594  $ uplo, n, imat, k, result( k )
595  END IF
596  nfail = nfail + 1
597  END IF
598  80 CONTINUE
599  nrun = nrun + 7 - k1
600  90 CONTINUE
601  100 CONTINUE
602  110 CONTINUE
603  120 CONTINUE
604  130 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 ', i1,
611  $ ', test(', i1, ')=', g12.5 )
612  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
613  $ ', type ', i1, ', test(', i1, ')=', g12.5 )
614  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
615  $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ') =',
616  $ g12.5 )
617  RETURN
618 *
619 * End of ZDRVPO
620 *
subroutine zpotrf(UPLO, N, A, LDA, INFO)
ZPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
Definition: zpotrf.f:102
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 zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
ZGET04
Definition: zget04.f:104
subroutine zlaqhe(UPLO, N, A, LDA, S, SCOND, AMAX, EQUED)
ZLAQHE scales a Hermitian matrix.
Definition: zlaqhe.f:136
subroutine zpot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZPOT02
Definition: zpot02.f:129
double precision function zlanhe(NORM, UPLO, N, A, LDA, WORK)
ZLANHE 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.
Definition: zlanhe.f:126
subroutine zlarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
ZLARHS
Definition: zlarhs.f:211
subroutine zposv(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
ZPOSV computes the solution to system of linear equations A * X = B for PO matrices ...
Definition: zposv.f:132
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:108
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine zlatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
ZLATB4
Definition: zlatb4.f:123
subroutine zlaipd(N, A, INDA, VINDA)
ZLAIPD
Definition: zlaipd.f:85
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine zpot01(UPLO, N, A, LDA, AFAC, LDAFAC, RWORK, RESID)
ZPOT01
Definition: zpot01.f:108
subroutine zpoequ(N, A, LDA, S, SCOND, AMAX, INFO)
ZPOEQU
Definition: zpoequ.f:115
double precision function dget06(RCOND, RCONDC)
DGET06
Definition: dget06.f:57
subroutine zerrvx(PATH, NUNIT)
ZERRVX
Definition: zerrvx.f:57
subroutine zlatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
ZLATMS
Definition: zlatms.f:334
subroutine zpotri(UPLO, N, A, LDA, INFO)
ZPOTRI
Definition: zpotri.f:97
subroutine zpot05(UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
ZPOT05
Definition: zpot05.f:167
subroutine zposvx(FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO)
ZPOSVX computes the solution to system of linear equations A * X = B for PO matrices ...
Definition: zposvx.f:308
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: