LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages

◆ sporfsx()

subroutine sporfsx ( character uplo,
character equed,
integer n,
integer nrhs,
real, dimension( lda, * ) a,
integer lda,
real, dimension( ldaf, * ) af,
integer ldaf,
real, dimension( * ) s,
real, dimension( ldb, * ) b,
integer ldb,
real, dimension( ldx, * ) x,
integer ldx,
real rcond,
real, dimension( * ) berr,
integer n_err_bnds,
real, dimension( nrhs, * ) err_bnds_norm,
real, dimension( nrhs, * ) err_bnds_comp,
integer nparams,
real, dimension( * ) params,
real, dimension( * ) work,
integer, dimension( * ) iwork,
integer info )

SPORFSX

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

Purpose:
!> !> SPORFSX improves the computed solution to a system of linear !> equations when the coefficient matrix is symmetric positive !> definite, and provides error bounds and backward error estimates !> for the solution. 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. !> !> The original system of linear equations may have been equilibrated !> before calling this routine, as described by arguments EQUED and S !> below. In this case, the solution and error bounds returned are !> for the original unequilibrated system. !>
!> Some optional parameters are bundled in the PARAMS array. These !> settings determine how refinement is performed, but often the !> defaults are acceptable. If the defaults are acceptable, users !> can pass NPARAMS = 0 which prevents the source code from accessing !> the PARAMS argument. !>
Parameters
[in]UPLO
!> UPLO is CHARACTER*1 !> = 'U': Upper triangle of A is stored; !> = 'L': Lower triangle of A is stored. !>
[in]EQUED
!> EQUED is CHARACTER*1 !> Specifies the form of equilibration that was done to A !> before calling this routine. This is needed to compute !> the solution and error bounds correctly. !> = 'N': No equilibration !> = 'Y': Both row and column equilibration, i.e., A has been !> replaced by diag(S) * A * diag(S). !> The right hand side B has been changed accordingly. !>
[in]N
!> N is INTEGER !> 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 matrices B and X. NRHS >= 0. !>
[in]A
!> A is REAL array, dimension (LDA,N) !> The symmetric matrix A. If UPLO = 'U', the leading N-by-N !> upper triangular part of A contains the upper triangular part !> of the matrix A, and the strictly lower triangular part of A !> is not referenced. If UPLO = 'L', the leading N-by-N lower !> triangular part of A contains the lower triangular part of !> the matrix A, and the strictly upper triangular part of A is !> not referenced. !>
[in]LDA
!> LDA is INTEGER !> The leading dimension of the array A. LDA >= max(1,N). !>
[in]AF
!> AF is REAL 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 SPOTRF. !>
[in]LDAF
!> LDAF is INTEGER !> The leading dimension of the array AF. LDAF >= max(1,N). !>
[in,out]S
!> S is REAL array, dimension (N) !> The scale factors for A. If EQUED = 'Y', A is multiplied on !> the left and right by diag(S). S is an input argument if FACT = !> 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED !> = 'Y', each element of S must be positive. If S is output, each !> element of S is a power of the radix. If S is input, each element !> of S 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 REAL 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]X
!> X is REAL array, dimension (LDX,NRHS) !> On entry, the solution matrix X, as computed by SGETRS. !> On exit, the improved solution matrix X. !>
[in]LDX
!> LDX is INTEGER !> The leading dimension of the array X. LDX >= max(1,N). !>
[out]RCOND
!> RCOND is REAL !> 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. !>
[out]BERR
!> BERR is REAL array, dimension (NRHS) !> Componentwise relative backward error. This is the !> componentwise relative backward error of each solution vector X(j) !> (i.e., the smallest relative change in any element of A or B that !> makes X(j) an exact solution). !>
[in]N_ERR_BNDS
!> N_ERR_BNDS is INTEGER !> Number of error bounds to return for each right hand side !> and each type (normwise or componentwise). See ERR_BNDS_NORM and !> ERR_BNDS_COMP below. !>
[out]ERR_BNDS_NORM
!> ERR_BNDS_NORM is REAL 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. !> !> See Lapack Working Note 165 for further details and extra !> cautions. !>
[out]ERR_BNDS_COMP
!> ERR_BNDS_COMP is REAL 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. !> !> See Lapack Working Note 165 for further details and extra !> cautions. !>
[in]NPARAMS
!> NPARAMS is INTEGER !> Specifies the number of parameters set in PARAMS. If <= 0, the !> PARAMS array is never referenced and default values are used. !>
[in,out]PARAMS
!> PARAMS is REAL array, dimension NPARAMS !> Specifies algorithm parameters. If an entry is < 0.0, then !> that entry will be filled with default value used for that !> parameter. Only positions up to NPARAMS are accessed; defaults !> are used for higher-numbered parameters. !> !> PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative !> refinement or not. !> Default: 1.0 !> = 0.0: No refinement is performed, and no error bounds are !> computed. !> = 1.0: Use the double-precision refinement algorithm, !> possibly with doubled-single computations if the !> compilation environment does not support DOUBLE !> PRECISION. !> (other values are reserved for future use) !> !> PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual !> computations allowed for refinement. !> Default: 10 !> 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. !> !> PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code !> will attempt to find a solution with small componentwise !> relative error in the double-precision algorithm. Positive !> is true, 0.0 is false. !> Default: 1.0 (attempt componentwise convergence) !>
[out]WORK
!> WORK is REAL array, dimension (4*N) !>
[out]IWORK
!> IWORK is INTEGER array, dimension (N) !>
[out]INFO
!> INFO is INTEGER !> = 0: Successful exit. The solution to every right-hand side is !> guaranteed. !> < 0: If INFO = -i, the i-th argument had an illegal value !> > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization !> has been completed, but the factor U is exactly singular, so !> the solution and error bounds could not be computed. RCOND = 0 !> is returned. !> = N+J: The solution corresponding to the Jth right-hand side is !> not guaranteed. The solutions corresponding to other right- !> hand sides K with K > J may not be guaranteed as well, but !> only the first such right-hand side is reported. If a small !> componentwise error is not requested (PARAMS(3) = 0.0) then !> the Jth right-hand side is the first with a normwise error !> bound that is not guaranteed (the smallest J such !> that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0) !> the Jth right-hand side is the first with either a normwise or !> componentwise error bound that is not guaranteed (the smallest !> J such that either ERR_BNDS_NORM(J,1) = 0.0 or !> ERR_BNDS_COMP(J,1) = 0.0). See the definition of !> ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information !> about all of the right-hand sides check ERR_BNDS_NORM or !> ERR_BNDS_COMP. !>
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 388 of file sporfsx.f.

393*
394* -- LAPACK computational routine --
395* -- LAPACK is a software package provided by Univ. of Tennessee, --
396* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
397*
398* .. Scalar Arguments ..
399 CHARACTER UPLO, EQUED
400 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
401 $ N_ERR_BNDS
402 REAL RCOND
403* ..
404* .. Array Arguments ..
405 INTEGER IWORK( * )
406 REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
407 $ X( LDX, * ), WORK( * )
408 REAL S( * ), PARAMS( * ), BERR( * ),
409 $ ERR_BNDS_NORM( NRHS, * ),
410 $ ERR_BNDS_COMP( NRHS, * )
411* ..
412*
413* ==================================================================
414*
415* .. Parameters ..
416 REAL ZERO, ONE
417 parameter( zero = 0.0e+0, one = 1.0e+0 )
418 REAL ITREF_DEFAULT, ITHRESH_DEFAULT,
419 $ COMPONENTWISE_DEFAULT
420 REAL RTHRESH_DEFAULT, DZTHRESH_DEFAULT
421 parameter( itref_default = 1.0 )
422 parameter( ithresh_default = 10.0 )
423 parameter( componentwise_default = 1.0 )
424 parameter( rthresh_default = 0.5 )
425 parameter( dzthresh_default = 0.25 )
426 INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
427 $ LA_LINRX_CWISE_I
428 parameter( la_linrx_itref_i = 1,
429 $ la_linrx_ithresh_i = 2 )
430 parameter( la_linrx_cwise_i = 3 )
431 INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
432 $ LA_LINRX_RCOND_I
433 parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
434 parameter( la_linrx_rcond_i = 3 )
435* ..
436* .. Local Scalars ..
437 CHARACTER(1) NORM
438 LOGICAL RCEQU
439 INTEGER J, PREC_TYPE, REF_TYPE
440 INTEGER N_NORMS
441 REAL ANORM, RCOND_TMP
442 REAL ILLRCOND_THRESH, ERR_LBND, CWISE_WRONG
443 LOGICAL IGNORE_CWISE
444 INTEGER ITHRESH
445 REAL RTHRESH, UNSTABLE_THRESH
446* ..
447* .. External Subroutines ..
449* ..
450* .. Intrinsic Functions ..
451 INTRINSIC max, sqrt
452* ..
453* .. External Functions ..
454 EXTERNAL lsame, ilaprec
455 EXTERNAL slamch, slansy, sla_porcond
456 REAL SLAMCH, SLANSY, SLA_PORCOND
457 LOGICAL LSAME
458 INTEGER ILAPREC
459* ..
460* .. Executable Statements ..
461*
462* Check the input parameters.
463*
464 info = 0
465 ref_type = int( itref_default )
466 IF ( nparams .GE. la_linrx_itref_i ) THEN
467 IF ( params( la_linrx_itref_i ) .LT. 0.0 ) THEN
468 params( la_linrx_itref_i ) = itref_default
469 ELSE
470 ref_type = params( la_linrx_itref_i )
471 END IF
472 END IF
473*
474* Set default parameters.
475*
476 illrcond_thresh = real( n ) * slamch( 'Epsilon' )
477 ithresh = int( ithresh_default )
478 rthresh = rthresh_default
479 unstable_thresh = dzthresh_default
480 ignore_cwise = componentwise_default .EQ. 0.0
481*
482 IF ( nparams.GE.la_linrx_ithresh_i ) THEN
483 IF ( params( la_linrx_ithresh_i ).LT.0.0 ) THEN
484 params( la_linrx_ithresh_i ) = ithresh
485 ELSE
486 ithresh = int( params( la_linrx_ithresh_i ) )
487 END IF
488 END IF
489 IF ( nparams.GE.la_linrx_cwise_i ) THEN
490 IF ( params( la_linrx_cwise_i ).LT.0.0 ) THEN
491 IF ( ignore_cwise ) THEN
492 params( la_linrx_cwise_i ) = 0.0
493 ELSE
494 params( la_linrx_cwise_i ) = 1.0
495 END IF
496 ELSE
497 ignore_cwise = params( la_linrx_cwise_i ) .EQ. 0.0
498 END IF
499 END IF
500 IF ( ref_type .EQ. 0 .OR. n_err_bnds .EQ. 0 ) THEN
501 n_norms = 0
502 ELSE IF ( ignore_cwise ) THEN
503 n_norms = 1
504 ELSE
505 n_norms = 2
506 END IF
507*
508 rcequ = lsame( equed, 'Y' )
509*
510* Test input parameters.
511*
512 IF (.NOT.lsame(uplo, 'U') .AND. .NOT.lsame(uplo, 'L')) THEN
513 info = -1
514 ELSE IF( .NOT.rcequ .AND. .NOT.lsame( equed, 'N' ) ) THEN
515 info = -2
516 ELSE IF( n.LT.0 ) THEN
517 info = -3
518 ELSE IF( nrhs.LT.0 ) THEN
519 info = -4
520 ELSE IF( lda.LT.max( 1, n ) ) THEN
521 info = -6
522 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
523 info = -8
524 ELSE IF( ldb.LT.max( 1, n ) ) THEN
525 info = -11
526 ELSE IF( ldx.LT.max( 1, n ) ) THEN
527 info = -13
528 END IF
529 IF( info.NE.0 ) THEN
530 CALL xerbla( 'SPORFSX', -info )
531 RETURN
532 END IF
533*
534* Quick return if possible.
535*
536 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
537 rcond = 1.0
538 DO j = 1, nrhs
539 berr( j ) = 0.0
540 IF ( n_err_bnds .GE. 1 ) THEN
541 err_bnds_norm( j, la_linrx_trust_i ) = 1.0
542 err_bnds_comp( j, la_linrx_trust_i ) = 1.0
543 END IF
544 IF ( n_err_bnds .GE. 2 ) THEN
545 err_bnds_norm( j, la_linrx_err_i ) = 0.0
546 err_bnds_comp( j, la_linrx_err_i ) = 0.0
547 END IF
548 IF ( n_err_bnds .GE. 3 ) THEN
549 err_bnds_norm( j, la_linrx_rcond_i ) = 1.0
550 err_bnds_comp( j, la_linrx_rcond_i ) = 1.0
551 END IF
552 END DO
553 RETURN
554 END IF
555*
556* Default to failure.
557*
558 rcond = 0.0
559 DO j = 1, nrhs
560 berr( j ) = 1.0
561 IF ( n_err_bnds .GE. 1 ) THEN
562 err_bnds_norm( j, la_linrx_trust_i ) = 1.0
563 err_bnds_comp( j, la_linrx_trust_i ) = 1.0
564 END IF
565 IF ( n_err_bnds .GE. 2 ) THEN
566 err_bnds_norm( j, la_linrx_err_i ) = 1.0
567 err_bnds_comp( j, la_linrx_err_i ) = 1.0
568 END IF
569 IF ( n_err_bnds .GE. 3 ) THEN
570 err_bnds_norm( j, la_linrx_rcond_i ) = 0.0
571 err_bnds_comp( j, la_linrx_rcond_i ) = 0.0
572 END IF
573 END DO
574*
575* Compute the norm of A and the reciprocal of the condition
576* number of A.
577*
578 norm = 'I'
579 anorm = slansy( norm, uplo, n, a, lda, work )
580 CALL spocon( uplo, n, af, ldaf, anorm, rcond, work,
581 $ iwork, info )
582*
583* Perform refinement on each right-hand side
584*
585 IF ( ref_type .NE. 0 ) THEN
586
587 prec_type = ilaprec( 'D' )
588
589 CALL sla_porfsx_extended( prec_type, uplo, n,
590 $ nrhs, a, lda, af, ldaf, rcequ, s, b,
591 $ ldb, x, ldx, berr, n_norms, err_bnds_norm, err_bnds_comp,
592 $ work( n+1 ), work( 1 ), work( 2*n+1 ), work( 1 ), rcond,
593 $ ithresh, rthresh, unstable_thresh, ignore_cwise,
594 $ info )
595 END IF
596
597 err_lbnd = max( 10.0, sqrt( real( n ) ) ) * slamch( 'Epsilon' )
598 IF ( n_err_bnds .GE. 1 .AND. n_norms .GE. 1 ) THEN
599*
600* Compute scaled normwise condition number cond(A*C).
601*
602 IF ( rcequ ) THEN
603 rcond_tmp = sla_porcond( uplo, n, a, lda, af, ldaf,
604 $ -1, s, info, work, iwork )
605 ELSE
606 rcond_tmp = sla_porcond( uplo, n, a, lda, af, ldaf,
607 $ 0, s, info, work, iwork )
608 END IF
609 DO j = 1, nrhs
610*
611* Cap the error at 1.0.
612*
613 IF ( n_err_bnds .GE. la_linrx_err_i
614 $ .AND. err_bnds_norm( j, la_linrx_err_i ) .GT. 1.0 )
615 $ err_bnds_norm( j, la_linrx_err_i ) = 1.0
616*
617* Threshold the error (see LAWN).
618*
619 IF ( rcond_tmp .LT. illrcond_thresh ) THEN
620 err_bnds_norm( j, la_linrx_err_i ) = 1.0
621 err_bnds_norm( j, la_linrx_trust_i ) = 0.0
622 IF ( info .LE. n ) info = n + j
623 ELSE IF ( err_bnds_norm( j, la_linrx_err_i ) .LT. err_lbnd )
624 $ THEN
625 err_bnds_norm( j, la_linrx_err_i ) = err_lbnd
626 err_bnds_norm( j, la_linrx_trust_i ) = 1.0
627 END IF
628*
629* Save the condition number.
630*
631 IF (n_err_bnds .GE. la_linrx_rcond_i) THEN
632 err_bnds_norm( j, la_linrx_rcond_i ) = rcond_tmp
633 END IF
634 END DO
635 END IF
636
637 IF ( n_err_bnds .GE. 1 .AND. n_norms .GE. 2 ) THEN
638*
639* Compute componentwise condition number cond(A*diag(Y(:,J))) for
640* each right-hand side using the current solution as an estimate of
641* the true solution. If the componentwise error estimate is too
642* large, then the solution is a lousy estimate of truth and the
643* estimated RCOND may be too optimistic. To avoid misleading users,
644* the inverse condition number is set to 0.0 when the estimated
645* cwise error is at least CWISE_WRONG.
646*
647 cwise_wrong = sqrt( slamch( 'Epsilon' ) )
648 DO j = 1, nrhs
649 IF (err_bnds_comp( j, la_linrx_err_i ) .LT. cwise_wrong )
650 $ THEN
651 rcond_tmp = sla_porcond( uplo, n, a, lda, af, ldaf, 1,
652 $ x( 1, j ), info, work, iwork )
653 ELSE
654 rcond_tmp = 0.0
655 END IF
656*
657* Cap the error at 1.0.
658*
659 IF ( n_err_bnds .GE. la_linrx_err_i
660 $ .AND. err_bnds_comp( j, la_linrx_err_i ) .GT. 1.0 )
661 $ err_bnds_comp( j, la_linrx_err_i ) = 1.0
662*
663* Threshold the error (see LAWN).
664*
665 IF ( rcond_tmp .LT. illrcond_thresh ) THEN
666 err_bnds_comp( j, la_linrx_err_i ) = 1.0
667 err_bnds_comp( j, la_linrx_trust_i ) = 0.0
668 IF ( params( la_linrx_cwise_i ) .EQ. 1.0
669 $ .AND. info.LT.n + j ) info = n + j
670 ELSE IF ( err_bnds_comp( j, la_linrx_err_i )
671 $ .LT. err_lbnd ) THEN
672 err_bnds_comp( j, la_linrx_err_i ) = err_lbnd
673 err_bnds_comp( j, la_linrx_trust_i ) = 1.0
674 END IF
675*
676* Save the condition number.
677*
678 IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
679 err_bnds_comp( j, la_linrx_rcond_i ) = rcond_tmp
680 END IF
681
682 END DO
683 END IF
684*
685 RETURN
686*
687* End of SPORFSX
688*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
integer function ilaprec(prec)
ILAPREC
Definition ilaprec.f:56
real function sla_porcond(uplo, n, a, lda, af, ldaf, cmode, c, info, work, iwork)
SLA_PORCOND estimates the Skeel condition number for a symmetric positive-definite matrix.
subroutine sla_porfsx_extended(prec_type, uplo, n, nrhs, a, lda, af, ldaf, colequ, c, b, ldb, y, ldy, berr_out, n_norms, err_bnds_norm, err_bnds_comp, res, ayb, dy, y_tail, rcond, ithresh, rthresh, dz_ub, ignore_cwise, info)
SLA_PORFSX_EXTENDED improves the computed solution to a system of linear equations for symmetric or H...
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slansy(norm, uplo, n, a, lda, work)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansy.f:120
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine spocon(uplo, n, a, lda, anorm, rcond, work, iwork, info)
SPOCON
Definition spocon.f:119
Here is the call graph for this function:
Here is the caller graph for this function: