LAPACK 3.12.1
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  boolean. Trust the answer if the
!>              reciprocal condition number is less than the threshold
!>              sqrt(n) * slamch('Epsilon').
!>
!>     err = 2  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 . 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  boolean. Trust the answer if the
!>              reciprocal condition number is less than the threshold
!>              sqrt(n) * slamch('Epsilon').
!>
!>     err = 2  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 . 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 378 of file zla_porfsx_extended.f.

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