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

Functions

subroutine dposv (uplo, n, nrhs, a, lda, b, ldb, info)
  DPOSV computes the solution to system of linear equations A * X = B for PO matrices
subroutine dposvx (fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info)
  DPOSVX computes the solution to system of linear equations A * X = B for PO matrices
subroutine dposvxx (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, iwork, info)
  DPOSVXX computes the solution to system of linear equations A * X = B for PO matrices
subroutine dsposv (uplo, n, nrhs, a, lda, b, ldb, x, ldx, work, swork, iter, info)
  DSPOSV computes the solution to system of linear equations A * X = B for PO matrices

Detailed Description

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

Function Documentation

◆ dposv()

subroutine dposv ( character uplo,
integer n,
integer nrhs,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( ldb, * ) b,
integer ldb,
integer info )

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

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

Purpose:
!>
!> DPOSV computes the solution to a real system of linear equations
!>    A * X = B,
!> where A is an N-by-N symmetric positive definite matrix and X and B
!> are N-by-NRHS matrices.
!>
!> The Cholesky decomposition is used to factor A 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.  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 DOUBLE PRECISION array, dimension (LDA,N)
!>          On entry, the symmetric 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**T*U or A = L*L**T.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in,out]B
!>          B is DOUBLE PRECISION 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 dposv.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 DOUBLE PRECISION A( LDA, * ), B( LDB, * )
141* ..
142*
143* =====================================================================
144*
145* .. External Functions ..
146 LOGICAL LSAME
147 EXTERNAL lsame
148* ..
149* .. External Subroutines ..
150 EXTERNAL dpotrf, dpotrs, 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( 'DPOSV ', -info )
173 RETURN
174 END IF
175*
176* Compute the Cholesky factorization A = U**T*U or A = L*L**T.
177*
178 CALL dpotrf( 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 dpotrs( uplo, n, nrhs, a, lda, b, ldb, info )
184*
185 END IF
186 RETURN
187*
188* End of DPOSV
189*
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
subroutine dpotrs(uplo, n, nrhs, a, lda, b, ldb, info)
DPOTRS
Definition dpotrs.f:110
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:107
#define max(a, b)
Definition macros.h:21

◆ dposvx()

subroutine dposvx ( character fact,
character uplo,
integer n,
integer nrhs,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( ldaf, * ) af,
integer ldaf,
character equed,
double precision, dimension( * ) s,
double precision, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( ldx, * ) x,
integer ldx,
double precision rcond,
double precision, dimension( * ) ferr,
double precision, dimension( * ) berr,
double precision, dimension( * ) work,
integer, dimension( * ) iwork,
integer info )

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

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

Purpose:
!>
!> DPOSVX uses the Cholesky factorization A = U**T*U or A = L*L**T to
!> compute the solution to a real system of linear equations
!>    A * X = B,
!> where A is an N-by-N symmetric 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**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.  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 DOUBLE PRECISION array, dimension (LDA,N)
!>          On entry, the symmetric 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 DOUBLE PRECISION 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':  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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (3*N)
!> 
[out]IWORK
!>          IWORK is INTEGER 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 304 of file dposvx.f.

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

◆ dposvxx()

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

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

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

Purpose:
!>
!>    DPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
!>    to compute the solution to a double precision system of linear equations
!>    A * X = B, where A is an N-by-N symmetric positive definite matrix
!>    and X and B are N-by-NRHS matrices.
!>
!>    If requested, both normwise and maximum componentwise error bounds
!>    are returned. DPOSVXX 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.
!>
!>    DPOSVXX 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
!>    DPOSVXX 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 DPOSVXX would itself produce.
!> 
Description:
!>
!>    The following steps are performed:
!>
!>    1. If FACT = 'E', double precision 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 DOUBLE PRECISION array, dimension (LDA,N)
!>     On entry, the symmetric 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION
!>     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 DOUBLE PRECISION
!>     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 DOUBLE PRECISION 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 DOUBLE PRECISION 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) * dlamch('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) * dlamch('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) * dlamch('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 DOUBLE PRECISION 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) * dlamch('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) * dlamch('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) * dlamch('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 DOUBLE PRECISION 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.0D+0
!>            = 0.0:  No refinement is performed, and no error bounds are
!>                    computed.
!>            = 1.0:  Use the extra-precise refinement algorithm.
!>              (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 DOUBLE PRECISION array, dimension (4*N)
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (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 490 of file dposvxx.f.

494*
495* -- LAPACK driver routine --
496* -- LAPACK is a software package provided by Univ. of Tennessee, --
497* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
498*
499* .. Scalar Arguments ..
500 CHARACTER EQUED, FACT, UPLO
501 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
502 $ N_ERR_BNDS
503 DOUBLE PRECISION RCOND, RPVGRW
504* ..
505* .. Array Arguments ..
506 INTEGER IWORK( * )
507 DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
508 $ X( LDX, * ), WORK( * )
509 DOUBLE PRECISION S( * ), PARAMS( * ), BERR( * ),
510 $ ERR_BNDS_NORM( NRHS, * ),
511 $ ERR_BNDS_COMP( NRHS, * )
512* ..
513*
514* ==================================================================
515*
516* .. Parameters ..
517 DOUBLE PRECISION ZERO, ONE
518 parameter( zero = 0.0d+0, one = 1.0d+0 )
519 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
520 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
521 INTEGER CMP_ERR_I, PIV_GROWTH_I
522 parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
523 $ berr_i = 3 )
524 parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
525 parameter( cmp_rcond_i = 7, cmp_err_i = 8,
526 $ piv_growth_i = 9 )
527* ..
528* .. Local Scalars ..
529 LOGICAL EQUIL, NOFACT, RCEQU
530 INTEGER INFEQU, J
531 DOUBLE PRECISION AMAX, BIGNUM, SMIN, SMAX,
532 $ SCOND, SMLNUM
533* ..
534* .. External Functions ..
535 EXTERNAL lsame, dlamch, dla_porpvgrw
536 LOGICAL LSAME
537 DOUBLE PRECISION DLAMCH, DLA_PORPVGRW
538* ..
539* .. External Subroutines ..
540 EXTERNAL dpoequb, dpotrf, dpotrs, dlacpy, dlaqsy,
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 = dlamch( '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 DPORFSX.
563*
564 rpvgrw = zero
565*
566* Test the input parameters. PARAMS is not tested until DPORFSX.
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( 'DPOSVXX', -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 dpoequb( n, a, lda, s, scond, amax, infequ )
620 IF( infequ.EQ.0 ) THEN
621*
622* Equilibrate the matrix.
623*
624 CALL dlaqsy( 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 dlascl2( n, nrhs, s, b, ldb )
632*
633 IF( nofact .OR. equil ) THEN
634*
635* Compute the Cholesky factorization of A.
636*
637 CALL dlacpy( uplo, n, n, a, lda, af, ldaf )
638 CALL dpotrf( uplo, n, af, ldaf, info )
639*
640* Return if INFO is non-zero.
641*
642 IF( info.NE.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 = dla_porpvgrw( uplo, info, a, lda, af, ldaf, work )
649 RETURN
650 ENDIF
651 END IF
652*
653* Compute the reciprocal growth factor RPVGRW.
654*
655 rpvgrw = dla_porpvgrw( uplo, n, a, lda, af, ldaf, work )
656*
657* Compute the solution matrix X.
658*
659 CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
660 CALL dpotrs( 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 dporfsx( 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, iwork, info )
668
669*
670* Scale solutions.
671*
672 IF ( rcequ ) THEN
673 CALL dlascl2 ( n, nrhs, s, x, ldx )
674 END IF
675*
676 RETURN
677*
678* End of DPOSVXX
679*
subroutine dlascl2(m, n, d, x, ldx)
DLASCL2 performs diagonal scaling on a vector.
Definition dlascl2.f:90
subroutine dporfsx(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, iwork, info)
DPORFSX
Definition dporfsx.f:394
subroutine dpoequb(n, a, lda, s, scond, amax, info)
DPOEQUB
Definition dpoequb.f:118
double precision function dla_porpvgrw(uplo, ncols, a, lda, af, ldaf, work)
DLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian...

◆ dsposv()

subroutine dsposv ( character uplo,
integer n,
integer nrhs,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( ldx, * ) x,
integer ldx,
double precision, dimension( n, * ) work,
real, dimension( * ) swork,
integer iter,
integer info )

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

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

Purpose:
!>
!> DSPOSV computes the solution to a real system of linear equations
!>    A * X = B,
!> where A is an N-by-N symmetric positive definite matrix and X and B
!> are N-by-NRHS matrices.
!>
!> DSPOSV first attempts to factorize the matrix in SINGLE PRECISION
!> and use this factorization within an iterative refinement procedure
!> to produce a solution with DOUBLE PRECISION normwise backward error
!> quality (see below). If the approach fails the method switches to a
!> DOUBLE PRECISION factorization and solve.
!>
!> The iterative refinement is not going to be a winning strategy if
!> the ratio SINGLE PRECISION performance over DOUBLE PRECISION
!> performance is too small. A reasonable strategy should take the
!> number of right-hand sides and the size of the matrix into account.
!> This might be done with a call to ILAENV in the future. Up to now, we
!> always try iterative refinement.
!>
!> The iterative refinement process is stopped if
!>     ITER > ITERMAX
!> or for all the RHS we have:
!>     RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX
!> where
!>     o ITER is the number of the current iteration in the iterative
!>       refinement process
!>     o RNRM is the infinity-norm of the residual
!>     o XNRM is the infinity-norm of the solution
!>     o ANRM is the infinity-operator-norm of the matrix A
!>     o EPS is the machine epsilon returned by DLAMCH('Epsilon')
!> The value ITERMAX and BWDMAX are fixed to 30 and 1.0D+00
!> respectively.
!> 
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 DOUBLE PRECISION array,
!>          dimension (LDA,N)
!>          On entry, the symmetric 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 iterative refinement has been successfully used
!>          (INFO = 0 and ITER >= 0, see description below), then A is
!>          unchanged, if double precision factorization has been used
!>          (INFO = 0 and ITER < 0, see description below), then the
!>          array A contains the factor U or L from the Cholesky
!>          factorization A = U**T*U or A = L*L**T.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in]B
!>          B is DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (LDX,NRHS)
!>          If INFO = 0, the N-by-NRHS solution matrix X.
!> 
[in]LDX
!>          LDX is INTEGER
!>          The leading dimension of the array X.  LDX >= max(1,N).
!> 
[out]WORK
!>          WORK is DOUBLE PRECISION array, dimension (N,NRHS)
!>          This array is used to hold the residual vectors.
!> 
[out]SWORK
!>          SWORK is REAL array, dimension (N*(N+NRHS))
!>          This array is used to use the single precision matrix and the
!>          right-hand sides or solutions in single precision.
!> 
[out]ITER
!>          ITER is INTEGER
!>          < 0: iterative refinement has failed, double precision
!>               factorization has been performed
!>               -1 : the routine fell back to full precision for
!>                    implementation- or machine-specific reasons
!>               -2 : narrowing the precision induced an overflow,
!>                    the routine fell back to full precision
!>               -3 : failure of SPOTRF
!>               -31: stop the iterative refinement after the 30th
!>                    iterations
!>          > 0: iterative refinement has been successfully used.
!>               Returns the number of iterations
!> 
[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 (DOUBLE
!>                PRECISION) 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 197 of file dsposv.f.

199*
200* -- LAPACK driver routine --
201* -- LAPACK is a software package provided by Univ. of Tennessee, --
202* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
203*
204* .. Scalar Arguments ..
205 CHARACTER UPLO
206 INTEGER INFO, ITER, LDA, LDB, LDX, N, NRHS
207* ..
208* .. Array Arguments ..
209 REAL SWORK( * )
210 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), WORK( N, * ),
211 $ X( LDX, * )
212* ..
213*
214* =====================================================================
215*
216* .. Parameters ..
217 LOGICAL DOITREF
218 parameter( doitref = .true. )
219*
220 INTEGER ITERMAX
221 parameter( itermax = 30 )
222*
223 DOUBLE PRECISION BWDMAX
224 parameter( bwdmax = 1.0e+00 )
225*
226 DOUBLE PRECISION NEGONE, ONE
227 parameter( negone = -1.0d+0, one = 1.0d+0 )
228*
229* .. Local Scalars ..
230 INTEGER I, IITER, PTSA, PTSX
231 DOUBLE PRECISION ANRM, CTE, EPS, RNRM, XNRM
232*
233* .. External Subroutines ..
234 EXTERNAL daxpy, dsymm, dlacpy, dlat2s, dlag2s, slag2d,
236* ..
237* .. External Functions ..
238 INTEGER IDAMAX
239 DOUBLE PRECISION DLAMCH, DLANSY
240 LOGICAL LSAME
241 EXTERNAL idamax, dlamch, dlansy, lsame
242* ..
243* .. Intrinsic Functions ..
244 INTRINSIC abs, dble, max, sqrt
245* ..
246* .. Executable Statements ..
247*
248 info = 0
249 iter = 0
250*
251* Test the input parameters.
252*
253 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
254 info = -1
255 ELSE IF( n.LT.0 ) THEN
256 info = -2
257 ELSE IF( nrhs.LT.0 ) THEN
258 info = -3
259 ELSE IF( lda.LT.max( 1, n ) ) THEN
260 info = -5
261 ELSE IF( ldb.LT.max( 1, n ) ) THEN
262 info = -7
263 ELSE IF( ldx.LT.max( 1, n ) ) THEN
264 info = -9
265 END IF
266 IF( info.NE.0 ) THEN
267 CALL xerbla( 'DSPOSV', -info )
268 RETURN
269 END IF
270*
271* Quick return if (N.EQ.0).
272*
273 IF( n.EQ.0 )
274 $ RETURN
275*
276* Skip single precision iterative refinement if a priori slower
277* than double precision factorization.
278*
279 IF( .NOT.doitref ) THEN
280 iter = -1
281 GO TO 40
282 END IF
283*
284* Compute some constants.
285*
286 anrm = dlansy( 'I', uplo, n, a, lda, work )
287 eps = dlamch( 'Epsilon' )
288 cte = anrm*eps*sqrt( dble( n ) )*bwdmax
289*
290* Set the indices PTSA, PTSX for referencing SA and SX in SWORK.
291*
292 ptsa = 1
293 ptsx = ptsa + n*n
294*
295* Convert B from double precision to single precision and store the
296* result in SX.
297*
298 CALL dlag2s( n, nrhs, b, ldb, swork( ptsx ), n, info )
299*
300 IF( info.NE.0 ) THEN
301 iter = -2
302 GO TO 40
303 END IF
304*
305* Convert A from double precision to single precision and store the
306* result in SA.
307*
308 CALL dlat2s( uplo, n, a, lda, swork( ptsa ), n, info )
309*
310 IF( info.NE.0 ) THEN
311 iter = -2
312 GO TO 40
313 END IF
314*
315* Compute the Cholesky factorization of SA.
316*
317 CALL spotrf( uplo, n, swork( ptsa ), n, info )
318*
319 IF( info.NE.0 ) THEN
320 iter = -3
321 GO TO 40
322 END IF
323*
324* Solve the system SA*SX = SB.
325*
326 CALL spotrs( uplo, n, nrhs, swork( ptsa ), n, swork( ptsx ), n,
327 $ info )
328*
329* Convert SX back to double precision
330*
331 CALL slag2d( n, nrhs, swork( ptsx ), n, x, ldx, info )
332*
333* Compute R = B - AX (R is WORK).
334*
335 CALL dlacpy( 'All', n, nrhs, b, ldb, work, n )
336*
337 CALL dsymm( 'Left', uplo, n, nrhs, negone, a, lda, x, ldx, one,
338 $ work, n )
339*
340* Check whether the NRHS normwise backward errors satisfy the
341* stopping criterion. If yes, set ITER=0 and return.
342*
343 DO i = 1, nrhs
344 xnrm = abs( x( idamax( n, x( 1, i ), 1 ), i ) )
345 rnrm = abs( work( idamax( n, work( 1, i ), 1 ), i ) )
346 IF( rnrm.GT.xnrm*cte )
347 $ GO TO 10
348 END DO
349*
350* If we are here, the NRHS normwise backward errors satisfy the
351* stopping criterion. We are good to exit.
352*
353 iter = 0
354 RETURN
355*
356 10 CONTINUE
357*
358 DO 30 iiter = 1, itermax
359*
360* Convert R (in WORK) from double precision to single precision
361* and store the result in SX.
362*
363 CALL dlag2s( n, nrhs, work, n, swork( ptsx ), n, info )
364*
365 IF( info.NE.0 ) THEN
366 iter = -2
367 GO TO 40
368 END IF
369*
370* Solve the system SA*SX = SR.
371*
372 CALL spotrs( uplo, n, nrhs, swork( ptsa ), n, swork( ptsx ), n,
373 $ info )
374*
375* Convert SX back to double precision and update the current
376* iterate.
377*
378 CALL slag2d( n, nrhs, swork( ptsx ), n, work, n, info )
379*
380 DO i = 1, nrhs
381 CALL daxpy( n, one, work( 1, i ), 1, x( 1, i ), 1 )
382 END DO
383*
384* Compute R = B - AX (R is WORK).
385*
386 CALL dlacpy( 'All', n, nrhs, b, ldb, work, n )
387*
388 CALL dsymm( 'L', uplo, n, nrhs, negone, a, lda, x, ldx, one,
389 $ work, n )
390*
391* Check whether the NRHS normwise backward errors satisfy the
392* stopping criterion. If yes, set ITER=IITER>0 and return.
393*
394 DO i = 1, nrhs
395 xnrm = abs( x( idamax( n, x( 1, i ), 1 ), i ) )
396 rnrm = abs( work( idamax( n, work( 1, i ), 1 ), i ) )
397 IF( rnrm.GT.xnrm*cte )
398 $ GO TO 20
399 END DO
400*
401* If we are here, the NRHS normwise backward errors satisfy the
402* stopping criterion, we are good to exit.
403*
404 iter = iiter
405*
406 RETURN
407*
408 20 CONTINUE
409*
410 30 CONTINUE
411*
412* If we are at this place of the code, this is because we have
413* performed ITER=ITERMAX iterations and never satisfied the
414* stopping criterion, set up the ITER flag accordingly and follow
415* up on double precision routine.
416*
417 iter = -itermax - 1
418*
419 40 CONTINUE
420*
421* Single-precision iterative refinement failed to converge to a
422* satisfactory solution, so we resort to double precision.
423*
424 CALL dpotrf( uplo, n, a, lda, info )
425*
426 IF( info.NE.0 )
427 $ RETURN
428*
429 CALL dlacpy( 'All', n, nrhs, b, ldb, x, ldx )
430 CALL dpotrs( uplo, n, nrhs, a, lda, x, ldx, info )
431*
432 RETURN
433*
434* End of DSPOSV
435*
subroutine slag2d(m, n, sa, ldsa, a, lda, info)
SLAG2D converts a single precision matrix to a double precision matrix.
Definition slag2d.f:104
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
subroutine dlat2s(uplo, n, a, lda, sa, ldsa, info)
DLAT2S converts a double-precision triangular matrix to a single-precision triangular matrix.
Definition dlat2s.f:111
subroutine dlag2s(m, n, a, lda, sa, ldsa, info)
DLAG2S converts a double precision matrix to a single precision matrix.
Definition dlag2s.f:108
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
subroutine dsymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
DSYMM
Definition dsymm.f:189
subroutine spotrs(uplo, n, nrhs, a, lda, b, ldb, info)
SPOTRS
Definition spotrs.f:110
subroutine spotrf(uplo, n, a, lda, info)
SPOTRF
Definition spotrf.f:107