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

◆ zppsvx()

subroutine zppsvx ( character fact,
character uplo,
integer n,
integer nrhs,
complex*16, dimension( * ) ap,
complex*16, dimension( * ) afp,
character equed,
double precision, dimension( * ) s,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( ldx, * ) x,
integer ldx,
double precision rcond,
double precision, dimension( * ) ferr,
double precision, dimension( * ) berr,
complex*16, dimension( * ) work,
double precision, dimension( * ) rwork,
integer info )

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

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

Purpose:
!> !> ZPPSVX 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 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 = '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, L is a lower triangular !> matrix, and **H indicates conjugate transpose. !> !> 3. If the leading principal minor of order i is not positive, !> 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, AFP contains the factored form of A. !> If EQUED = 'Y', the matrix A has been equilibrated !> with scaling factors given by S. AP and AFP will not !> be modified. !> = 'N': The matrix A will be copied to AFP and factored. !> = 'E': The matrix A will be equilibrated if necessary, then !> 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,out]AP
!> AP is COMPLEX*16 array, dimension (N*(N+1)/2) !> On entry, the upper or lower triangle of the Hermitian matrix !> A, packed columnwise in a linear array, except if FACT = 'F' !> and EQUED = 'Y', then A must contain the equilibrated matrix !> diag(S)*A*diag(S). 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)*(2n-j)/2) = A(i,j) for j<=i<=n. !> See below for further details. 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,out]AFP
!> AFP is COMPLEX*16 array, dimension (N*(N+1)/2) !> If FACT = 'F', then AFP 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 AFP is the factored !> form of the equilibrated matrix A. !> !> If FACT = 'N', then AFP 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 AFP 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 AP for the form of the !> equilibrated matrix). !>
[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 DOUBLE PRECISION 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*16 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*16 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 DOUBLE PRECISION !> 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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*16 array, dimension (2*N) !>
[out]RWORK
!> RWORK is DOUBLE PRECISION 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 principal minor of order i of A !> is not positive, 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.
Further Details:
!> !> The packed storage scheme is illustrated by the following example !> when N = 4, UPLO = 'U': !> !> Two-dimensional storage of the Hermitian matrix A: !> !> a11 a12 a13 a14 !> a22 a23 a24 !> a33 a34 (aij = conjg(aji)) !> a44 !> !> Packed storage of the upper triangle of A: !> !> AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ] !>

Definition at line 307 of file zppsvx.f.

310*
311* -- LAPACK driver routine --
312* -- LAPACK is a software package provided by Univ. of Tennessee, --
313* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
314*
315* .. Scalar Arguments ..
316 CHARACTER EQUED, FACT, UPLO
317 INTEGER INFO, LDB, LDX, N, NRHS
318 DOUBLE PRECISION RCOND
319* ..
320* .. Array Arguments ..
321 DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * ), S( * )
322 COMPLEX*16 AFP( * ), AP( * ), B( LDB, * ), WORK( * ),
323 $ X( LDX, * )
324* ..
325*
326* =====================================================================
327*
328* .. Parameters ..
329 DOUBLE PRECISION ZERO, ONE
330 parameter( zero = 0.0d+0, one = 1.0d+0 )
331* ..
332* .. Local Scalars ..
333 LOGICAL EQUIL, NOFACT, RCEQU
334 INTEGER I, INFEQU, J
335 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
336* ..
337* .. External Functions ..
338 LOGICAL LSAME
339 DOUBLE PRECISION DLAMCH, ZLANHP
340 EXTERNAL lsame, dlamch, zlanhp
341* ..
342* .. External Subroutines ..
343 EXTERNAL xerbla, zcopy, zlacpy, zlaqhp, zppcon,
344 $ zppequ,
346* ..
347* .. Intrinsic Functions ..
348 INTRINSIC max, min
349* ..
350* .. Executable Statements ..
351*
352 info = 0
353 nofact = lsame( fact, 'N' )
354 equil = lsame( fact, 'E' )
355 IF( nofact .OR. equil ) THEN
356 equed = 'N'
357 rcequ = .false.
358 ELSE
359 rcequ = lsame( equed, 'Y' )
360 smlnum = dlamch( 'Safe minimum' )
361 bignum = one / smlnum
362 END IF
363*
364* Test the input parameters.
365*
366 IF( .NOT.nofact .AND.
367 $ .NOT.equil .AND.
368 $ .NOT.lsame( fact, 'F' ) )
369 $ THEN
370 info = -1
371 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
372 $ .NOT.lsame( uplo, 'L' ) )
373 $ THEN
374 info = -2
375 ELSE IF( n.LT.0 ) THEN
376 info = -3
377 ELSE IF( nrhs.LT.0 ) THEN
378 info = -4
379 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
380 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
381 info = -7
382 ELSE
383 IF( rcequ ) THEN
384 smin = bignum
385 smax = zero
386 DO 10 j = 1, n
387 smin = min( smin, s( j ) )
388 smax = max( smax, s( j ) )
389 10 CONTINUE
390 IF( smin.LE.zero ) THEN
391 info = -8
392 ELSE IF( n.GT.0 ) THEN
393 scond = max( smin, smlnum ) / min( smax, bignum )
394 ELSE
395 scond = one
396 END IF
397 END IF
398 IF( info.EQ.0 ) THEN
399 IF( ldb.LT.max( 1, n ) ) THEN
400 info = -10
401 ELSE IF( ldx.LT.max( 1, n ) ) THEN
402 info = -12
403 END IF
404 END IF
405 END IF
406*
407 IF( info.NE.0 ) THEN
408 CALL xerbla( 'ZPPSVX', -info )
409 RETURN
410 END IF
411*
412 IF( equil ) THEN
413*
414* Compute row and column scalings to equilibrate the matrix A.
415*
416 CALL zppequ( uplo, n, ap, s, scond, amax, infequ )
417 IF( infequ.EQ.0 ) THEN
418*
419* Equilibrate the matrix.
420*
421 CALL zlaqhp( uplo, n, ap, s, scond, amax, equed )
422 rcequ = lsame( equed, 'Y' )
423 END IF
424 END IF
425*
426* Scale the right-hand side.
427*
428 IF( rcequ ) THEN
429 DO 30 j = 1, nrhs
430 DO 20 i = 1, n
431 b( i, j ) = s( i )*b( i, j )
432 20 CONTINUE
433 30 CONTINUE
434 END IF
435*
436 IF( nofact .OR. equil ) THEN
437*
438* Compute the Cholesky factorization A = U**H * U or A = L * L**H.
439*
440 CALL zcopy( n*( n+1 ) / 2, ap, 1, afp, 1 )
441 CALL zpptrf( uplo, n, afp, info )
442*
443* Return if INFO is non-zero.
444*
445 IF( info.GT.0 )THEN
446 rcond = zero
447 RETURN
448 END IF
449 END IF
450*
451* Compute the norm of the matrix A.
452*
453 anorm = zlanhp( 'I', uplo, n, ap, rwork )
454*
455* Compute the reciprocal of the condition number of A.
456*
457 CALL zppcon( uplo, n, afp, anorm, rcond, work, rwork, info )
458*
459* Compute the solution matrix X.
460*
461 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
462 CALL zpptrs( uplo, n, nrhs, afp, x, ldx, info )
463*
464* Use iterative refinement to improve the computed solution and
465* compute error bounds and backward error estimates for it.
466*
467 CALL zpprfs( uplo, n, nrhs, ap, afp, b, ldb, x, ldx, ferr,
468 $ berr,
469 $ work, rwork, info )
470*
471* Transform the solution matrix X to a solution of the original
472* system.
473*
474 IF( rcequ ) THEN
475 DO 50 j = 1, nrhs
476 DO 40 i = 1, n
477 x( i, j ) = s( i )*x( i, j )
478 40 CONTINUE
479 50 CONTINUE
480 DO 60 j = 1, nrhs
481 ferr( j ) = ferr( j ) / scond
482 60 CONTINUE
483 END IF
484*
485* Set INFO = N+1 if the matrix is singular to working precision.
486*
487 IF( rcond.LT.dlamch( 'Epsilon' ) )
488 $ info = n + 1
489*
490 RETURN
491*
492* End of ZPPSVX
493*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlanhp(norm, uplo, n, ap, work)
ZLANHP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhp.f:115
subroutine zlaqhp(uplo, n, ap, s, scond, amax, equed)
ZLAQHP scales a Hermitian matrix stored in packed form.
Definition zlaqhp.f:124
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zppcon(uplo, n, ap, anorm, rcond, work, rwork, info)
ZPPCON
Definition zppcon.f:117
subroutine zppequ(uplo, n, ap, s, scond, amax, info)
ZPPEQU
Definition zppequ.f:115
subroutine zpprfs(uplo, n, nrhs, ap, afp, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZPPRFS
Definition zpprfs.f:170
subroutine zpptrf(uplo, n, ap, info)
ZPPTRF
Definition zpptrf.f:117
subroutine zpptrs(uplo, n, nrhs, ap, b, ldb, info)
ZPPTRS
Definition zpptrs.f:106
Here is the call graph for this function:
Here is the caller graph for this function: