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

◆ 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 303 of file zposvx.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 DOUBLE PRECISION RCOND
315* ..
316* .. Array Arguments ..
317 DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * ), S( * )
318 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
319 $ WORK( * ), X( LDX, * )
320* ..
321*
322* =====================================================================
323*
324* .. Parameters ..
325 DOUBLE PRECISION ZERO, ONE
326 parameter( zero = 0.0d+0, one = 1.0d+0 )
327* ..
328* .. Local Scalars ..
329 LOGICAL EQUIL, NOFACT, RCEQU
330 INTEGER I, INFEQU, J
331 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
332* ..
333* .. External Functions ..
334 LOGICAL LSAME
335 DOUBLE PRECISION DLAMCH, ZLANHE
336 EXTERNAL lsame, dlamch, zlanhe
337* ..
338* .. External Subroutines ..
339 EXTERNAL xerbla, zlacpy, zlaqhe, zpocon, zpoequ, 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. .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( 'ZPOSVX', -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 zpoequ( n, a, lda, s, scond, amax, infequ )
413 IF( infequ.EQ.0 ) THEN
414*
415* Equilibrate the matrix.
416*
417 CALL zlaqhe( 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 zlacpy( uplo, n, n, a, lda, af, ldaf )
437 CALL zpotrf( 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 = zlanhe( '1', uplo, n, a, lda, rwork )
450*
451* Compute the reciprocal of the condition number of A.
452*
453 CALL zpocon( uplo, n, af, ldaf, anorm, rcond, work, rwork, info )
454*
455* Compute the solution matrix X.
456*
457 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
458 CALL zpotrs( 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 zporfs( 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.dlamch( 'Epsilon' ) )
483 $ info = n + 1
484*
485 RETURN
486*
487* End of ZPOSVX
488*
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:103
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:124
subroutine zlaqhe(uplo, n, a, lda, s, scond, amax, equed)
ZLAQHE scales a Hermitian matrix.
Definition zlaqhe.f:134
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:121
subroutine zpoequ(n, a, lda, s, scond, amax, info)
ZPOEQU
Definition zpoequ.f:113
subroutine zporfs(uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZPORFS
Definition zporfs.f:183
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:107
subroutine zpotrs(uplo, n, nrhs, a, lda, b, ldb, info)
ZPOTRS
Definition zpotrs.f:110
Here is the call graph for this function:
Here is the caller graph for this function: