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

◆ zposvx()

subroutine zposvx ( character fact,
character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldaf, * ) af,
integer ldaf,
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 )

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

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

Purpose:
!> !> ZPOSVX 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 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, 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*16 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*16 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 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 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*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.

Definition at line 301 of file zposvx.f.

305*
306* -- LAPACK driver routine --
307* -- LAPACK is a software package provided by Univ. of Tennessee, --
308* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
309*
310* .. Scalar Arguments ..
311 CHARACTER EQUED, FACT, UPLO
312 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
313 DOUBLE PRECISION RCOND
314* ..
315* .. Array Arguments ..
316 DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * ), S( * )
317 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
318 $ WORK( * ), X( LDX, * )
319* ..
320*
321* =====================================================================
322*
323* .. Parameters ..
324 DOUBLE PRECISION ZERO, ONE
325 parameter( zero = 0.0d+0, one = 1.0d+0 )
326* ..
327* .. Local Scalars ..
328 LOGICAL EQUIL, NOFACT, RCEQU
329 INTEGER I, INFEQU, J
330 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
331* ..
332* .. External Functions ..
333 LOGICAL LSAME
334 DOUBLE PRECISION DLAMCH, ZLANHE
335 EXTERNAL lsame, dlamch, zlanhe
336* ..
337* .. External Subroutines ..
338 EXTERNAL xerbla, zlacpy, zlaqhe, zpocon, zpoequ,
339 $ zporfs,
340 $ zpotrf, zpotrs
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 = dlamch( 'Safe minimum' )
356 bignum = one / smlnum
357 END IF
358*
359* Test the input parameters.
360*
361 IF( .NOT.nofact .AND.
362 $ .NOT.equil .AND.
363 $ .NOT.lsame( fact, 'F' ) )
364 $ THEN
365 info = -1
366 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
367 $ .NOT.lsame( uplo, 'L' ) )
368 $ THEN
369 info = -2
370 ELSE IF( n.LT.0 ) THEN
371 info = -3
372 ELSE IF( nrhs.LT.0 ) THEN
373 info = -4
374 ELSE IF( lda.LT.max( 1, n ) ) THEN
375 info = -6
376 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
377 info = -8
378 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
379 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
380 info = -9
381 ELSE
382 IF( rcequ ) THEN
383 smin = bignum
384 smax = zero
385 DO 10 j = 1, n
386 smin = min( smin, s( j ) )
387 smax = max( smax, s( j ) )
388 10 CONTINUE
389 IF( smin.LE.zero ) THEN
390 info = -10
391 ELSE IF( n.GT.0 ) THEN
392 scond = max( smin, smlnum ) / min( smax, bignum )
393 ELSE
394 scond = one
395 END IF
396 END IF
397 IF( info.EQ.0 ) THEN
398 IF( ldb.LT.max( 1, n ) ) THEN
399 info = -12
400 ELSE IF( ldx.LT.max( 1, n ) ) THEN
401 info = -14
402 END IF
403 END IF
404 END IF
405*
406 IF( info.NE.0 ) THEN
407 CALL xerbla( 'ZPOSVX', -info )
408 RETURN
409 END IF
410*
411 IF( equil ) THEN
412*
413* Compute row and column scalings to equilibrate the matrix A.
414*
415 CALL zpoequ( n, a, lda, s, scond, amax, infequ )
416 IF( infequ.EQ.0 ) THEN
417*
418* Equilibrate the matrix.
419*
420 CALL zlaqhe( uplo, n, a, lda, s, scond, amax, equed )
421 rcequ = lsame( equed, 'Y' )
422 END IF
423 END IF
424*
425* Scale the right hand side.
426*
427 IF( rcequ ) THEN
428 DO 30 j = 1, nrhs
429 DO 20 i = 1, n
430 b( i, j ) = s( i )*b( i, j )
431 20 CONTINUE
432 30 CONTINUE
433 END IF
434*
435 IF( nofact .OR. equil ) THEN
436*
437* Compute the Cholesky factorization A = U**H *U or A = L*L**H.
438*
439 CALL zlacpy( uplo, n, n, a, lda, af, ldaf )
440 CALL zpotrf( uplo, n, af, ldaf, info )
441*
442* Return if INFO is non-zero.
443*
444 IF( info.GT.0 )THEN
445 rcond = zero
446 RETURN
447 END IF
448 END IF
449*
450* Compute the norm of the matrix A.
451*
452 anorm = zlanhe( '1', uplo, n, a, lda, rwork )
453*
454* Compute the reciprocal of the condition number of A.
455*
456 CALL zpocon( uplo, n, af, ldaf, anorm, rcond, work, rwork,
457 $ info )
458*
459* Compute the solution matrix X.
460*
461 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
462 CALL zpotrs( uplo, n, nrhs, af, ldaf, 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 zporfs( uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx,
468 $ ferr, berr, work, rwork, info )
469*
470* Transform the solution matrix X to a solution of the original
471* system.
472*
473 IF( rcequ ) THEN
474 DO 50 j = 1, nrhs
475 DO 40 i = 1, n
476 x( i, j ) = s( i )*x( i, j )
477 40 CONTINUE
478 50 CONTINUE
479 DO 60 j = 1, nrhs
480 ferr( j ) = ferr( j ) / scond
481 60 CONTINUE
482 END IF
483*
484* Set INFO = N+1 if the matrix is singular to working precision.
485*
486 IF( rcond.LT.dlamch( 'Epsilon' ) )
487 $ info = n + 1
488*
489 RETURN
490*
491* End of ZPOSVX
492*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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 zlanhe(norm, uplo, n, a, lda, work)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhe.f:122
subroutine zlaqhe(uplo, n, a, lda, s, scond, amax, equed)
ZLAQHE scales a Hermitian matrix.
Definition zlaqhe.f:132
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zpocon(uplo, n, a, lda, anorm, rcond, work, rwork, info)
ZPOCON
Definition zpocon.f:119
subroutine zpoequ(n, a, lda, s, scond, amax, info)
ZPOEQU
Definition zpoequ.f:111
subroutine zporfs(uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZPORFS
Definition zporfs.f:181
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:105
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: