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

Functions

subroutine zhesv (uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
  ZHESV computes the solution to system of linear equations A * X = B for HE matrices
subroutine zhesv_aa (uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
  ZHESV_AA computes the solution to system of linear equations A * X = B for HE matrices
subroutine zhesv_aa_2stage (uplo, n, nrhs, a, lda, tb, ltb, ipiv, ipiv2, b, ldb, work, lwork, info)
  ZHESV_AA_2STAGE computes the solution to system of linear equations A * X = B for HE matrices
subroutine zhesv_rk (uplo, n, nrhs, a, lda, e, ipiv, b, ldb, work, lwork, info)
  ZHESV_RK computes the solution to system of linear equations A * X = B for SY matrices
subroutine zhesv_rook (uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
 ZHESV_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using the bounded Bunch-Kaufman ("rook") diagonal pivoting method
subroutine zhesvx (fact, uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, lwork, rwork, info)
  ZHESVX computes the solution to system of linear equations A * X = B for HE matrices
subroutine zhesvxx (fact, uplo, n, nrhs, a, lda, af, ldaf, ipiv, equed, s, b, ldb, x, ldx, rcond, rpvgrw, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, rwork, info)
  ZHESVXX computes the solution to system of linear equations A * X = B for HE matrices

Detailed Description

This is the group of complex16 solve driver functions for HE matrices

Function Documentation

◆ zhesv()

subroutine zhesv ( character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
integer, dimension( * ) ipiv,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( * ) work,
integer lwork,
integer info )

ZHESV computes the solution to system of linear equations A * X = B for HE matrices

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

Purpose:
!>
!> ZHESV computes the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian matrix and X and B are N-by-NRHS
!> matrices.
!>
!> The diagonal pivoting method is used to factor A as
!>    A = U * D * U**H,  if UPLO = 'U', or
!>    A = L * D * L**H,  if UPLO = 'L',
!> where U (or L) is a product of permutation and unit upper (lower)
!> triangular matrices, and D is Hermitian and block diagonal with
!> 1-by-1 and 2-by-2 diagonal blocks.  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 COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the Hermitian 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 block diagonal matrix D and the
!>          multipliers used to obtain the factor U or L from the
!>          factorization A = U*D*U**H or A = L*D*L**H as computed by
!>          ZHETRF.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          Details of the interchanges and the block structure of D, as
!>          determined by ZHETRF.  If IPIV(k) > 0, then rows and columns
!>          k and IPIV(k) were interchanged, and D(k,k) is a 1-by-1
!>          diagonal block.  If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0,
!>          then rows and columns k-1 and -IPIV(k) were interchanged and
!>          D(k-1:k,k-1:k) is a 2-by-2 diagonal block.  If UPLO = 'L' and
!>          IPIV(k) = IPIV(k+1) < 0, then rows and columns k+1 and
!>          -IPIV(k) were interchanged and D(k:k+1,k:k+1) is a 2-by-2
!>          diagonal block.
!> 
[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 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]WORK
!>          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of WORK.  LWORK >= 1, and for best performance
!>          LWORK >= max(1,N*NB), where NB is the optimal blocksize for
!>          ZHETRF.
!>          for LWORK < N, TRS will be done with Level BLAS 2
!>          for LWORK >= N, TRS will be done with Level BLAS 3
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!>          > 0: if INFO = i, D(i,i) is exactly zero.  The factorization
!>               has been completed, but the block diagonal matrix D is
!>               exactly singular, so the solution could not be computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 169 of file zhesv.f.

171*
172* -- LAPACK driver routine --
173* -- LAPACK is a software package provided by Univ. of Tennessee, --
174* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176* .. Scalar Arguments ..
177 CHARACTER UPLO
178 INTEGER INFO, LDA, LDB, LWORK, N, NRHS
179* ..
180* .. Array Arguments ..
181 INTEGER IPIV( * )
182 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
183* ..
184*
185* =====================================================================
186*
187* .. Local Scalars ..
188 LOGICAL LQUERY
189 INTEGER LWKOPT, NB
190* ..
191* .. External Functions ..
192 LOGICAL LSAME
193 INTEGER ILAENV
194 EXTERNAL lsame, ilaenv
195* ..
196* .. External Subroutines ..
197 EXTERNAL xerbla, zhetrf, zhetrs, zhetrs2
198* ..
199* .. Intrinsic Functions ..
200 INTRINSIC max
201* ..
202* .. Executable Statements ..
203*
204* Test the input parameters.
205*
206 info = 0
207 lquery = ( lwork.EQ.-1 )
208 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
209 info = -1
210 ELSE IF( n.LT.0 ) THEN
211 info = -2
212 ELSE IF( nrhs.LT.0 ) THEN
213 info = -3
214 ELSE IF( lda.LT.max( 1, n ) ) THEN
215 info = -5
216 ELSE IF( ldb.LT.max( 1, n ) ) THEN
217 info = -8
218 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
219 info = -10
220 END IF
221*
222 IF( info.EQ.0 ) THEN
223 IF( n.EQ.0 ) THEN
224 lwkopt = 1
225 ELSE
226 nb = ilaenv( 1, 'ZHETRF', uplo, n, -1, -1, -1 )
227 lwkopt = n*nb
228 END IF
229 work( 1 ) = lwkopt
230 END IF
231*
232 IF( info.NE.0 ) THEN
233 CALL xerbla( 'ZHESV ', -info )
234 RETURN
235 ELSE IF( lquery ) THEN
236 RETURN
237 END IF
238*
239* Compute the factorization A = U*D*U**H or A = L*D*L**H.
240*
241 CALL zhetrf( uplo, n, a, lda, ipiv, work, lwork, info )
242 IF( info.EQ.0 ) THEN
243*
244* Solve the system A*X = B, overwriting B with X.
245*
246 IF ( lwork.LT.n ) THEN
247*
248* Solve with TRS ( Use Level BLAS 2)
249*
250 CALL zhetrs( uplo, n, nrhs, a, lda, ipiv, b, ldb, info )
251*
252 ELSE
253*
254* Solve with TRS2 ( Use Level BLAS 3)
255*
256 CALL zhetrs2( uplo,n,nrhs,a,lda,ipiv,b,ldb,work,info )
257*
258 END IF
259*
260 END IF
261*
262 work( 1 ) = lwkopt
263*
264 RETURN
265*
266* End of ZHESV
267*
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
subroutine zhetrf(uplo, n, a, lda, ipiv, work, lwork, info)
ZHETRF
Definition zhetrf.f:177
subroutine zhetrs(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
ZHETRS
Definition zhetrs.f:120
subroutine zhetrs2(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, info)
ZHETRS2
Definition zhetrs2.f:127
#define max(a, b)
Definition macros.h:21

◆ zhesv_aa()

subroutine zhesv_aa ( character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
integer, dimension( * ) ipiv,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( * ) work,
integer lwork,
integer info )

ZHESV_AA computes the solution to system of linear equations A * X = B for HE matrices

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

Purpose:
!>
!> ZHESV_AA computes the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian matrix and X and B are N-by-NRHS
!> matrices.
!>
!> Aasen's algorithm is used to factor A as
!>    A = U**H * T * U,  if UPLO = 'U', or
!>    A = L * T * L**H,  if UPLO = 'L',
!> where U (or L) is a product of permutation and unit upper (lower)
!> triangular matrices, and T is Hermitian and tridiagonal. 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 COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the Hermitian 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 tridiagonal matrix T and the
!>          multipliers used to obtain the factor U or L from the
!>          factorization A = U**H*T*U or A = L*T*L**H as computed by
!>          ZHETRF_AA.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          On exit, it contains the details of the interchanges, i.e.,
!>          the row and column k of A were interchanged with the
!>          row and column IPIV(k).
!> 
[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 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]WORK
!>          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of WORK.  LWORK >= MAX(1,2*N,3*N-2), and for best 
!>          performance LWORK >= max(1,N*NB), where NB is the optimal
!>          blocksize for ZHETRF.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!>          > 0: if INFO = i, D(i,i) is exactly zero.  The factorization
!>               has been completed, but the block diagonal matrix D is
!>               exactly singular, so the solution could not be computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 160 of file zhesv_aa.f.

162*
163* -- LAPACK driver routine --
164* -- LAPACK is a software package provided by Univ. of Tennessee, --
165* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
166*
167* .. Scalar Arguments ..
168 CHARACTER UPLO
169 INTEGER INFO, LDA, LDB, LWORK, N, NRHS
170* ..
171* .. Array Arguments ..
172 INTEGER IPIV( * )
173 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
174* ..
175*
176* =====================================================================
177*
178* .. Local Scalars ..
179 LOGICAL LQUERY
180 INTEGER LWKOPT, LWKOPT_HETRF, LWKOPT_HETRS
181* ..
182* .. External Functions ..
183 LOGICAL LSAME
184 INTEGER ILAENV
185 EXTERNAL lsame, ilaenv
186* ..
187* .. External Subroutines ..
188 EXTERNAL xerbla, zhetrf_aa, zhetrs_aa
189* ..
190* .. Intrinsic Functions ..
191 INTRINSIC max
192* ..
193* .. Executable Statements ..
194*
195* Test the input parameters.
196*
197 info = 0
198 lquery = ( lwork.EQ.-1 )
199 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
200 info = -1
201 ELSE IF( n.LT.0 ) THEN
202 info = -2
203 ELSE IF( nrhs.LT.0 ) THEN
204 info = -3
205 ELSE IF( lda.LT.max( 1, n ) ) THEN
206 info = -5
207 ELSE IF( ldb.LT.max( 1, n ) ) THEN
208 info = -8
209 ELSE IF( lwork.LT.max(2*n, 3*n-2) .AND. .NOT.lquery ) THEN
210 info = -10
211 END IF
212*
213 IF( info.EQ.0 ) THEN
214 CALL zhetrf_aa( uplo, n, a, lda, ipiv, work, -1, info )
215 lwkopt_hetrf = int( work(1) )
216 CALL zhetrs_aa( uplo, n, nrhs, a, lda, ipiv, b, ldb, work,
217 $ -1, info )
218 lwkopt_hetrs = int( work(1) )
219 lwkopt = max( lwkopt_hetrf, lwkopt_hetrs )
220 work( 1 ) = lwkopt
221 END IF
222*
223 IF( info.NE.0 ) THEN
224 CALL xerbla( 'ZHESV_AA ', -info )
225 RETURN
226 ELSE IF( lquery ) THEN
227 RETURN
228 END IF
229*
230* Compute the factorization A = U**H*T*U or A = L*T*L**H.
231*
232 CALL zhetrf_aa( uplo, n, a, lda, ipiv, work, lwork, info )
233 IF( info.EQ.0 ) THEN
234*
235* Solve the system A*X = B, overwriting B with X.
236*
237 CALL zhetrs_aa( uplo, n, nrhs, a, lda, ipiv, b, ldb, work,
238 $ lwork, info )
239*
240 END IF
241*
242 work( 1 ) = lwkopt
243*
244 RETURN
245*
246* End of ZHESV_AA
247*
subroutine zhetrf_aa(uplo, n, a, lda, ipiv, work, lwork, info)
ZHETRF_AA
Definition zhetrf_aa.f:132
subroutine zhetrs_aa(uplo, n, nrhs, a, lda, ipiv, b, ldb, work, lwork, info)
ZHETRS_AA
Definition zhetrs_aa.f:132

◆ zhesv_aa_2stage()

subroutine zhesv_aa_2stage ( character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( * ) tb,
integer ltb,
integer, dimension( * ) ipiv,
integer, dimension( * ) ipiv2,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( * ) work,
integer lwork,
integer info )

ZHESV_AA_2STAGE computes the solution to system of linear equations A * X = B for HE matrices

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

Purpose:
!>
!> ZHESV_AA_2STAGE computes the solution to a complex system of 
!> linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian matrix and X and B are N-by-NRHS
!> matrices.
!>
!> Aasen's 2-stage algorithm is used to factor A as
!>    A = U**H * T * U,  if UPLO = 'U', or
!>    A = L * T * L**H,  if UPLO = 'L',
!> where U (or L) is a product of permutation and unit upper (lower)
!> triangular matrices, and T is Hermitian and band. The matrix T is
!> then LU-factored with partial pivoting. The factored form of A
!> is then used to solve the system of equations A * X = B.
!>
!> This is the blocked version of the algorithm, calling Level 3 BLAS.
!> 
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 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 COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the hermitian 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, L is stored below (or above) the subdiaonal blocks,
!>          when UPLO  is 'L' (or 'U').
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]TB
!>          TB is COMPLEX*16 array, dimension (LTB)
!>          On exit, details of the LU factorization of the band matrix.
!> 
[in]LTB
!>          LTB is INTEGER
!>          The size of the array TB. LTB >= 4*N, internally
!>          used to select NB such that LTB >= (3*NB+1)*N.
!>
!>          If LTB = -1, then a workspace query is assumed; the
!>          routine only calculates the optimal size of LTB, 
!>          returns this value as the first entry of TB, and
!>          no error message related to LTB is issued by XERBLA.
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          On exit, it contains the details of the interchanges, i.e.,
!>          the row and column k of A were interchanged with the
!>          row and column IPIV(k).
!> 
[out]IPIV2
!>          IPIV2 is INTEGER array, dimension (N)
!>          On exit, it contains the details of the interchanges, i.e.,
!>          the row and column k of T were interchanged with the
!>          row and column IPIV(k).
!> 
[in,out]B
!>          B is COMPLEX*16 array, dimension (LDB,NRHS)
!>          On entry, the right hand side matrix B.
!>          On exit, the solution matrix X.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]WORK
!>          WORK is COMPLEX*16 workspace of size LWORK
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The size of WORK. LWORK >= N, internally used to select NB
!>          such that LWORK >= N*NB.
!>
!>          If LWORK = -1, then a workspace query is assumed; the
!>          routine only calculates the optimal size of the WORK array,
!>          returns this value as the first entry of the WORK array, and
!>          no error message related to LWORK is issued by XERBLA.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value.
!>          > 0:  if INFO = i, band LU factorization failed on i-th column
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 184 of file zhesv_aa_2stage.f.

187*
188* -- LAPACK driver routine --
189* -- LAPACK is a software package provided by Univ. of Tennessee, --
190* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191*
192 IMPLICIT NONE
193*
194* .. Scalar Arguments ..
195 CHARACTER UPLO
196 INTEGER N, NRHS, LDA, LDB, LTB, LWORK, INFO
197* ..
198* .. Array Arguments ..
199 INTEGER IPIV( * ), IPIV2( * )
200 COMPLEX*16 A( LDA, * ), B( LDB, * ), TB( * ), WORK( * )
201* ..
202*
203* =====================================================================
204* .. Parameters ..
205 COMPLEX*16 ZERO, ONE
206 parameter( zero = ( 0.0d+0, 0.0d+0 ),
207 $ one = ( 1.0d+0, 0.0d+0 ) )
208*
209* .. Local Scalars ..
210 LOGICAL UPPER, TQUERY, WQUERY
211 INTEGER LWKOPT
212* ..
213* .. External Functions ..
214 LOGICAL LSAME
215 INTEGER ILAENV
216 EXTERNAL lsame, ilaenv
217* ..
218* .. External Subroutines ..
220* ..
221* .. Intrinsic Functions ..
222 INTRINSIC max
223* ..
224* .. Executable Statements ..
225*
226* Test the input parameters.
227*
228 info = 0
229 upper = lsame( uplo, 'U' )
230 wquery = ( lwork.EQ.-1 )
231 tquery = ( ltb.EQ.-1 )
232 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
233 info = -1
234 ELSE IF( n.LT.0 ) THEN
235 info = -2
236 ELSE IF( nrhs.LT.0 ) THEN
237 info = -3
238 ELSE IF( lda.LT.max( 1, n ) ) THEN
239 info = -5
240 ELSE IF( ltb.LT.( 4*n ) .AND. .NOT.tquery ) THEN
241 info = -7
242 ELSE IF( ldb.LT.max( 1, n ) ) THEN
243 info = -11
244 ELSE IF( lwork.LT.n .AND. .NOT.wquery ) THEN
245 info = -13
246 END IF
247*
248 IF( info.EQ.0 ) THEN
249 CALL zhetrf_aa_2stage( uplo, n, a, lda, tb, -1, ipiv,
250 $ ipiv2, work, -1, info )
251 lwkopt = int( work(1) )
252 END IF
253*
254 IF( info.NE.0 ) THEN
255 CALL xerbla( 'ZHESV_AA_2STAGE', -info )
256 RETURN
257 ELSE IF( wquery .OR. tquery ) THEN
258 RETURN
259 END IF
260*
261* Compute the factorization A = U**H*T*U or A = L*T*L**H.
262*
263 CALL zhetrf_aa_2stage( uplo, n, a, lda, tb, ltb, ipiv, ipiv2,
264 $ work, lwork, info )
265 IF( info.EQ.0 ) THEN
266*
267* Solve the system A*X = B, overwriting B with X.
268*
269 CALL zhetrs_aa_2stage( uplo, n, nrhs, a, lda, tb, ltb, ipiv,
270 $ ipiv2, b, ldb, info )
271*
272 END IF
273*
274 work( 1 ) = lwkopt
275*
276 RETURN
277*
278* End of ZHESV_AA_2STAGE
279*
subroutine zhetrs_aa_2stage(uplo, n, nrhs, a, lda, tb, ltb, ipiv, ipiv2, b, ldb, info)
ZHETRS_AA_2STAGE
subroutine zhetrf_aa_2stage(uplo, n, a, lda, tb, ltb, ipiv, ipiv2, work, lwork, info)
ZHETRF_AA_2STAGE

◆ zhesv_rk()

subroutine zhesv_rk ( character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( * ) e,
integer, dimension( * ) ipiv,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( * ) work,
integer lwork,
integer info )

ZHESV_RK computes the solution to system of linear equations A * X = B for SY matrices

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

Purpose:
!> ZHESV_RK computes the solution to a complex system of linear
!> equations A * X = B, where A is an N-by-N Hermitian matrix
!> and X and B are N-by-NRHS matrices.
!>
!> The bounded Bunch-Kaufman (rook) diagonal pivoting method is used
!> to factor A as
!>    A = P*U*D*(U**H)*(P**T),  if UPLO = 'U', or
!>    A = P*L*D*(L**H)*(P**T),  if UPLO = 'L',
!> where U (or L) is unit upper (or lower) triangular matrix,
!> U**H (or L**H) is the conjugate of U (or L), P is a permutation
!> matrix, P**T is the transpose of P, and D is Hermitian and block
!> diagonal with 1-by-1 and 2-by-2 diagonal blocks.
!>
!> ZHETRF_RK is called to compute the factorization of a complex
!> Hermitian matrix.  The factored form of A is then used to solve
!> the system of equations A * X = B by calling BLAS3 routine ZHETRS_3.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>          Specifies whether the upper or lower triangular part of the
!>          Hermitian matrix A is stored:
!>          = '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 COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the Hermitian 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, diagonal of the block diagonal
!>          matrix D and factors U or L  as computed by ZHETRF_RK:
!>            a) ONLY diagonal elements of the Hermitian block diagonal
!>               matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
!>               (superdiagonal (or subdiagonal) elements of D
!>                are stored on exit in array E), and
!>            b) If UPLO = 'U': factor U in the superdiagonal part of A.
!>               If UPLO = 'L': factor L in the subdiagonal part of A.
!>
!>          For more info see the description of ZHETRF_RK routine.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]E
!>          E is COMPLEX*16 array, dimension (N)
!>          On exit, contains the output computed by the factorization
!>          routine ZHETRF_RK, i.e. the superdiagonal (or subdiagonal)
!>          elements of the Hermitian block diagonal matrix D
!>          with 1-by-1 or 2-by-2 diagonal blocks, where
!>          If UPLO = 'U': E(i) = D(i-1,i), i=2:N, E(1) is set to 0;
!>          If UPLO = 'L': E(i) = D(i+1,i), i=1:N-1, E(N) is set to 0.
!>
!>          NOTE: For 1-by-1 diagonal block D(k), where
!>          1 <= k <= N, the element E(k) is set to 0 in both
!>          UPLO = 'U' or UPLO = 'L' cases.
!>
!>          For more info see the description of ZHETRF_RK routine.
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          Details of the interchanges and the block structure of D,
!>          as determined by ZHETRF_RK.
!>
!>          For more info see the description of ZHETRF_RK routine.
!> 
[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 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]WORK
!>          WORK is COMPLEX*16 array, dimension ( MAX(1,LWORK) ).
!>          Work array used in the factorization stage.
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of WORK.  LWORK >= 1. For best performance
!>          of factorization stage LWORK >= max(1,N*NB), where NB is
!>          the optimal blocksize for ZHETRF_RK.
!>
!>          If LWORK = -1, then a workspace query is assumed;
!>          the routine only calculates the optimal size of the WORK
!>          array for factorization stage, returns this value as
!>          the first entry of the WORK array, and no error message
!>          related to LWORK is issued by XERBLA.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>
!>          < 0: If INFO = -k, the k-th argument had an illegal value
!>
!>          > 0: If INFO = k, the matrix A is singular, because:
!>                 If UPLO = 'U': column k in the upper
!>                 triangular part of A contains all zeros.
!>                 If UPLO = 'L': column k in the lower
!>                 triangular part of A contains all zeros.
!>
!>               Therefore D(k,k) is exactly zero, and superdiagonal
!>               elements of column k of U (or subdiagonal elements of
!>               column k of L ) are all zeros. The factorization has
!>               been completed, but the block diagonal matrix D is
!>               exactly singular, and division by zero will occur if
!>               it is used to solve a system of equations.
!>
!>               NOTE: INFO only stores the first occurrence of
!>               a singularity, any subsequent occurrence of singularity
!>               is not stored in INFO even though the factorization
!>               always completes.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
!>
!>  December 2016,  Igor Kozachenko,
!>                  Computer Science Division,
!>                  University of California, Berkeley
!>
!>  September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
!>                  School of Mathematics,
!>                  University of Manchester
!>
!> 

Definition at line 226 of file zhesv_rk.f.

228*
229* -- LAPACK driver routine --
230* -- LAPACK is a software package provided by Univ. of Tennessee, --
231* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
232*
233* .. Scalar Arguments ..
234 CHARACTER UPLO
235 INTEGER INFO, LDA, LDB, LWORK, N, NRHS
236* ..
237* .. Array Arguments ..
238 INTEGER IPIV( * )
239 COMPLEX*16 A( LDA, * ), B( LDB, * ), E( * ), WORK( * )
240* ..
241*
242* =====================================================================
243*
244* .. Local Scalars ..
245 LOGICAL LQUERY
246 INTEGER LWKOPT
247* ..
248* .. External Functions ..
249 LOGICAL LSAME
250 EXTERNAL lsame
251* ..
252* .. External Subroutines ..
253 EXTERNAL xerbla, zhetrf_rk, zhetrs_3
254* ..
255* .. Intrinsic Functions ..
256 INTRINSIC max
257* ..
258* .. Executable Statements ..
259*
260* Test the input parameters.
261*
262 info = 0
263 lquery = ( lwork.EQ.-1 )
264 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
265 info = -1
266 ELSE IF( n.LT.0 ) THEN
267 info = -2
268 ELSE IF( nrhs.LT.0 ) THEN
269 info = -3
270 ELSE IF( lda.LT.max( 1, n ) ) THEN
271 info = -5
272 ELSE IF( ldb.LT.max( 1, n ) ) THEN
273 info = -9
274 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
275 info = -11
276 END IF
277*
278 IF( info.EQ.0 ) THEN
279 IF( n.EQ.0 ) THEN
280 lwkopt = 1
281 ELSE
282 CALL zhetrf_rk( uplo, n, a, lda, e, ipiv, work, -1, info )
283 lwkopt = dble( work(1) )
284 END IF
285 work( 1 ) = lwkopt
286 END IF
287*
288 IF( info.NE.0 ) THEN
289 CALL xerbla( 'ZHESV_RK ', -info )
290 RETURN
291 ELSE IF( lquery ) THEN
292 RETURN
293 END IF
294*
295* Compute the factorization A = P*U*D*(U**H)*(P**T) or
296* A = P*U*D*(U**H)*(P**T).
297*
298 CALL zhetrf_rk( uplo, n, a, lda, e, ipiv, work, lwork, info )
299*
300 IF( info.EQ.0 ) THEN
301*
302* Solve the system A*X = B with BLAS3 solver, overwriting B with X.
303*
304 CALL zhetrs_3( uplo, n, nrhs, a, lda, e, ipiv, b, ldb, info )
305*
306 END IF
307*
308 work( 1 ) = lwkopt
309*
310 RETURN
311*
312* End of ZHESV_RK
313*
subroutine zhetrf_rk(uplo, n, a, lda, e, ipiv, work, lwork, info)
ZHETRF_RK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bunch...
Definition zhetrf_rk.f:259
subroutine zhetrs_3(uplo, n, nrhs, a, lda, e, ipiv, b, ldb, info)
ZHETRS_3
Definition zhetrs_3.f:165

◆ zhesv_rook()

subroutine zhesv_rook ( character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
integer, dimension( * ) ipiv,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( * ) work,
integer lwork,
integer info )

ZHESV_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using the bounded Bunch-Kaufman ("rook") diagonal pivoting method

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

Purpose:
!>
!> ZHESV_ROOK computes the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian matrix and X and B are N-by-NRHS
!> matrices.
!>
!> The bounded Bunch-Kaufman () diagonal pivoting method is used
!> to factor A as
!>    A = U * D * U**T,  if UPLO = 'U', or
!>    A = L * D * L**T,  if UPLO = 'L',
!> where U (or L) is a product of permutation and unit upper (lower)
!> triangular matrices, and D is Hermitian and block diagonal with
!> 1-by-1 and 2-by-2 diagonal blocks.
!>
!> ZHETRF_ROOK is called to compute the factorization of a complex
!> Hermition matrix A using the bounded Bunch-Kaufman () diagonal
!> pivoting method.
!>
!> The factored form of A is then used to solve the system
!> of equations A * X = B by calling ZHETRS_ROOK (uses BLAS 2).
!> 
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 COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the Hermitian 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 block diagonal matrix D and the
!>          multipliers used to obtain the factor U or L from the
!>          factorization A = U*D*U**H or A = L*D*L**H as computed by
!>          ZHETRF_ROOK.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          Details of the interchanges and the block structure of D.
!>
!>          If UPLO = 'U':
!>             Only the last KB elements of IPIV are set.
!>
!>             If IPIV(k) > 0, then rows and columns k and IPIV(k) were
!>             interchanged and D(k,k) is a 1-by-1 diagonal block.
!>
!>             If IPIV(k) < 0 and IPIV(k-1) < 0, then rows and
!>             columns k and -IPIV(k) were interchanged and rows and
!>             columns k-1 and -IPIV(k-1) were inerchaged,
!>             D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
!>
!>          If UPLO = 'L':
!>             Only the first KB elements of IPIV are set.
!>
!>             If IPIV(k) > 0, then rows and columns k and IPIV(k)
!>             were interchanged and D(k,k) is a 1-by-1 diagonal block.
!>
!>             If IPIV(k) < 0 and IPIV(k+1) < 0, then rows and
!>             columns k and -IPIV(k) were interchanged and rows and
!>             columns k+1 and -IPIV(k+1) were inerchaged,
!>             D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
!> 
[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 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]WORK
!>          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of WORK.  LWORK >= 1, and for best performance
!>          LWORK >= max(1,N*NB), where NB is the optimal blocksize for
!>          ZHETRF_ROOK.
!>          for LWORK < N, TRS will be done with Level BLAS 2
!>          for LWORK >= N, TRS will be done with Level BLAS 3
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!>          > 0: if INFO = i, D(i,i) is exactly zero.  The factorization
!>               has been completed, but the block diagonal matrix D is
!>               exactly singular, so the solution could not be computed.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
!>
!>  November 2013,  Igor Kozachenko,
!>                  Computer Science Division,
!>                  University of California, Berkeley
!>
!>  September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
!>                  School of Mathematics,
!>                  University of Manchester
!>
!> 

Definition at line 203 of file zhesv_rook.f.

205*
206* -- LAPACK driver routine --
207* -- LAPACK is a software package provided by Univ. of Tennessee, --
208* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
209*
210* .. Scalar Arguments ..
211 CHARACTER UPLO
212 INTEGER INFO, LDA, LDB, LWORK, N, NRHS
213* ..
214* .. Array Arguments ..
215 INTEGER IPIV( * )
216 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
217* ..
218*
219* =====================================================================
220*
221* .. Local Scalars ..
222 LOGICAL LQUERY
223 INTEGER LWKOPT, NB
224* ..
225* .. External Functions ..
226 LOGICAL LSAME
227 INTEGER ILAENV
228 EXTERNAL lsame, ilaenv
229* ..
230* .. External Subroutines ..
232* ..
233* .. Intrinsic Functions ..
234 INTRINSIC max
235* ..
236* .. Executable Statements ..
237*
238* Test the input parameters.
239*
240 info = 0
241 lquery = ( lwork.EQ.-1 )
242 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
243 info = -1
244 ELSE IF( n.LT.0 ) THEN
245 info = -2
246 ELSE IF( nrhs.LT.0 ) THEN
247 info = -3
248 ELSE IF( lda.LT.max( 1, n ) ) THEN
249 info = -5
250 ELSE IF( ldb.LT.max( 1, n ) ) THEN
251 info = -8
252 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
253 info = -10
254 END IF
255*
256 IF( info.EQ.0 ) THEN
257 IF( n.EQ.0 ) THEN
258 lwkopt = 1
259 ELSE
260 nb = ilaenv( 1, 'ZHETRF_ROOK', uplo, n, -1, -1, -1 )
261 lwkopt = n*nb
262 END IF
263 work( 1 ) = lwkopt
264 END IF
265*
266 IF( info.NE.0 ) THEN
267 CALL xerbla( 'ZHESV_ROOK ', -info )
268 RETURN
269 ELSE IF( lquery ) THEN
270 RETURN
271 END IF
272*
273* Compute the factorization A = U*D*U**H or A = L*D*L**H.
274*
275 CALL zhetrf_rook( uplo, n, a, lda, ipiv, work, lwork, info )
276 IF( info.EQ.0 ) THEN
277*
278* Solve the system A*X = B, overwriting B with X.
279*
280* Solve with TRS ( Use Level BLAS 2)
281*
282 CALL zhetrs_rook( uplo, n, nrhs, a, lda, ipiv, b, ldb, info )
283*
284 END IF
285*
286 work( 1 ) = lwkopt
287*
288 RETURN
289*
290* End of ZHESV_ROOK
291*
subroutine zhetrs_rook(uplo, n, nrhs, a, lda, ipiv, b, ldb, info)
ZHETRS_ROOK computes the solution to a system of linear equations A * X = B for HE matrices using fac...
subroutine zhetrf_rook(uplo, n, a, lda, ipiv, work, lwork, info)
ZHETRF_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bun...

◆ zhesvx()

subroutine zhesvx ( character fact,
character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldaf, * ) af,
integer ldaf,
integer, dimension( * ) ipiv,
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,
integer lwork,
double precision, dimension( * ) rwork,
integer info )

ZHESVX computes the solution to system of linear equations A * X = B for HE matrices

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

Purpose:
!>
!> ZHESVX uses the diagonal pivoting factorization to compute the
!> solution to a complex system of linear equations A * X = B,
!> where A is an N-by-N Hermitian 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 = 'N', the diagonal pivoting method is used to factor A.
!>    The form of the factorization is
!>       A = U * D * U**H,  if UPLO = 'U', or
!>       A = L * D * L**H,  if UPLO = 'L',
!>    where U (or L) is a product of permutation and unit upper (lower)
!>    triangular matrices, and D is Hermitian and block diagonal with
!>    1-by-1 and 2-by-2 diagonal blocks.
!>
!> 2. If some D(i,i)=0, so that D is exactly singular, then the routine
!>    returns with INFO = i. Otherwise, the factored form of A is used
!>    to estimate the condition number of the matrix A.  If the
!>    reciprocal of the condition number is less than machine precision,
!>    INFO = N+1 is returned as a warning, but the routine still goes on
!>    to solve for X and compute error bounds as described below.
!>
!> 3. The system of equations is solved for X using the factored form
!>    of A.
!>
!> 4. Iterative refinement is applied to improve the computed solution
!>    matrix and calculate error bounds and backward error estimates
!>    for it.
!> 
Parameters
[in]FACT
!>          FACT is CHARACTER*1
!>          Specifies whether or not the factored form of A has been
!>          supplied on entry.
!>          = 'F':  On entry, AF and IPIV contain the factored form
!>                  of A.  A, AF and IPIV will not be modified.
!>          = 'N':  The matrix A will be 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]A
!>          A is COMPLEX*16 array, dimension (LDA,N)
!>          The Hermitian 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.
!> 
[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 block diagonal matrix D and the multipliers used
!>          to obtain the factor U or L from the factorization
!>          A = U*D*U**H or A = L*D*L**H as computed by ZHETRF.
!>
!>          If FACT = 'N', then AF is an output argument and on exit
!>          returns the block diagonal matrix D and the multipliers used
!>          to obtain the factor U or L from the factorization
!>          A = U*D*U**H or A = L*D*L**H.
!> 
[in]LDAF
!>          LDAF is INTEGER
!>          The leading dimension of the array AF.  LDAF >= max(1,N).
!> 
[in,out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          If FACT = 'F', then IPIV is an input argument and on entry
!>          contains details of the interchanges and the block structure
!>          of D, as determined by ZHETRF.
!>          If IPIV(k) > 0, then rows and columns k and IPIV(k) were
!>          interchanged and D(k,k) is a 1-by-1 diagonal block.
!>          If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
!>          columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
!>          is a 2-by-2 diagonal block.  If UPLO = 'L' and IPIV(k) =
!>          IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
!>          interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
!>
!>          If FACT = 'N', then IPIV is an output argument and on exit
!>          contains details of the interchanges and the block structure
!>          of D, as determined by ZHETRF.
!> 
[in]B
!>          B is COMPLEX*16 array, dimension (LDB,NRHS)
!>          The N-by-NRHS right hand side matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]X
!>          X is COMPLEX*16 array, dimension (LDX,NRHS)
!>          If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X.
!> 
[in]LDX
!>          LDX is INTEGER
!>          The leading dimension of the array X.  LDX >= max(1,N).
!> 
[out]RCOND
!>          RCOND is DOUBLE PRECISION
!>          The estimate of the reciprocal condition number of the matrix
!>          A.  If RCOND is less than the machine precision (in
!>          particular, if RCOND = 0), the matrix is singular to working
!>          precision.  This condition is indicated by a return code of
!>          INFO > 0.
!> 
[out]FERR
!>          FERR is 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 (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of WORK.  LWORK >= max(1,2*N), and for best
!>          performance, when FACT = 'N', LWORK >= max(1,2*N,N*NB), where
!>          NB is the optimal blocksize for ZHETRF.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[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:  D(i,i) is exactly zero.  The factorization
!>                       has been completed but the factor D is exactly
!>                       singular, so the solution and error bounds could
!>                       not be computed. RCOND = 0 is returned.
!>                = N+1: D is nonsingular, but RCOND is less than machine
!>                       precision, meaning that the matrix is singular
!>                       to working precision.  Nevertheless, the
!>                       solution and error bounds are computed because
!>                       there are a number of situations where the
!>                       computed solution can be more accurate than the
!>                       value of RCOND would suggest.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 282 of file zhesvx.f.

285*
286* -- LAPACK driver routine --
287* -- LAPACK is a software package provided by Univ. of Tennessee, --
288* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
289*
290* .. Scalar Arguments ..
291 CHARACTER FACT, UPLO
292 INTEGER INFO, LDA, LDAF, LDB, LDX, LWORK, N, NRHS
293 DOUBLE PRECISION RCOND
294* ..
295* .. Array Arguments ..
296 INTEGER IPIV( * )
297 DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * )
298 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
299 $ WORK( * ), X( LDX, * )
300* ..
301*
302* =====================================================================
303*
304* .. Parameters ..
305 DOUBLE PRECISION ZERO
306 parameter( zero = 0.0d+0 )
307* ..
308* .. Local Scalars ..
309 LOGICAL LQUERY, NOFACT
310 INTEGER LWKOPT, NB
311 DOUBLE PRECISION ANORM
312* ..
313* .. External Functions ..
314 LOGICAL LSAME
315 INTEGER ILAENV
316 DOUBLE PRECISION DLAMCH, ZLANHE
317 EXTERNAL lsame, ilaenv, dlamch, zlanhe
318* ..
319* .. External Subroutines ..
320 EXTERNAL xerbla, zhecon, zherfs, zhetrf, zhetrs, zlacpy
321* ..
322* .. Intrinsic Functions ..
323 INTRINSIC max
324* ..
325* .. Executable Statements ..
326*
327* Test the input parameters.
328*
329 info = 0
330 nofact = lsame( fact, 'N' )
331 lquery = ( lwork.EQ.-1 )
332 IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
333 info = -1
334 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) )
335 $ THEN
336 info = -2
337 ELSE IF( n.LT.0 ) THEN
338 info = -3
339 ELSE IF( nrhs.LT.0 ) THEN
340 info = -4
341 ELSE IF( lda.LT.max( 1, n ) ) THEN
342 info = -6
343 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
344 info = -8
345 ELSE IF( ldb.LT.max( 1, n ) ) THEN
346 info = -11
347 ELSE IF( ldx.LT.max( 1, n ) ) THEN
348 info = -13
349 ELSE IF( lwork.LT.max( 1, 2*n ) .AND. .NOT.lquery ) THEN
350 info = -18
351 END IF
352*
353 IF( info.EQ.0 ) THEN
354 lwkopt = max( 1, 2*n )
355 IF( nofact ) THEN
356 nb = ilaenv( 1, 'ZHETRF', uplo, n, -1, -1, -1 )
357 lwkopt = max( lwkopt, n*nb )
358 END IF
359 work( 1 ) = lwkopt
360 END IF
361*
362 IF( info.NE.0 ) THEN
363 CALL xerbla( 'ZHESVX', -info )
364 RETURN
365 ELSE IF( lquery ) THEN
366 RETURN
367 END IF
368*
369 IF( nofact ) THEN
370*
371* Compute the factorization A = U*D*U**H or A = L*D*L**H.
372*
373 CALL zlacpy( uplo, n, n, a, lda, af, ldaf )
374 CALL zhetrf( uplo, n, af, ldaf, ipiv, work, lwork, info )
375*
376* Return if INFO is non-zero.
377*
378 IF( info.GT.0 )THEN
379 rcond = zero
380 RETURN
381 END IF
382 END IF
383*
384* Compute the norm of the matrix A.
385*
386 anorm = zlanhe( 'I', uplo, n, a, lda, rwork )
387*
388* Compute the reciprocal of the condition number of A.
389*
390 CALL zhecon( uplo, n, af, ldaf, ipiv, anorm, rcond, work, info )
391*
392* Compute the solution vectors X.
393*
394 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
395 CALL zhetrs( uplo, n, nrhs, af, ldaf, ipiv, x, ldx, info )
396*
397* Use iterative refinement to improve the computed solutions and
398* compute error bounds and backward error estimates for them.
399*
400 CALL zherfs( uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x,
401 $ ldx, ferr, berr, work, rwork, info )
402*
403* Set INFO = N+1 if the matrix is singular to working precision.
404*
405 IF( rcond.LT.dlamch( 'Epsilon' ) )
406 $ info = n + 1
407*
408 work( 1 ) = lwkopt
409*
410 RETURN
411*
412* End of ZHESVX
413*
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 zhecon(uplo, n, a, lda, ipiv, anorm, rcond, work, info)
ZHECON
Definition zhecon.f:125
subroutine zherfs(uplo, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZHERFS
Definition zherfs.f:192
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

◆ zhesvxx()

subroutine zhesvxx ( character fact,
character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldaf, * ) af,
integer ldaf,
integer, dimension( * ) ipiv,
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 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,
complex*16, dimension( * ) work,
double precision, dimension( * ) rwork,
integer info )

ZHESVXX computes the solution to system of linear equations A * X = B for HE matrices

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

Purpose:
!>
!>    ZHESVXX uses the diagonal pivoting factorization to compute the
!>    solution to a complex*16 system of linear equations A * X = B, where
!>    A is an N-by-N Hermitian matrix and X and B are N-by-NRHS
!>    matrices.
!>
!>    If requested, both normwise and maximum componentwise error bounds
!>    are returned. ZHESVXX 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.
!>
!>    ZHESVXX 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
!>    ZHESVXX 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 ZHESVXX 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 LU decomposition is used to factor
!>    the matrix A (after equilibration if FACT = 'E') as
!>
!>       A = U * D * U**T,  if UPLO = 'U', or
!>       A = L * D * L**T,  if UPLO = 'L',
!>
!>    where U (or L) is a product of permutation and unit upper (lower)
!>    triangular matrices, and D is Hermitian and block diagonal with
!>    1-by-1 and 2-by-2 diagonal blocks.
!>
!>    3. If some D(i,i)=0, so that D is exactly singular, then the
!>    routine returns with INFO = i. Otherwise, the factored form of A
!>    is used to estimate the condition number of the matrix A (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(R) 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 and IPIV contain the factored form of A.
!>               If EQUED is not 'N', the matrix A has been
!>               equilibrated with scaling factors given by S.
!>               A, AF, and IPIV 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 COMPLEX*16 array, dimension (LDA,N)
!>     The Hermitian 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 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 block diagonal matrix D and the multipliers
!>     used to obtain the factor U or L from the factorization A =
!>     U*D*U**H or A = L*D*L**H as computed by ZHETRF.
!>
!>     If FACT = 'N', then AF is an output argument and on exit
!>     returns the block diagonal matrix D and the multipliers
!>     used to obtain the factor U or L from the factorization A =
!>     U*D*U**H or A = L*D*L**H.
!> 
[in]LDAF
!>          LDAF is INTEGER
!>     The leading dimension of the array AF.  LDAF >= max(1,N).
!> 
[in,out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>     If FACT = 'F', then IPIV is an input argument and on entry
!>     contains details of the interchanges and the block
!>     structure of D, as determined by ZHETRF.  If IPIV(k) > 0,
!>     then rows and columns k and IPIV(k) were interchanged and
!>     D(k,k) is a 1-by-1 diagonal block.  If UPLO = 'U' and
!>     IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and
!>     -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2
!>     diagonal block.  If UPLO = 'L' and IPIV(k) = IPIV(k+1) < 0,
!>     then rows and columns k+1 and -IPIV(k) were interchanged
!>     and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
!>
!>     If FACT = 'N', then IPIV is an output argument and on exit
!>     contains details of the interchanges and the block
!>     structure of D, as determined by ZHETRF.
!> 
[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 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 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, 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 COMPLEX*16 array, dimension (5*N)
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (2*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 502 of file zhesvxx.f.

506*
507* -- LAPACK driver routine --
508* -- LAPACK is a software package provided by Univ. of Tennessee, --
509* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
510*
511* .. Scalar Arguments ..
512 CHARACTER EQUED, FACT, UPLO
513 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
514 $ N_ERR_BNDS
515 DOUBLE PRECISION RCOND, RPVGRW
516* ..
517* .. Array Arguments ..
518 INTEGER IPIV( * )
519 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
520 $ WORK( * ), X( LDX, * )
521 DOUBLE PRECISION S( * ), PARAMS( * ), BERR( * ), RWORK( * ),
522 $ ERR_BNDS_NORM( NRHS, * ),
523 $ ERR_BNDS_COMP( NRHS, * )
524* ..
525*
526* ==================================================================
527*
528* .. Parameters ..
529 DOUBLE PRECISION ZERO, ONE
530 parameter( zero = 0.0d+0, one = 1.0d+0 )
531 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
532 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
533 INTEGER CMP_ERR_I, PIV_GROWTH_I
534 parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
535 $ berr_i = 3 )
536 parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
537 parameter( cmp_rcond_i = 7, cmp_err_i = 8,
538 $ piv_growth_i = 9 )
539* ..
540* .. Local Scalars ..
541 LOGICAL EQUIL, NOFACT, RCEQU
542 INTEGER INFEQU, J
543 DOUBLE PRECISION AMAX, BIGNUM, SMIN, SMAX, SCOND, SMLNUM
544* ..
545* .. External Functions ..
546 EXTERNAL lsame, dlamch, zla_herpvgrw
547 LOGICAL LSAME
548 DOUBLE PRECISION DLAMCH, ZLA_HERPVGRW
549* ..
550* .. External Subroutines ..
551 EXTERNAL zheequb, zhetrf, zhetrs, zlacpy,
553* ..
554* .. Intrinsic Functions ..
555 INTRINSIC max, min
556* ..
557* .. Executable Statements ..
558*
559 info = 0
560 nofact = lsame( fact, 'N' )
561 equil = lsame( fact, 'E' )
562 smlnum = dlamch( 'Safe minimum' )
563 bignum = one / smlnum
564 IF( nofact .OR. equil ) THEN
565 equed = 'N'
566 rcequ = .false.
567 ELSE
568 rcequ = lsame( equed, 'Y' )
569 ENDIF
570*
571* Default is failure. If an input parameter is wrong or
572* factorization fails, make everything look horrible. Only the
573* pivot growth is set here, the rest is initialized in ZHERFSX.
574*
575 rpvgrw = zero
576*
577* Test the input parameters. PARAMS is not tested until ZHERFSX.
578*
579 IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.
580 $ lsame( fact, 'F' ) ) THEN
581 info = -1
582 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
583 $ .NOT.lsame( uplo, 'L' ) ) THEN
584 info = -2
585 ELSE IF( n.LT.0 ) THEN
586 info = -3
587 ELSE IF( nrhs.LT.0 ) THEN
588 info = -4
589 ELSE IF( lda.LT.max( 1, n ) ) THEN
590 info = -6
591 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
592 info = -8
593 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
594 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
595 info = -9
596 ELSE
597 IF ( rcequ ) THEN
598 smin = bignum
599 smax = zero
600 DO 10 j = 1, n
601 smin = min( smin, s( j ) )
602 smax = max( smax, s( j ) )
603 10 CONTINUE
604 IF( smin.LE.zero ) THEN
605 info = -10
606 ELSE IF( n.GT.0 ) THEN
607 scond = max( smin, smlnum ) / min( smax, bignum )
608 ELSE
609 scond = one
610 END IF
611 END IF
612 IF( info.EQ.0 ) THEN
613 IF( ldb.LT.max( 1, n ) ) THEN
614 info = -12
615 ELSE IF( ldx.LT.max( 1, n ) ) THEN
616 info = -14
617 END IF
618 END IF
619 END IF
620*
621 IF( info.NE.0 ) THEN
622 CALL xerbla( 'ZHESVXX', -info )
623 RETURN
624 END IF
625*
626 IF( equil ) THEN
627*
628* Compute row and column scalings to equilibrate the matrix A.
629*
630 CALL zheequb( uplo, n, a, lda, s, scond, amax, work, infequ )
631 IF( infequ.EQ.0 ) THEN
632*
633* Equilibrate the matrix.
634*
635 CALL zlaqhe( uplo, n, a, lda, s, scond, amax, equed )
636 rcequ = lsame( equed, 'Y' )
637 END IF
638 END IF
639*
640* Scale the right-hand side.
641*
642 IF( rcequ ) CALL zlascl2( n, nrhs, s, b, ldb )
643*
644 IF( nofact .OR. equil ) THEN
645*
646* Compute the LDL^H or UDU^H factorization of A.
647*
648 CALL zlacpy( uplo, n, n, a, lda, af, ldaf )
649 CALL zhetrf( uplo, n, af, ldaf, ipiv, work, 5*max(1,n), info )
650*
651* Return if INFO is non-zero.
652*
653 IF( info.GT.0 ) THEN
654*
655* Pivot in column INFO is exactly 0
656* Compute the reciprocal pivot growth factor of the
657* leading rank-deficient INFO columns of A.
658*
659 IF( n.GT.0 )
660 $ rpvgrw = zla_herpvgrw( uplo, n, info, a, lda, af, ldaf,
661 $ ipiv, rwork )
662 RETURN
663 END IF
664 END IF
665*
666* Compute the reciprocal pivot growth factor RPVGRW.
667*
668 IF( n.GT.0 )
669 $ rpvgrw = zla_herpvgrw( uplo, n, info, a, lda, af, ldaf, ipiv,
670 $ rwork )
671*
672* Compute the solution matrix X.
673*
674 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
675 CALL zhetrs( uplo, n, nrhs, af, ldaf, ipiv, x, ldx, info )
676*
677* Use iterative refinement to improve the computed solution and
678* compute error bounds and backward error estimates for it.
679*
680 CALL zherfsx( uplo, equed, n, nrhs, a, lda, af, ldaf, ipiv,
681 $ s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm,
682 $ err_bnds_comp, nparams, params, work, rwork, info )
683*
684* Scale solutions.
685*
686 IF ( rcequ ) THEN
687 CALL zlascl2 ( n, nrhs, s, x, ldx )
688 END IF
689*
690 RETURN
691*
692* End of ZHESVXX
693*
subroutine zlaqhe(uplo, n, a, lda, s, scond, amax, equed)
ZLAQHE scales a Hermitian matrix.
Definition zlaqhe.f:134
double precision function zla_herpvgrw(uplo, n, info, a, lda, af, ldaf, ipiv, work)
ZLA_HERPVGRW
subroutine zheequb(uplo, n, a, lda, s, scond, amax, work, info)
ZHEEQUB
Definition zheequb.f:132
subroutine zherfsx(uplo, equed, n, nrhs, a, lda, af, ldaf, ipiv, s, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, rwork, info)
ZHERFSX
Definition zherfsx.f:401
subroutine zlascl2(m, n, d, x, ldx)
ZLASCL2 performs diagonal scaling on a vector.
Definition zlascl2.f:91
#define min(a, b)
Definition macros.h:20