OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches

Functions

subroutine cposv (uplo, n, nrhs, a, lda, b, ldb, info)
  CPOSV computes the solution to system of linear equations A * X = B for PO matrices
subroutine cposvx (fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
  CPOSVX computes the solution to system of linear equations A * X = B for PO matrices
subroutine cposvxx (fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, rwork, info)
  CPOSVXX computes the solution to system of linear equations A * X = B for PO matrices

Detailed Description

This is the group of complex solve driver functions for PO matrices

Function Documentation

◆ cposv()

subroutine cposv ( character uplo,
integer n,
integer nrhs,
complex, dimension( lda, * ) a,
integer lda,
complex, dimension( ldb, * ) b,
integer ldb,
integer info )

CPOSV computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
!>
!> CPOSV computes the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian positive definite matrix and X and B
!> are N-by-NRHS matrices.
!>
!> The Cholesky decomposition is used to factor A as
!>    A = U**H* U,  if UPLO = 'U', or
!>    A = L * L**H,  if UPLO = 'L',
!> where U is an upper triangular matrix and  L is a lower triangular
!> matrix.  The factored form of A is then used to solve the system of
!> equations A * X = B.
!> 
Parameters
[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.  NRHS >= 0.
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA,N)
!>          On entry, the Hermitian 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.
!>
!>          On exit, if INFO = 0, the factor U or L from the Cholesky
!>          factorization A = U**H*U or A = L*L**H.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in,out]B
!>          B is COMPLEX array, dimension (LDB,NRHS)
!>          On entry, the N-by-NRHS right hand side matrix B.
!>          On exit, if INFO = 0, the N-by-NRHS solution matrix X.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!>          > 0:  if INFO = i, the leading minor of order i of A is not
!>                positive definite, so the factorization could not be
!>                completed, and the solution has not been computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 129 of file cposv.f.

130*
131* -- LAPACK driver routine --
132* -- LAPACK is a software package provided by Univ. of Tennessee, --
133* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
134*
135* .. Scalar Arguments ..
136 CHARACTER UPLO
137 INTEGER INFO, LDA, LDB, N, NRHS
138* ..
139* .. Array Arguments ..
140 COMPLEX A( LDA, * ), B( LDB, * )
141* ..
142*
143* =====================================================================
144*
145* .. External Functions ..
146 LOGICAL LSAME
147 EXTERNAL lsame
148* ..
149* .. External Subroutines ..
150 EXTERNAL cpotrf, cpotrs, xerbla
151* ..
152* .. Intrinsic Functions ..
153 INTRINSIC max
154* ..
155* .. Executable Statements ..
156*
157* Test the input parameters.
158*
159 info = 0
160 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
161 info = -1
162 ELSE IF( n.LT.0 ) THEN
163 info = -2
164 ELSE IF( nrhs.LT.0 ) THEN
165 info = -3
166 ELSE IF( lda.LT.max( 1, n ) ) THEN
167 info = -5
168 ELSE IF( ldb.LT.max( 1, n ) ) THEN
169 info = -7
170 END IF
171 IF( info.NE.0 ) THEN
172 CALL xerbla( 'CPOSV ', -info )
173 RETURN
174 END IF
175*
176* Compute the Cholesky factorization A = U**H*U or A = L*L**H.
177*
178 CALL cpotrf( uplo, n, a, lda, info )
179 IF( info.EQ.0 ) THEN
180*
181* Solve the system A*X = B, overwriting B with X.
182*
183 CALL cpotrs( uplo, n, nrhs, a, lda, b, ldb, info )
184*
185 END IF
186 RETURN
187*
188* End of CPOSV
189*
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
subroutine cpotrs(uplo, n, nrhs, a, lda, b, ldb, info)
CPOTRS
Definition cpotrs.f:110
subroutine cpotrf(uplo, n, a, lda, info)
CPOTRF
Definition cpotrf.f:107
#define max(a, b)
Definition macros.h:21

◆ cposvx()

subroutine cposvx ( character fact,
character uplo,
integer n,
integer nrhs,
complex, dimension( lda, * ) a,
integer lda,
complex, dimension( ldaf, * ) af,
integer ldaf,
character equed,
real, dimension( * ) s,
complex, dimension( ldb, * ) b,
integer ldb,
complex, dimension( ldx, * ) x,
integer ldx,
real rcond,
real, dimension( * ) ferr,
real, dimension( * ) berr,
complex, dimension( * ) work,
real, dimension( * ) rwork,
integer info )

CPOSVX computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
!>
!> CPOSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to
!> compute the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian positive definite matrix and X and B
!> are N-by-NRHS matrices.
!>
!> Error bounds on the solution and a condition estimate are also
!> provided.
!> 
Description:
!>
!> The following steps are performed:
!>
!> 1. If FACT = 'E', real scaling factors are computed to equilibrate
!>    the system:
!>       diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B
!>    Whether or not the system will be equilibrated depends on the
!>    scaling of the matrix A, but if equilibration is used, A is
!>    overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
!>
!> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
!>    factor the matrix A (after equilibration if FACT = 'E') as
!>       A = U**H* U,  if UPLO = 'U', or
!>       A = L * L**H,  if UPLO = 'L',
!>    where U is an upper triangular matrix and L is a lower triangular
!>    matrix.
!>
!> 3. If the leading i-by-i principal minor is not positive definite,
!>    then the routine returns with INFO = i. Otherwise, the factored
!>    form of A is used to estimate the condition number of the matrix
!>    A.  If the reciprocal of the condition number is less than machine
!>    precision, INFO = N+1 is returned as a warning, but the routine
!>    still goes on to solve for X and compute error bounds as
!>    described below.
!>
!> 4. The system of equations is solved for X using the factored form
!>    of A.
!>
!> 5. Iterative refinement is applied to improve the computed solution
!>    matrix and calculate error bounds and backward error estimates
!>    for it.
!>
!> 6. If equilibration was used, the matrix X is premultiplied by
!>    diag(S) so that it solves the original system before
!>    equilibration.
!> 
Parameters
[in]FACT
!>          FACT is CHARACTER*1
!>          Specifies whether or not the factored form of the matrix A is
!>          supplied on entry, and if not, whether the matrix A should be
!>          equilibrated before it is factored.
!>          = 'F':  On entry, AF contains the factored form of A.
!>                  If EQUED = 'Y', the matrix A has been equilibrated
!>                  with scaling factors given by S.  A and AF will not
!>                  be modified.
!>          = 'N':  The matrix A will be copied to AF and factored.
!>          = 'E':  The matrix A will be equilibrated if necessary, then
!>                  copied to AF and factored.
!> 
[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 matrices B and X.  NRHS >= 0.
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA,N)
!>          On entry, the Hermitian matrix A, except if FACT = 'F' and
!>          EQUED = 'Y', then A must contain the equilibrated matrix
!>          diag(S)*A*diag(S).  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.  A is not modified if
!>          FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit.
!>
!>          On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
!>          diag(S)*A*diag(S).
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in,out]AF
!>          AF is COMPLEX array, dimension (LDAF,N)
!>          If FACT = 'F', then AF is an input argument and on entry
!>          contains the triangular factor U or L from the Cholesky
!>          factorization A = U**H*U or A = L*L**H, in the same storage
!>          format as A.  If EQUED .ne. 'N', then AF is the factored form
!>          of the equilibrated matrix diag(S)*A*diag(S).
!>
!>          If FACT = 'N', then AF is an output argument and on exit
!>          returns the triangular factor U or L from the Cholesky
!>          factorization A = U**H*U or A = L*L**H of the original
!>          matrix A.
!>
!>          If FACT = 'E', then AF is an output argument and on exit
!>          returns the triangular factor U or L from the Cholesky
!>          factorization A = U**H*U or A = L*L**H of the equilibrated
!>          matrix A (see the description of A for the form of the
!>          equilibrated matrix).
!> 
[in]LDAF
!>          LDAF is INTEGER
!>          The leading dimension of the array AF.  LDAF >= max(1,N).
!> 
[in,out]EQUED
!>          EQUED is CHARACTER*1
!>          Specifies the form of equilibration that was done.
!>          = 'N':  No equilibration (always true if FACT = 'N').
!>          = 'Y':  Equilibration was done, i.e., A has been replaced by
!>                  diag(S) * A * diag(S).
!>          EQUED is an input argument if FACT = 'F'; otherwise, it is an
!>          output argument.
!> 
[in,out]S
!>          S is REAL array, dimension (N)
!>          The scale factors for A; not accessed if EQUED = 'N'.  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.
!> 
[in,out]B
!>          B is COMPLEX array, dimension (LDB,NRHS)
!>          On entry, the N-by-NRHS righthand side matrix B.
!>          On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
!>          B is overwritten by diag(S) * B.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]X
!>          X is COMPLEX array, dimension (LDX,NRHS)
!>          If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
!>          the original system of equations.  Note that if EQUED = 'Y',
!>          A and B are modified on exit, and the solution to the
!>          equilibrated system is inv(diag(S))*X.
!> 
[in]LDX
!>          LDX is INTEGER
!>          The leading dimension of the array X.  LDX >= max(1,N).
!> 
[out]RCOND
!>          RCOND is REAL
!>          The estimate of the reciprocal condition number of the matrix
!>          A after equilibration (if done).  If RCOND is less than the
!>          machine precision (in particular, if RCOND = 0), the matrix
!>          is singular to working precision.  This condition is
!>          indicated by a return code of INFO > 0.
!> 
[out]FERR
!>          FERR is REAL array, dimension (NRHS)
!>          The estimated forward error bound for each solution vector
!>          X(j) (the j-th column of the solution matrix X).
!>          If XTRUE is the true solution corresponding to X(j), FERR(j)
!>          is an estimated upper bound for the magnitude of the largest
!>          element in (X(j) - XTRUE) divided by the magnitude of the
!>          largest element in X(j).  The estimate is as reliable as
!>          the estimate for RCOND, and is almost always a slight
!>          overestimate of the true error.
!> 
[out]BERR
!>          BERR is REAL array, dimension (NRHS)
!>          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).
!> 
[out]WORK
!>          WORK is COMPLEX array, dimension (2*N)
!> 
[out]RWORK
!>          RWORK is REAL array, dimension (N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!>          > 0: if INFO = i, and i is
!>                <= N:  the leading minor of order i of A is
!>                       not positive definite, so the factorization
!>                       could not be completed, and the solution has not
!>                       been computed. RCOND = 0 is returned.
!>                = N+1: U is nonsingular, but RCOND is less than machine
!>                       precision, meaning that the matrix is singular
!>                       to working precision.  Nevertheless, the
!>                       solution and error bounds are computed because
!>                       there are a number of situations where the
!>                       computed solution can be more accurate than the
!>                       value of RCOND would suggest.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 303 of file cposvx.f.

306*
307* -- LAPACK driver routine --
308* -- LAPACK is a software package provided by Univ. of Tennessee, --
309* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
310*
311* .. Scalar Arguments ..
312 CHARACTER EQUED, FACT, UPLO
313 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
314 REAL RCOND
315* ..
316* .. Array Arguments ..
317 REAL BERR( * ), FERR( * ), RWORK( * ), S( * )
318 COMPLEX A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
319 $ WORK( * ), X( LDX, * )
320* ..
321*
322* =====================================================================
323*
324* .. Parameters ..
325 REAL ZERO, ONE
326 parameter( zero = 0.0e+0, one = 1.0e+0 )
327* ..
328* .. Local Scalars ..
329 LOGICAL EQUIL, NOFACT, RCEQU
330 INTEGER I, INFEQU, J
331 REAL AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
332* ..
333* .. External Functions ..
334 LOGICAL LSAME
335 REAL CLANHE, SLAMCH
336 EXTERNAL lsame, clanhe, slamch
337* ..
338* .. External Subroutines ..
339 EXTERNAL clacpy, claqhe, cpocon, cpoequ, cporfs, cpotrf,
340 $ cpotrs, xerbla
341* ..
342* .. Intrinsic Functions ..
343 INTRINSIC max, min
344* ..
345* .. Executable Statements ..
346*
347 info = 0
348 nofact = lsame( fact, 'N' )
349 equil = lsame( fact, 'E' )
350 IF( nofact .OR. equil ) THEN
351 equed = 'N'
352 rcequ = .false.
353 ELSE
354 rcequ = lsame( equed, 'Y' )
355 smlnum = slamch( 'Safe minimum' )
356 bignum = one / smlnum
357 END IF
358*
359* Test the input parameters.
360*
361 IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
362 $ THEN
363 info = -1
364 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) )
365 $ THEN
366 info = -2
367 ELSE IF( n.LT.0 ) THEN
368 info = -3
369 ELSE IF( nrhs.LT.0 ) THEN
370 info = -4
371 ELSE IF( lda.LT.max( 1, n ) ) THEN
372 info = -6
373 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
374 info = -8
375 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
376 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
377 info = -9
378 ELSE
379 IF( rcequ ) THEN
380 smin = bignum
381 smax = zero
382 DO 10 j = 1, n
383 smin = min( smin, s( j ) )
384 smax = max( smax, s( j ) )
385 10 CONTINUE
386 IF( smin.LE.zero ) THEN
387 info = -10
388 ELSE IF( n.GT.0 ) THEN
389 scond = max( smin, smlnum ) / min( smax, bignum )
390 ELSE
391 scond = one
392 END IF
393 END IF
394 IF( info.EQ.0 ) THEN
395 IF( ldb.LT.max( 1, n ) ) THEN
396 info = -12
397 ELSE IF( ldx.LT.max( 1, n ) ) THEN
398 info = -14
399 END IF
400 END IF
401 END IF
402*
403 IF( info.NE.0 ) THEN
404 CALL xerbla( 'CPOSVX', -info )
405 RETURN
406 END IF
407*
408 IF( equil ) THEN
409*
410* Compute row and column scalings to equilibrate the matrix A.
411*
412 CALL cpoequ( n, a, lda, s, scond, amax, infequ )
413 IF( infequ.EQ.0 ) THEN
414*
415* Equilibrate the matrix.
416*
417 CALL claqhe( uplo, n, a, lda, s, scond, amax, equed )
418 rcequ = lsame( equed, 'Y' )
419 END IF
420 END IF
421*
422* Scale the right hand side.
423*
424 IF( rcequ ) THEN
425 DO 30 j = 1, nrhs
426 DO 20 i = 1, n
427 b( i, j ) = s( i )*b( i, j )
428 20 CONTINUE
429 30 CONTINUE
430 END IF
431*
432 IF( nofact .OR. equil ) THEN
433*
434* Compute the Cholesky factorization A = U**H *U or A = L*L**H.
435*
436 CALL clacpy( uplo, n, n, a, lda, af, ldaf )
437 CALL cpotrf( uplo, n, af, ldaf, info )
438*
439* Return if INFO is non-zero.
440*
441 IF( info.GT.0 )THEN
442 rcond = zero
443 RETURN
444 END IF
445 END IF
446*
447* Compute the norm of the matrix A.
448*
449 anorm = clanhe( '1', uplo, n, a, lda, rwork )
450*
451* Compute the reciprocal of the condition number of A.
452*
453 CALL cpocon( uplo, n, af, ldaf, anorm, rcond, work, rwork, info )
454*
455* Compute the solution matrix X.
456*
457 CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
458 CALL cpotrs( uplo, n, nrhs, af, ldaf, x, ldx, info )
459*
460* Use iterative refinement to improve the computed solution and
461* compute error bounds and backward error estimates for it.
462*
463 CALL cporfs( uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx,
464 $ ferr, berr, work, rwork, info )
465*
466* Transform the solution matrix X to a solution of the original
467* system.
468*
469 IF( rcequ ) THEN
470 DO 50 j = 1, nrhs
471 DO 40 i = 1, n
472 x( i, j ) = s( i )*x( i, j )
473 40 CONTINUE
474 50 CONTINUE
475 DO 60 j = 1, nrhs
476 ferr( j ) = ferr( j ) / scond
477 60 CONTINUE
478 END IF
479*
480* Set INFO = N+1 if the matrix is singular to working precision.
481*
482 IF( rcond.LT.slamch( 'Epsilon' ) )
483 $ info = n + 1
484*
485 RETURN
486*
487* End of CPOSVX
488*
real function clanhe(norm, uplo, n, a, lda, work)
CLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clanhe.f:124
subroutine claqhe(uplo, n, a, lda, s, scond, amax, equed)
CLAQHE scales a Hermitian matrix.
Definition claqhe.f:134
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
subroutine cporfs(uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CPORFS
Definition cporfs.f:183
subroutine cpocon(uplo, n, a, lda, anorm, rcond, work, rwork, info)
CPOCON
Definition cpocon.f:121
subroutine cpoequ(n, a, lda, s, scond, amax, info)
CPOEQU
Definition cpoequ.f:113
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
#define min(a, b)
Definition macros.h:20

◆ cposvxx()

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

CPOSVXX computes the solution to system of linear equations A * X = B for PO matrices

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

Purpose:
!>
!>    CPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
!>    to compute the solution to a complex system of linear equations
!>    A * X = B, where A is an N-by-N Hermitian positive definite matrix
!>    and X and B are N-by-NRHS matrices.
!>
!>    If requested, both normwise and maximum componentwise error bounds
!>    are returned. CPOSVXX will return a solution with a tiny
!>    guaranteed error (O(eps) where eps is the working machine
!>    precision) unless the matrix is very ill-conditioned, in which
!>    case a warning is returned. Relevant condition numbers also are
!>    calculated and returned.
!>
!>    CPOSVXX accepts user-provided factorizations and equilibration
!>    factors; see the definitions of the FACT and EQUED options.
!>    Solving with refinement and using a factorization from a previous
!>    CPOSVXX call will also produce a solution with either O(eps)
!>    errors or warnings, but we cannot make that claim for general
!>    user-provided factorizations and equilibration factors if they
!>    differ from what CPOSVXX would itself produce.
!> 
Description:
!>
!>    The following steps are performed:
!>
!>    1. If FACT = 'E', real scaling factors are computed to equilibrate
!>    the system:
!>
!>      diag(S)*A*diag(S)     *inv(diag(S))*X = diag(S)*B
!>
!>    Whether or not the system will be equilibrated depends on the
!>    scaling of the matrix A, but if equilibration is used, A is
!>    overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
!>
!>    2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
!>    factor the matrix A (after equilibration if FACT = 'E') as
!>       A = U**T* U,  if UPLO = 'U', or
!>       A = L * L**T,  if UPLO = 'L',
!>    where U is an upper triangular matrix and L is a lower triangular
!>    matrix.
!>
!>    3. If the leading i-by-i principal minor is not positive definite,
!>    then the routine returns with INFO = i. Otherwise, the factored
!>    form of A is used to estimate the condition number of the matrix
!>    A (see argument RCOND).  If the reciprocal of the condition number
!>    is less than machine precision, the routine still goes on to solve
!>    for X and compute error bounds as described below.
!>
!>    4. The system of equations is solved for X using the factored form
!>    of A.
!>
!>    5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
!>    the routine will use iterative refinement to try to get a small
!>    error and error bounds.  Refinement calculates the residual to at
!>    least twice the working precision.
!>
!>    6. If equilibration was used, the matrix X is premultiplied by
!>    diag(S) so that it solves the original system before
!>    equilibration.
!> 
!>     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]FACT
!>          FACT is CHARACTER*1
!>     Specifies whether or not the factored form of the matrix A is
!>     supplied on entry, and if not, whether the matrix A should be
!>     equilibrated before it is factored.
!>       = 'F':  On entry, AF contains the factored form of A.
!>               If EQUED is not 'N', the matrix A has been
!>               equilibrated with scaling factors given by S.
!>               A and AF are not modified.
!>       = 'N':  The matrix A will be copied to AF and factored.
!>       = 'E':  The matrix A will be equilibrated if necessary, then
!>               copied to AF and factored.
!> 
[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 matrices B and X.  NRHS >= 0.
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA,N)
!>     On entry, the Hermitian matrix A, except if FACT = 'F' and EQUED =
!>     'Y', then A must contain the equilibrated matrix
!>     diag(S)*A*diag(S).  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.  A is
!>     not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
!>     'N' on exit.
!>
!>     On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
!>     diag(S)*A*diag(S).
!> 
[in]LDA
!>          LDA is INTEGER
!>     The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in,out]AF
!>          AF is COMPLEX array, dimension (LDAF,N)
!>     If FACT = 'F', then AF is an input argument and on entry
!>     contains the triangular factor U or L from the Cholesky
!>     factorization A = U**T*U or A = L*L**T, in the same storage
!>     format as A.  If EQUED .ne. 'N', then AF is the factored
!>     form of the equilibrated matrix diag(S)*A*diag(S).
!>
!>     If FACT = 'N', then AF is an output argument and on exit
!>     returns the triangular factor U or L from the Cholesky
!>     factorization A = U**T*U or A = L*L**T of the original
!>     matrix A.
!>
!>     If FACT = 'E', then AF is an output argument and on exit
!>     returns the triangular factor U or L from the Cholesky
!>     factorization A = U**T*U or A = L*L**T of the equilibrated
!>     matrix A (see the description of A for the form of the
!>     equilibrated matrix).
!> 
[in]LDAF
!>          LDAF is INTEGER
!>     The leading dimension of the array AF.  LDAF >= max(1,N).
!> 
[in,out]EQUED
!>          EQUED is CHARACTER*1
!>     Specifies the form of equilibration that was done.
!>       = 'N':  No equilibration (always true if FACT = 'N').
!>       = 'Y':  Both row and column equilibration, i.e., A has been
!>               replaced by diag(S) * A * diag(S).
!>     EQUED is an input argument if FACT = 'F'; otherwise, it is an
!>     output argument.
!> 
[in,out]S
!>          S is REAL array, dimension (N)
!>     The row 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,out]B
!>          B is COMPLEX array, dimension (LDB,NRHS)
!>     On entry, the N-by-NRHS right hand side matrix B.
!>     On exit,
!>     if EQUED = 'N', B is not modified;
!>     if EQUED = 'Y', B is overwritten by diag(S)*B;
!> 
[in]LDB
!>          LDB is INTEGER
!>     The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]X
!>          X is COMPLEX array, dimension (LDX,NRHS)
!>     If INFO = 0, the N-by-NRHS solution matrix X to the original
!>     system of equations.  Note that A and B are modified on exit if
!>     EQUED .ne. 'N', and the solution to the equilibrated system is
!>     inv(diag(S))*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]RPVGRW
!>          RPVGRW is REAL
!>     Reciprocal pivot growth.  On exit, this contains the reciprocal
!>     pivot growth factor norm(A)/norm(U). The 
!>     norm is used.  If this is much less than 1, then the stability of
!>     the LU factorization of the (equilibrated) matrix A could be poor.
!>     This also means that the solution X, estimated condition numbers,
!>     and error bounds could be unreliable. If factorization fails with
!>     0<INFO<=N, then this contains the reciprocal pivot growth factor
!>     for the leading INFO columns of A.
!> 
[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 COMPLEX array, dimension (2*N)
!> 
[out]RWORK
!>          RWORK is REAL array, dimension (2*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 492 of file cposvxx.f.

496*
497* -- LAPACK driver routine --
498* -- LAPACK is a software package provided by Univ. of Tennessee, --
499* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
500*
501* .. Scalar Arguments ..
502 CHARACTER EQUED, FACT, UPLO
503 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
504 $ N_ERR_BNDS
505 REAL RCOND, RPVGRW
506* ..
507* .. Array Arguments ..
508 COMPLEX A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
509 $ WORK( * ), X( LDX, * )
510 REAL S( * ), PARAMS( * ), BERR( * ), RWORK( * ),
511 $ ERR_BNDS_NORM( NRHS, * ),
512 $ ERR_BNDS_COMP( NRHS, * )
513* ..
514*
515* ==================================================================
516*
517* .. Parameters ..
518 REAL ZERO, ONE
519 parameter( zero = 0.0e+0, one = 1.0e+0 )
520 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
521 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
522 INTEGER CMP_ERR_I, PIV_GROWTH_I
523 parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
524 $ berr_i = 3 )
525 parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
526 parameter( cmp_rcond_i = 7, cmp_err_i = 8,
527 $ piv_growth_i = 9 )
528* ..
529* .. Local Scalars ..
530 LOGICAL EQUIL, NOFACT, RCEQU
531 INTEGER INFEQU, J
532 REAL AMAX, BIGNUM, SMIN, SMAX, SCOND, SMLNUM
533* ..
534* .. External Functions ..
535 EXTERNAL lsame, slamch, cla_porpvgrw
536 LOGICAL LSAME
537 REAL SLAMCH, CLA_PORPVGRW
538* ..
539* .. External Subroutines ..
540 EXTERNAL cpoequb, cpotrf, cpotrs, clacpy,
542* ..
543* .. Intrinsic Functions ..
544 INTRINSIC max, min
545* ..
546* .. Executable Statements ..
547*
548 info = 0
549 nofact = lsame( fact, 'N' )
550 equil = lsame( fact, 'E' )
551 smlnum = slamch( 'Safe minimum' )
552 bignum = one / smlnum
553 IF( nofact .OR. equil ) THEN
554 equed = 'N'
555 rcequ = .false.
556 ELSE
557 rcequ = lsame( equed, 'Y' )
558 ENDIF
559*
560* Default is failure. If an input parameter is wrong or
561* factorization fails, make everything look horrible. Only the
562* pivot growth is set here, the rest is initialized in CPORFSX.
563*
564 rpvgrw = zero
565*
566* Test the input parameters. PARAMS is not tested until CPORFSX.
567*
568 IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.
569 $ lsame( fact, 'F' ) ) THEN
570 info = -1
571 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
572 $ .NOT.lsame( uplo, 'L' ) ) THEN
573 info = -2
574 ELSE IF( n.LT.0 ) THEN
575 info = -3
576 ELSE IF( nrhs.LT.0 ) THEN
577 info = -4
578 ELSE IF( lda.LT.max( 1, n ) ) THEN
579 info = -6
580 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
581 info = -8
582 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
583 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
584 info = -9
585 ELSE
586 IF ( rcequ ) THEN
587 smin = bignum
588 smax = zero
589 DO 10 j = 1, n
590 smin = min( smin, s( j ) )
591 smax = max( smax, s( j ) )
592 10 CONTINUE
593 IF( smin.LE.zero ) THEN
594 info = -10
595 ELSE IF( n.GT.0 ) THEN
596 scond = max( smin, smlnum ) / min( smax, bignum )
597 ELSE
598 scond = one
599 END IF
600 END IF
601 IF( info.EQ.0 ) THEN
602 IF( ldb.LT.max( 1, n ) ) THEN
603 info = -12
604 ELSE IF( ldx.LT.max( 1, n ) ) THEN
605 info = -14
606 END IF
607 END IF
608 END IF
609*
610 IF( info.NE.0 ) THEN
611 CALL xerbla( 'CPOSVXX', -info )
612 RETURN
613 END IF
614*
615 IF( equil ) THEN
616*
617* Compute row and column scalings to equilibrate the matrix A.
618*
619 CALL cpoequb( n, a, lda, s, scond, amax, infequ )
620 IF( infequ.EQ.0 ) THEN
621*
622* Equilibrate the matrix.
623*
624 CALL claqhe( uplo, n, a, lda, s, scond, amax, equed )
625 rcequ = lsame( equed, 'Y' )
626 END IF
627 END IF
628*
629* Scale the right-hand side.
630*
631 IF( rcequ ) CALL clascl2( n, nrhs, s, b, ldb )
632*
633 IF( nofact .OR. equil ) THEN
634*
635* Compute the Cholesky factorization of A.
636*
637 CALL clacpy( uplo, n, n, a, lda, af, ldaf )
638 CALL cpotrf( uplo, n, af, ldaf, info )
639*
640* Return if INFO is non-zero.
641*
642 IF( info.GT.0 ) THEN
643*
644* Pivot in column INFO is exactly 0
645* Compute the reciprocal pivot growth factor of the
646* leading rank-deficient INFO columns of A.
647*
648 rpvgrw = cla_porpvgrw( uplo, n, a, lda, af, ldaf, rwork )
649 RETURN
650 END IF
651 END IF
652*
653* Compute the reciprocal pivot growth factor RPVGRW.
654*
655 rpvgrw = cla_porpvgrw( uplo, n, a, lda, af, ldaf, rwork )
656*
657* Compute the solution matrix X.
658*
659 CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
660 CALL cpotrs( uplo, n, nrhs, af, ldaf, x, ldx, info )
661*
662* Use iterative refinement to improve the computed solution and
663* compute error bounds and backward error estimates for it.
664*
665 CALL cporfsx( uplo, equed, n, nrhs, a, lda, af, ldaf,
666 $ s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm,
667 $ err_bnds_comp, nparams, params, work, rwork, info )
668
669*
670* Scale solutions.
671*
672 IF ( rcequ ) THEN
673 CALL clascl2( n, nrhs, s, x, ldx )
674 END IF
675*
676 RETURN
677*
678* End of CPOSVXX
679*
subroutine clascl2(m, n, d, x, ldx)
CLASCL2 performs diagonal scaling on a vector.
Definition clascl2.f:91
subroutine cpoequb(n, a, lda, s, scond, amax, info)
CPOEQUB
Definition cpoequb.f:119
subroutine cporfsx(uplo, equed, n, nrhs, a, lda, af, ldaf, s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, rwork, info)
CPORFSX
Definition cporfsx.f:393
real function cla_porpvgrw(uplo, ncols, a, lda, af, ldaf, work)
CLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian...