LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zla_porfsx_extended()

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 responsible 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' or '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 < 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 define 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.

Definition at line 380 of file zla_porfsx_extended.f.

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