LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zla_porfsx_extended ( integer  PREC_TYPE,
character  UPLO,
integer  N,
integer  NRHS,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldaf, * )  AF,
integer  LDAF,
logical  COLEQU,
double precision, dimension( * )  C,
complex*16, dimension( ldb, * )  B,
integer  LDB,
complex*16, dimension( ldy, * )  Y,
integer  LDY,
double precision, dimension( * )  BERR_OUT,
integer  N_NORMS,
double precision, dimension( nrhs, * )  ERR_BNDS_NORM,
double precision, dimension( nrhs, * )  ERR_BNDS_COMP,
complex*16, dimension( * )  RES,
double precision, dimension( * )  AYB,
complex*16, dimension( * )  DY,
complex*16, dimension( * )  Y_TAIL,
double precision  RCOND,
integer  ITHRESH,
double precision  RTHRESH,
double precision  DZ_UB,
logical  IGNORE_CWISE,
integer  INFO 
)

ZLA_PORFSX_EXTENDED improves the computed solution to a system of linear equations for symmetric or Hermitian positive-definite matrices by performing extra-precise iterative refinement and provides error bounds and backward error estimates for the solution.

Download ZLA_PORFSX_EXTENDED + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZLA_PORFSX_EXTENDED improves the computed solution to a system of
 linear equations by performing extra-precise iterative refinement
 and provides error bounds and backward error estimates for the solution.
 This subroutine is called by ZPORFSX to perform iterative refinement.
 In addition to normwise error bound, the code provides maximum
 componentwise error bound if possible. See comments for ERR_BNDS_NORM
 and ERR_BNDS_COMP for details of the error bounds. Note that this
 subroutine is only resonsible for setting the second fields of
 ERR_BNDS_NORM and ERR_BNDS_COMP.
Parameters
[in]PREC_TYPE
          PREC_TYPE is INTEGER
     Specifies the intermediate precision to be used in refinement.
     The value is defined by ILAPREC(P) where P is a CHARACTER and
     P    = 'S':  Single
          = 'D':  Double
          = 'I':  Indigenous
          = 'X', 'E':  Extra
[in]UPLO
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
     The number of linear equations, i.e., the order of the
     matrix A.  N >= 0.
[in]NRHS
          NRHS is INTEGER
     The number of right-hand-sides, i.e., the number of columns of the
     matrix B.
[in]A
          A is COMPLEX*16 array, dimension (LDA,N)
     On entry, the N-by-N matrix A.
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in]AF
          AF is COMPLEX*16 array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by ZPOTRF.
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in]COLEQU
          COLEQU is LOGICAL
     If .TRUE. then column equilibration was done to A before calling
     this routine. This is needed to compute the solution and error
     bounds correctly.
[in]C
          C is DOUBLE PRECISION array, dimension (N)
     The column scale factors for A. If COLEQU = .FALSE., C
     is not accessed. If C is input, each element of C should be a power
     of the radix to ensure a reliable solution and error estimates.
     Scaling by powers of the radix does not cause rounding errors unless
     the result underflows or overflows. Rounding errors during scaling
     lead to refining with a matrix that is not equivalent to the
     input matrix, producing error estimates that may not be
     reliable.
[in]B
          B is COMPLEX*16 array, dimension (LDB,NRHS)
     The right-hand-side matrix B.
[in]LDB
          LDB is INTEGER
     The leading dimension of the array B.  LDB >= max(1,N).
[in,out]Y
          Y is COMPLEX*16 array, dimension
                    (LDY,NRHS)
     On entry, the solution matrix X, as computed by ZPOTRS.
     On exit, the improved solution matrix Y.
[in]LDY
          LDY is INTEGER
     The leading dimension of the array Y.  LDY >= max(1,N).
[out]BERR_OUT
          BERR_OUT is DOUBLE PRECISION array, dimension (NRHS)
     On exit, BERR_OUT(j) contains the componentwise relative backward
     error for right-hand-side j from the formula
         max(i) ( abs(RES(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
     where abs(Z) is the componentwise absolute value of the matrix
     or vector Z. This is computed by ZLA_LIN_BERR.
[in]N_NORMS
          N_NORMS is INTEGER
     Determines which error bounds to return (see ERR_BNDS_NORM
     and ERR_BNDS_COMP).
     If N_NORMS >= 1 return normwise error bounds.
     If N_NORMS >= 2 return componentwise error bounds.
[in,out]ERR_BNDS_NORM
          ERR_BNDS_NORM is DOUBLE PRECISION array, dimension
                    (NRHS, N_ERR_BNDS)
     For each right-hand side, this array contains information about
     various error bounds and condition numbers corresponding to the
     normwise relative error, which is defined as follows:

     Normwise relative error in the ith solution vector:
             max_j (abs(XTRUE(j,i) - X(j,i)))
            ------------------------------
                  max_j abs(X(j,i))

     The array is indexed by the type of error information as described
     below. There currently are up to three pieces of information
     returned.

     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
     right-hand side.

     The second index in ERR_BNDS_NORM(:,err) contains the following
     three fields:
     err = 1 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" error bound: The estimated forward error,
              almost certainly within a factor of 10 of the true error
              so long as the next entry is greater than the threshold
              sqrt(n) * slamch('Epsilon'). This error bound should only
              be trusted if the previous boolean is true.

     err = 3  Reciprocal condition number: Estimated normwise
              reciprocal condition number.  Compared with the threshold
              sqrt(n) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". These reciprocal condition
              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
              appropriately scaled matrix Z.
              Let Z = S*A, where S scales each row by a power of the
              radix so all absolute row sums of Z are approximately 1.

     This subroutine is only responsible for setting the second field
     above.
     See Lapack Working Note 165 for further details and extra
     cautions.
[in,out]ERR_BNDS_COMP
          ERR_BNDS_COMP is DOUBLE PRECISION array, dimension
                    (NRHS, N_ERR_BNDS)
     For each right-hand side, this array contains information about
     various error bounds and condition numbers corresponding to the
     componentwise relative error, which is defined as follows:

     Componentwise relative error in the ith solution vector:
                    abs(XTRUE(j,i) - X(j,i))
             max_j ----------------------
                         abs(X(j,i))

     The array is indexed by the right-hand side i (on which the
     componentwise relative error depends), and the type of error
     information as described below. There currently are up to three
     pieces of information returned for each right-hand side. If
     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS .LT. 3, then at most
     the first (:,N_ERR_BNDS) entries are returned.

     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
     right-hand side.

     The second index in ERR_BNDS_COMP(:,err) contains the following
     three fields:
     err = 1 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" error bound: The estimated forward error,
              almost certainly within a factor of 10 of the true error
              so long as the next entry is greater than the threshold
              sqrt(n) * slamch('Epsilon'). This error bound should only
              be trusted if the previous boolean is true.

     err = 3  Reciprocal condition number: Estimated componentwise
              reciprocal condition number.  Compared with the threshold
              sqrt(n) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". These reciprocal condition
              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
              appropriately scaled matrix Z.
              Let Z = S*(A*diag(x)), where x is the solution for the
              current right-hand side and S scales each row of
              A*diag(x) by a power of the radix so all absolute row
              sums of Z are approximately 1.

     This subroutine is only responsible for setting the second field
     above.
     See Lapack Working Note 165 for further details and extra
     cautions.
[in]RES
          RES is COMPLEX*16 array, dimension (N)
     Workspace to hold the intermediate residual.
[in]AYB
          AYB is DOUBLE PRECISION array, dimension (N)
     Workspace.
[in]DY
          DY is COMPLEX*16 PRECISION array, dimension (N)
     Workspace to hold the intermediate solution.
[in]Y_TAIL
          Y_TAIL is COMPLEX*16 array, dimension (N)
     Workspace to hold the trailing bits of the intermediate solution.
[in]RCOND
          RCOND is DOUBLE PRECISION
     Reciprocal scaled condition number.  This is an estimate of the
     reciprocal Skeel condition number of the matrix A after
     equilibration (if done).  If this is less than the machine
     precision (in particular, if it is zero), the matrix is singular
     to working precision.  Note that the error may still be small even
     if this number is very small and the matrix appears ill-
     conditioned.
[in]ITHRESH
          ITHRESH is INTEGER
     The maximum number of residual computations allowed for
     refinement. The default is 10. For 'aggressive' set to 100 to
     permit convergence using approximate factorizations or
     factorizations other than LU. If the factorization uses a
     technique other than Gaussian elimination, the guarantees in
     ERR_BNDS_NORM and ERR_BNDS_COMP may no longer be trustworthy.
[in]RTHRESH
          RTHRESH is DOUBLE PRECISION
     Determines when to stop refinement if the error estimate stops
     decreasing. Refinement will stop when the next solution no longer
     satisfies norm(dx_{i+1}) < RTHRESH * norm(dx_i) where norm(Z) is
     the infinity norm of Z. RTHRESH satisfies 0 < RTHRESH <= 1. The
     default value is 0.5. For 'aggressive' set to 0.9 to permit
     convergence on extremely ill-conditioned matrices. See LAWN 165
     for more details.
[in]DZ_UB
          DZ_UB is DOUBLE PRECISION
     Determines when to start considering componentwise convergence.
     Componentwise convergence is only considered after each component
     of the solution Y is stable, which we definte as the relative
     change in each component being less than DZ_UB. The default value
     is 0.25, requiring the first bit to be stable. See LAWN 165 for
     more details.
[in]IGNORE_CWISE
          IGNORE_CWISE is LOGICAL
     If .TRUE. then ignore componentwise convergence. Default value
     is .FALSE..
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
       < 0:  if INFO = -i, the ith argument to ZPOTRS had an illegal
             value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 392 of file zla_porfsx_extended.f.

392 *
393 * -- LAPACK computational routine (version 3.4.2) --
394 * -- LAPACK is a software package provided by Univ. of Tennessee, --
395 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
396 * September 2012
397 *
398 * .. Scalar Arguments ..
399  INTEGER info, lda, ldaf, ldb, ldy, n, nrhs, prec_type,
400  $ n_norms, ithresh
401  CHARACTER uplo
402  LOGICAL colequ, ignore_cwise
403  DOUBLE PRECISION rthresh, dz_ub
404 * ..
405 * .. Array Arguments ..
406  COMPLEX*16 a( lda, * ), af( ldaf, * ), b( ldb, * ),
407  $ y( ldy, * ), res( * ), dy( * ), y_tail( * )
408  DOUBLE PRECISION c( * ), ayb( * ), rcond, berr_out( * ),
409  $ err_bnds_norm( nrhs, * ),
410  $ err_bnds_comp( nrhs, * )
411 * ..
412 *
413 * =====================================================================
414 *
415 * .. Local Scalars ..
416  INTEGER uplo2, cnt, i, j, x_state, z_state,
417  $ y_prec_state
418  DOUBLE PRECISION yk, dyk, ymin, normy, normx, normdx, dxrat,
419  $ dzrat, prevnormdx, prev_dz_z, dxratmax,
420  $ dzratmax, dx_x, dz_z, final_dx_x, final_dz_z,
421  $ eps, hugeval, incr_thresh
422  LOGICAL incr_prec
423  COMPLEX*16 zdum
424 * ..
425 * .. Parameters ..
426  INTEGER unstable_state, working_state, conv_state,
427  $ noprog_state, base_residual, extra_residual,
428  $ extra_y
429  parameter ( unstable_state = 0, working_state = 1,
430  $ conv_state = 2, noprog_state = 3 )
431  parameter ( base_residual = 0, extra_residual = 1,
432  $ extra_y = 2 )
433  INTEGER final_nrm_err_i, final_cmp_err_i, berr_i
434  INTEGER rcond_i, nrm_rcond_i, nrm_err_i, cmp_rcond_i
435  INTEGER cmp_err_i, piv_growth_i
436  parameter ( final_nrm_err_i = 1, final_cmp_err_i = 2,
437  $ berr_i = 3 )
438  parameter ( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
439  parameter ( cmp_rcond_i = 7, cmp_err_i = 8,
440  $ piv_growth_i = 9 )
441  INTEGER la_linrx_itref_i, la_linrx_ithresh_i,
442  $ la_linrx_cwise_i
443  parameter ( la_linrx_itref_i = 1,
444  $ la_linrx_ithresh_i = 2 )
445  parameter ( la_linrx_cwise_i = 3 )
446  INTEGER la_linrx_trust_i, la_linrx_err_i,
447  $ la_linrx_rcond_i
448  parameter ( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
449  parameter ( la_linrx_rcond_i = 3 )
450 * ..
451 * .. External Functions ..
452  LOGICAL lsame
453  EXTERNAL ilauplo
454  INTEGER ilauplo
455 * ..
456 * .. External Subroutines ..
457  EXTERNAL zaxpy, zcopy, zpotrs, zhemv, blas_zhemv_x,
458  $ blas_zhemv2_x, zla_heamv, zla_wwaddw,
460  DOUBLE PRECISION dlamch
461 * ..
462 * .. Intrinsic Functions ..
463  INTRINSIC abs, dble, dimag, max, min
464 * ..
465 * .. Statement Functions ..
466  DOUBLE PRECISION cabs1
467 * ..
468 * .. Statement Function Definitions ..
469  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
470 * ..
471 * .. Executable Statements ..
472 *
473  IF (info.NE.0) RETURN
474  eps = dlamch( 'Epsilon' )
475  hugeval = dlamch( 'Overflow' )
476 * Force HUGEVAL to Inf
477  hugeval = hugeval * hugeval
478 * Using HUGEVAL may lead to spurious underflows.
479  incr_thresh = dble(n) * eps
480 
481  IF (lsame(uplo, 'L')) THEN
482  uplo2 = ilauplo( 'L' )
483  ELSE
484  uplo2 = ilauplo( 'U' )
485  ENDIF
486 
487  DO j = 1, nrhs
488  y_prec_state = extra_residual
489  IF (y_prec_state .EQ. extra_y) THEN
490  DO i = 1, n
491  y_tail( i ) = 0.0d+0
492  END DO
493  END IF
494 
495  dxrat = 0.0d+0
496  dxratmax = 0.0d+0
497  dzrat = 0.0d+0
498  dzratmax = 0.0d+0
499  final_dx_x = hugeval
500  final_dz_z = hugeval
501  prevnormdx = hugeval
502  prev_dz_z = hugeval
503  dz_z = hugeval
504  dx_x = hugeval
505 
506  x_state = working_state
507  z_state = unstable_state
508  incr_prec = .false.
509 
510  DO cnt = 1, ithresh
511 *
512 * Compute residual RES = B_s - op(A_s) * Y,
513 * op(A) = A, A**T, or A**H depending on TRANS (and type).
514 *
515  CALL zcopy( n, b( 1, j ), 1, res, 1 )
516  IF (y_prec_state .EQ. base_residual) THEN
517  CALL zhemv(uplo, n, dcmplx(-1.0d+0), a, lda, y(1,j), 1,
518  $ dcmplx(1.0d+0), res, 1)
519  ELSE IF (y_prec_state .EQ. extra_residual) THEN
520  CALL blas_zhemv_x(uplo2, n, dcmplx(-1.0d+0), a, lda,
521  $ y( 1, j ), 1, dcmplx(1.0d+0), res, 1, prec_type)
522  ELSE
523  CALL blas_zhemv2_x(uplo2, n, dcmplx(-1.0d+0), a, lda,
524  $ y(1, j), y_tail, 1, dcmplx(1.0d+0), res, 1,
525  $ prec_type)
526  END IF
527 
528 ! XXX: RES is no longer needed.
529  CALL zcopy( n, res, 1, dy, 1 )
530  CALL zpotrs( uplo, n, 1, af, ldaf, dy, n, info)
531 *
532 * Calculate relative changes DX_X, DZ_Z and ratios DXRAT, DZRAT.
533 *
534  normx = 0.0d+0
535  normy = 0.0d+0
536  normdx = 0.0d+0
537  dz_z = 0.0d+0
538  ymin = hugeval
539 
540  DO i = 1, n
541  yk = cabs1(y(i, j))
542  dyk = cabs1(dy(i))
543 
544  IF (yk .NE. 0.0d+0) THEN
545  dz_z = max( dz_z, dyk / yk )
546  ELSE IF (dyk .NE. 0.0d+0) THEN
547  dz_z = hugeval
548  END IF
549 
550  ymin = min( ymin, yk )
551 
552  normy = max( normy, yk )
553 
554  IF ( colequ ) THEN
555  normx = max(normx, yk * c(i))
556  normdx = max(normdx, dyk * c(i))
557  ELSE
558  normx = normy
559  normdx = max(normdx, dyk)
560  END IF
561  END DO
562 
563  IF (normx .NE. 0.0d+0) THEN
564  dx_x = normdx / normx
565  ELSE IF (normdx .EQ. 0.0d+0) THEN
566  dx_x = 0.0d+0
567  ELSE
568  dx_x = hugeval
569  END IF
570 
571  dxrat = normdx / prevnormdx
572  dzrat = dz_z / prev_dz_z
573 *
574 * Check termination criteria.
575 *
576  IF (ymin*rcond .LT. incr_thresh*normy
577  $ .AND. y_prec_state .LT. extra_y)
578  $ incr_prec = .true.
579 
580  IF (x_state .EQ. noprog_state .AND. dxrat .LE. rthresh)
581  $ x_state = working_state
582  IF (x_state .EQ. working_state) THEN
583  IF (dx_x .LE. eps) THEN
584  x_state = conv_state
585  ELSE IF (dxrat .GT. rthresh) THEN
586  IF (y_prec_state .NE. extra_y) THEN
587  incr_prec = .true.
588  ELSE
589  x_state = noprog_state
590  END IF
591  ELSE
592  IF (dxrat .GT. dxratmax) dxratmax = dxrat
593  END IF
594  IF (x_state .GT. working_state) final_dx_x = dx_x
595  END IF
596 
597  IF (z_state .EQ. unstable_state .AND. dz_z .LE. dz_ub)
598  $ z_state = working_state
599  IF (z_state .EQ. noprog_state .AND. dzrat .LE. rthresh)
600  $ z_state = working_state
601  IF (z_state .EQ. working_state) THEN
602  IF (dz_z .LE. eps) THEN
603  z_state = conv_state
604  ELSE IF (dz_z .GT. dz_ub) THEN
605  z_state = unstable_state
606  dzratmax = 0.0d+0
607  final_dz_z = hugeval
608  ELSE IF (dzrat .GT. rthresh) THEN
609  IF (y_prec_state .NE. extra_y) THEN
610  incr_prec = .true.
611  ELSE
612  z_state = noprog_state
613  END IF
614  ELSE
615  IF (dzrat .GT. dzratmax) dzratmax = dzrat
616  END IF
617  IF (z_state .GT. working_state) final_dz_z = dz_z
618  END IF
619 
620  IF ( x_state.NE.working_state.AND.
621  $ (ignore_cwise.OR.z_state.NE.working_state) )
622  $ GOTO 666
623 
624  IF (incr_prec) THEN
625  incr_prec = .false.
626  y_prec_state = y_prec_state + 1
627  DO i = 1, n
628  y_tail( i ) = 0.0d+0
629  END DO
630  END IF
631 
632  prevnormdx = normdx
633  prev_dz_z = dz_z
634 *
635 * Update soluton.
636 *
637  IF (y_prec_state .LT. extra_y) THEN
638  CALL zaxpy( n, dcmplx(1.0d+0), dy, 1, y(1,j), 1 )
639  ELSE
640  CALL zla_wwaddw(n, y(1,j), y_tail, dy)
641  END IF
642 
643  END DO
644 * Target of "IF (Z_STOP .AND. X_STOP)". Sun's f77 won't EXIT.
645  666 CONTINUE
646 *
647 * Set final_* when cnt hits ithresh.
648 *
649  IF (x_state .EQ. working_state) final_dx_x = dx_x
650  IF (z_state .EQ. working_state) final_dz_z = dz_z
651 *
652 * Compute error bounds.
653 *
654  IF (n_norms .GE. 1) THEN
655  err_bnds_norm( j, la_linrx_err_i ) =
656  $ final_dx_x / (1 - dxratmax)
657  END IF
658  IF (n_norms .GE. 2) THEN
659  err_bnds_comp( j, la_linrx_err_i ) =
660  $ final_dz_z / (1 - dzratmax)
661  END IF
662 *
663 * Compute componentwise relative backward error from formula
664 * max(i) ( abs(R(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
665 * where abs(Z) is the componentwise absolute value of the matrix
666 * or vector Z.
667 *
668 * Compute residual RES = B_s - op(A_s) * Y,
669 * op(A) = A, A**T, or A**H depending on TRANS (and type).
670 *
671  CALL zcopy( n, b( 1, j ), 1, res, 1 )
672  CALL zhemv(uplo, n, dcmplx(-1.0d+0), a, lda, y(1,j), 1,
673  $ dcmplx(1.0d+0), res, 1)
674 
675  DO i = 1, n
676  ayb( i ) = cabs1( b( i, j ) )
677  END DO
678 *
679 * Compute abs(op(A_s))*abs(Y) + abs(B_s).
680 *
681  CALL zla_heamv (uplo2, n, 1.0d+0,
682  $ a, lda, y(1, j), 1, 1.0d+0, ayb, 1)
683 
684  CALL zla_lin_berr (n, n, 1, res, ayb, berr_out(j))
685 *
686 * End of loop for each RHS.
687 *
688  END DO
689 *
690  RETURN
subroutine zhemv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZHEMV
Definition: zhemv.f:156
subroutine zla_heamv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZLA_HEAMV computes a matrix-vector product using a Hermitian indefinite matrix to calculate error bou...
Definition: zla_heamv.f:180
subroutine zla_lin_berr(N, NZ, NRHS, RES, AYB, BERR)
ZLA_LIN_BERR computes a component-wise relative backward error.
Definition: zla_lin_berr.f:103
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zla_wwaddw(N, X, Y, W)
ZLA_WWADDW adds a vector into a doubled-single vector.
Definition: zla_wwaddw.f:83
integer function ilauplo(UPLO)
ILAUPLO
Definition: ilauplo.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:53
subroutine zpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
ZPOTRS
Definition: zpotrs.f:112

Here is the call graph for this function:

Here is the caller graph for this function: