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

◆ cspsvx()

subroutine cspsvx ( character fact,
character uplo,
integer n,
integer nrhs,
complex, dimension( * ) ap,
complex, dimension( * ) afp,
integer, dimension( * ) ipiv,
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 )

CSPSVX computes the solution to system of linear equations A * X = B for OTHER matrices

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

Purpose:
!>
!> CSPSVX uses the diagonal pivoting factorization A = U*D*U**T or
!> A = L*D*L**T to compute the solution to a complex system of linear
!> equations A * X = B, where A is an N-by-N symmetric matrix stored
!> in packed format 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 = 'N', the diagonal pivoting method is used to factor A as
!>       A = U * D * U**T,  if UPLO = 'U', or
!>       A = L * D * L**T,  if UPLO = 'L',
!>    where U (or L) is a product of permutation and unit upper (lower)
!>    triangular matrices and D is symmetric and block diagonal with
!>    1-by-1 and 2-by-2 diagonal blocks.
!>
!> 2. If some D(i,i)=0, so that D is exactly singular, 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.
!>
!> 3. The system of equations is solved for X using the factored form
!>    of A.
!>
!> 4. Iterative refinement is applied to improve the computed solution
!>    matrix and calculate error bounds and backward error estimates
!>    for it.
!> 
Parameters
[in]FACT
!>          FACT is CHARACTER*1
!>          Specifies whether or not the factored form of A has been
!>          supplied on entry.
!>          = 'F':  On entry, AFP and IPIV contain the factored form
!>                  of A.  AP, AFP and IPIV will not be modified.
!>          = 'N':  The matrix A will be copied to AFP 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]AP
!>          AP is COMPLEX array, dimension (N*(N+1)/2)
!>          The upper or lower triangle of the symmetric matrix A, packed
!>          columnwise in a linear array.  The j-th column of A is stored
!>          in the array AP as follows:
!>          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
!>          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.
!>          See below for further details.
!> 
[in,out]AFP
!>          AFP is COMPLEX array, dimension (N*(N+1)/2)
!>          If FACT = 'F', then AFP is an input argument and on entry
!>          contains the block diagonal matrix D and the multipliers used
!>          to obtain the factor U or L from the factorization
!>          A = U*D*U**T or A = L*D*L**T as computed by CSPTRF, stored as
!>          a packed triangular matrix in the same storage format as A.
!>
!>          If FACT = 'N', then AFP is an output argument and on exit
!>          contains the block diagonal matrix D and the multipliers used
!>          to obtain the factor U or L from the factorization
!>          A = U*D*U**T or A = L*D*L**T as computed by CSPTRF, stored as
!>          a packed triangular matrix in the same storage format as A.
!> 
[in,out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          If FACT = 'F', then IPIV is an input argument and on entry
!>          contains details of the interchanges and the block structure
!>          of D, as determined by CSPTRF.
!>          If IPIV(k) > 0, then rows and columns k and IPIV(k) were
!>          interchanged and D(k,k) is a 1-by-1 diagonal block.
!>          If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
!>          columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
!>          is a 2-by-2 diagonal block.  If UPLO = 'L' and IPIV(k) =
!>          IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
!>          interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
!>
!>          If FACT = 'N', then IPIV is an output argument and on exit
!>          contains details of the interchanges and the block structure
!>          of D, as determined by CSPTRF.
!> 
[in]B
!>          B is COMPLEX array, dimension (LDB,NRHS)
!>          The N-by-NRHS right hand side matrix 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.
!> 
[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.  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:  D(i,i) is exactly zero.  The factorization
!>                       has been completed but the factor D is exactly
!>                       singular, so the solution and error bounds could
!>                       not be computed. RCOND = 0 is returned.
!>                = N+1: D 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.
Further Details:
!>
!>  The packed storage scheme is illustrated by the following example
!>  when N = 4, UPLO = 'U':
!>
!>  Two-dimensional storage of the symmetric matrix A:
!>
!>     a11 a12 a13 a14
!>         a22 a23 a24
!>             a33 a34     (aij = aji)
!>                 a44
!>
!>  Packed storage of the upper triangle of A:
!>
!>  AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]
!> 

Definition at line 273 of file cspsvx.f.

276*
277* -- LAPACK driver routine --
278* -- LAPACK is a software package provided by Univ. of Tennessee, --
279* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
280*
281* .. Scalar Arguments ..
282 CHARACTER FACT, UPLO
283 INTEGER INFO, LDB, LDX, N, NRHS
284 REAL RCOND
285* ..
286* .. Array Arguments ..
287 INTEGER IPIV( * )
288 REAL BERR( * ), FERR( * ), RWORK( * )
289 COMPLEX AFP( * ), AP( * ), B( LDB, * ), WORK( * ),
290 $ X( LDX, * )
291* ..
292*
293* =====================================================================
294*
295* .. Parameters ..
296 REAL ZERO
297 parameter( zero = 0.0e+0 )
298* ..
299* .. Local Scalars ..
300 LOGICAL NOFACT
301 REAL ANORM
302* ..
303* .. External Functions ..
304 LOGICAL LSAME
305 REAL CLANSP, SLAMCH
306 EXTERNAL lsame, clansp, slamch
307* ..
308* .. External Subroutines ..
309 EXTERNAL ccopy, clacpy, cspcon, csprfs, csptrf,
310 $ csptrs,
311 $ xerbla
312* ..
313* .. Intrinsic Functions ..
314 INTRINSIC max
315* ..
316* .. Executable Statements ..
317*
318* Test the input parameters.
319*
320 info = 0
321 nofact = lsame( fact, 'N' )
322 IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
323 info = -1
324 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
325 $ .NOT.lsame( uplo, 'L' ) )
326 $ THEN
327 info = -2
328 ELSE IF( n.LT.0 ) THEN
329 info = -3
330 ELSE IF( nrhs.LT.0 ) THEN
331 info = -4
332 ELSE IF( ldb.LT.max( 1, n ) ) THEN
333 info = -9
334 ELSE IF( ldx.LT.max( 1, n ) ) THEN
335 info = -11
336 END IF
337 IF( info.NE.0 ) THEN
338 CALL xerbla( 'CSPSVX', -info )
339 RETURN
340 END IF
341*
342 IF( nofact ) THEN
343*
344* Compute the factorization A = U*D*U**T or A = L*D*L**T.
345*
346 CALL ccopy( n*( n+1 ) / 2, ap, 1, afp, 1 )
347 CALL csptrf( uplo, n, afp, ipiv, info )
348*
349* Return if INFO is non-zero.
350*
351 IF( info.GT.0 )THEN
352 rcond = zero
353 RETURN
354 END IF
355 END IF
356*
357* Compute the norm of the matrix A.
358*
359 anorm = clansp( 'I', uplo, n, ap, rwork )
360*
361* Compute the reciprocal of the condition number of A.
362*
363 CALL cspcon( uplo, n, afp, ipiv, anorm, rcond, work, info )
364*
365* Compute the solution vectors X.
366*
367 CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
368 CALL csptrs( uplo, n, nrhs, afp, ipiv, x, ldx, info )
369*
370* Use iterative refinement to improve the computed solutions and
371* compute error bounds and backward error estimates for them.
372*
373 CALL csprfs( uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx,
374 $ ferr,
375 $ berr, work, rwork, info )
376*
377* Set INFO = N+1 if the matrix is singular to working precision.
378*
379 IF( rcond.LT.slamch( 'Epsilon' ) )
380 $ info = n + 1
381*
382 RETURN
383*
384* End of CSPSVX
385*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine cspcon(uplo, n, ap, ipiv, anorm, rcond, work, info)
CSPCON
Definition cspcon.f:117
subroutine csprfs(uplo, n, nrhs, ap, afp, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CSPRFS
Definition csprfs.f:179
subroutine csptrf(uplo, n, ap, ipiv, info)
CSPTRF
Definition csptrf.f:156
subroutine csptrs(uplo, n, nrhs, ap, ipiv, b, ldb, info)
CSPTRS
Definition csptrs.f:113
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:101
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clansp(norm, uplo, n, ap, work)
CLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clansp.f:113
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: