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

Functions

subroutine zcposv (uplo, n, nrhs, a, lda, b, ldb, x, ldx, work, swork, rwork, iter, info)
  ZCPOSV computes the solution to system of linear equations A * X = B for PO matrices
subroutine zposv (uplo, n, nrhs, a, lda, b, ldb, info)
  ZPOSV computes the solution to system of linear equations A * X = B for PO matrices
subroutine zposvx (fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
  ZPOSVX computes the solution to system of linear equations A * X = B for PO matrices
subroutine zposvxx (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, rwork, info)
  ZPOSVXX computes the solution to system of linear equations A * X = B for PO matrices

Detailed Description

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

Function Documentation

◆ zcposv()

subroutine zcposv ( character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( ldx, * ) x,
integer ldx,
complex*16, dimension( n, * ) work,
complex, dimension( * ) swork,
double precision, dimension( * ) rwork,
integer iter,
integer info )

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

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

Purpose:
!>
!> ZCPOSV computes the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian positive definite matrix and X and B
!> are N-by-NRHS matrices.
!>
!> ZCPOSV first attempts to factorize the matrix in COMPLEX and use this
!> factorization within an iterative refinement procedure to produce a
!> solution with COMPLEX*16 normwise backward error quality (see below).
!> If the approach fails the method switches to a COMPLEX*16
!> factorization and solve.
!>
!> The iterative refinement is not going to be a winning strategy if
!> the ratio COMPLEX performance over COMPLEX*16 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 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.
!>
!>          Note that the imaginary parts of the diagonal
!>          elements need not be set and are assumed to be zero.
!>
!>          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**H*U or A = L*L**H.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[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, 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 COMPLEX*16 array, dimension (N,NRHS)
!>          This array is used to hold the residual vectors.
!> 
[out]SWORK
!>          SWORK is COMPLEX 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]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (N)
!> 
[out]ITER
!>          ITER is INTEGER
!>          < 0: iterative refinement has failed, COMPLEX*16
!>               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 CPOTRF
!>               -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
!>                (COMPLEX*16) 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 207 of file zcposv.f.

209*
210* -- LAPACK driver routine --
211* -- LAPACK is a software package provided by Univ. of Tennessee, --
212* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
213*
214* .. Scalar Arguments ..
215 CHARACTER UPLO
216 INTEGER INFO, ITER, LDA, LDB, LDX, N, NRHS
217* ..
218* .. Array Arguments ..
219 DOUBLE PRECISION RWORK( * )
220 COMPLEX SWORK( * )
221 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( N, * ),
222 $ X( LDX, * )
223* ..
224*
225* =====================================================================
226*
227* .. Parameters ..
228 LOGICAL DOITREF
229 parameter( doitref = .true. )
230*
231 INTEGER ITERMAX
232 parameter( itermax = 30 )
233*
234 DOUBLE PRECISION BWDMAX
235 parameter( bwdmax = 1.0e+00 )
236*
237 COMPLEX*16 NEGONE, ONE
238 parameter( negone = ( -1.0d+00, 0.0d+00 ),
239 $ one = ( 1.0d+00, 0.0d+00 ) )
240*
241* .. Local Scalars ..
242 INTEGER I, IITER, PTSA, PTSX
243 DOUBLE PRECISION ANRM, CTE, EPS, RNRM, XNRM
244 COMPLEX*16 ZDUM
245*
246* .. External Subroutines ..
247 EXTERNAL zaxpy, zhemm, zlacpy, zlat2c, zlag2c, clag2z,
249* ..
250* .. External Functions ..
251 INTEGER IZAMAX
252 DOUBLE PRECISION DLAMCH, ZLANHE
253 LOGICAL LSAME
254 EXTERNAL izamax, dlamch, zlanhe, lsame
255* ..
256* .. Intrinsic Functions ..
257 INTRINSIC abs, dble, max, sqrt
258* .. Statement Functions ..
259 DOUBLE PRECISION CABS1
260* ..
261* .. Statement Function definitions ..
262 cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
263* ..
264* .. Executable Statements ..
265*
266 info = 0
267 iter = 0
268*
269* Test the input parameters.
270*
271 IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
272 info = -1
273 ELSE IF( n.LT.0 ) THEN
274 info = -2
275 ELSE IF( nrhs.LT.0 ) THEN
276 info = -3
277 ELSE IF( lda.LT.max( 1, n ) ) THEN
278 info = -5
279 ELSE IF( ldb.LT.max( 1, n ) ) THEN
280 info = -7
281 ELSE IF( ldx.LT.max( 1, n ) ) THEN
282 info = -9
283 END IF
284 IF( info.NE.0 ) THEN
285 CALL xerbla( 'ZCPOSV', -info )
286 RETURN
287 END IF
288*
289* Quick return if (N.EQ.0).
290*
291 IF( n.EQ.0 )
292 $ RETURN
293*
294* Skip single precision iterative refinement if a priori slower
295* than double precision factorization.
296*
297 IF( .NOT.doitref ) THEN
298 iter = -1
299 GO TO 40
300 END IF
301*
302* Compute some constants.
303*
304 anrm = zlanhe( 'I', uplo, n, a, lda, rwork )
305 eps = dlamch( 'Epsilon' )
306 cte = anrm*eps*sqrt( dble( n ) )*bwdmax
307*
308* Set the indices PTSA, PTSX for referencing SA and SX in SWORK.
309*
310 ptsa = 1
311 ptsx = ptsa + n*n
312*
313* Convert B from double precision to single precision and store the
314* result in SX.
315*
316 CALL zlag2c( n, nrhs, b, ldb, swork( ptsx ), n, info )
317*
318 IF( info.NE.0 ) THEN
319 iter = -2
320 GO TO 40
321 END IF
322*
323* Convert A from double precision to single precision and store the
324* result in SA.
325*
326 CALL zlat2c( uplo, n, a, lda, swork( ptsa ), n, info )
327*
328 IF( info.NE.0 ) THEN
329 iter = -2
330 GO TO 40
331 END IF
332*
333* Compute the Cholesky factorization of SA.
334*
335 CALL cpotrf( uplo, n, swork( ptsa ), n, info )
336*
337 IF( info.NE.0 ) THEN
338 iter = -3
339 GO TO 40
340 END IF
341*
342* Solve the system SA*SX = SB.
343*
344 CALL cpotrs( uplo, n, nrhs, swork( ptsa ), n, swork( ptsx ), n,
345 $ info )
346*
347* Convert SX back to COMPLEX*16
348*
349 CALL clag2z( n, nrhs, swork( ptsx ), n, x, ldx, info )
350*
351* Compute R = B - AX (R is WORK).
352*
353 CALL zlacpy( 'All', n, nrhs, b, ldb, work, n )
354*
355 CALL zhemm( 'Left', uplo, n, nrhs, negone, a, lda, x, ldx, one,
356 $ work, n )
357*
358* Check whether the NRHS normwise backward errors satisfy the
359* stopping criterion. If yes, set ITER=0 and return.
360*
361 DO i = 1, nrhs
362 xnrm = cabs1( x( izamax( n, x( 1, i ), 1 ), i ) )
363 rnrm = cabs1( work( izamax( n, work( 1, i ), 1 ), i ) )
364 IF( rnrm.GT.xnrm*cte )
365 $ GO TO 10
366 END DO
367*
368* If we are here, the NRHS normwise backward errors satisfy the
369* stopping criterion. We are good to exit.
370*
371 iter = 0
372 RETURN
373*
374 10 CONTINUE
375*
376 DO 30 iiter = 1, itermax
377*
378* Convert R (in WORK) from double precision to single precision
379* and store the result in SX.
380*
381 CALL zlag2c( n, nrhs, work, n, swork( ptsx ), n, info )
382*
383 IF( info.NE.0 ) THEN
384 iter = -2
385 GO TO 40
386 END IF
387*
388* Solve the system SA*SX = SR.
389*
390 CALL cpotrs( uplo, n, nrhs, swork( ptsa ), n, swork( ptsx ), n,
391 $ info )
392*
393* Convert SX back to double precision and update the current
394* iterate.
395*
396 CALL clag2z( n, nrhs, swork( ptsx ), n, work, n, info )
397*
398 DO i = 1, nrhs
399 CALL zaxpy( n, one, work( 1, i ), 1, x( 1, i ), 1 )
400 END DO
401*
402* Compute R = B - AX (R is WORK).
403*
404 CALL zlacpy( 'All', n, nrhs, b, ldb, work, n )
405*
406 CALL zhemm( 'L', uplo, n, nrhs, negone, a, lda, x, ldx, one,
407 $ work, n )
408*
409* Check whether the NRHS normwise backward errors satisfy the
410* stopping criterion. If yes, set ITER=IITER>0 and return.
411*
412 DO i = 1, nrhs
413 xnrm = cabs1( x( izamax( n, x( 1, i ), 1 ), i ) )
414 rnrm = cabs1( work( izamax( n, work( 1, i ), 1 ), i ) )
415 IF( rnrm.GT.xnrm*cte )
416 $ GO TO 20
417 END DO
418*
419* If we are here, the NRHS normwise backward errors satisfy the
420* stopping criterion, we are good to exit.
421*
422 iter = iiter
423*
424 RETURN
425*
426 20 CONTINUE
427*
428 30 CONTINUE
429*
430* If we are at this place of the code, this is because we have
431* performed ITER=ITERMAX iterations and never satisfied the
432* stopping criterion, set up the ITER flag accordingly and follow
433* up on double precision routine.
434*
435 iter = -itermax - 1
436*
437 40 CONTINUE
438*
439* Single-precision iterative refinement failed to converge to a
440* satisfactory solution, so we resort to double precision.
441*
442 CALL zpotrf( uplo, n, a, lda, info )
443*
444 IF( info.NE.0 )
445 $ RETURN
446*
447 CALL zlacpy( 'All', n, nrhs, b, ldb, x, ldx )
448 CALL zpotrs( uplo, n, nrhs, a, lda, x, ldx, info )
449*
450 RETURN
451*
452* End of ZCPOSV
453*
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
integer function izamax(n, zx, incx)
IZAMAX
Definition izamax.f:71
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 zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine clag2z(m, n, sa, ldsa, a, lda, info)
CLAG2Z converts a complex single precision matrix to a complex double precision matrix.
Definition clag2z.f:103
subroutine zlat2c(uplo, n, a, lda, sa, ldsa, info)
ZLAT2C converts a double complex triangular matrix to a complex triangular matrix.
Definition zlat2c.f:111
subroutine zlag2c(m, n, a, lda, sa, ldsa, info)
ZLAG2C converts a complex double precision matrix to a complex single precision matrix.
Definition zlag2c.f:107
subroutine zpotrs(uplo, n, nrhs, a, lda, b, ldb, info)
ZPOTRS
Definition zpotrs.f:110
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
subroutine zhemm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
ZHEMM
Definition zhemm.f:191
subroutine cpotrs(uplo, n, nrhs, a, lda, b, ldb, info)
CPOTRS
Definition cpotrs.f:110
subroutine cpotrf(uplo, n, a, lda, info)
CPOTRF
Definition cpotrf.f:107
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
Definition zpotrf.f:102
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
#define max(a, b)
Definition macros.h:21

◆ zposv()

subroutine zposv ( character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldb, * ) b,
integer ldb,
integer info )

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

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

Purpose:
!>
!> ZPOSV computes the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian positive definite matrix and X and B
!> are N-by-NRHS matrices.
!>
!> The Cholesky decomposition is used to factor A as
!>    A = U**H* U,  if UPLO = 'U', or
!>    A = L * L**H,  if UPLO = 'L',
!> where U is an upper triangular matrix and  L is a lower triangular
!> matrix.  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 factor U or L from the Cholesky
!>          factorization A = U**H *U or A = L*L**H.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[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]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 zposv.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 COMPLEX*16 A( LDA, * ), B( LDB, * )
141* ..
142*
143* =====================================================================
144*
145* .. External Functions ..
146 LOGICAL LSAME
147 EXTERNAL lsame
148* ..
149* .. External Subroutines ..
150 EXTERNAL xerbla, zpotrf, zpotrs
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( 'ZPOSV ', -info )
173 RETURN
174 END IF
175*
176* Compute the Cholesky factorization A = U**H *U or A = L*L**H.
177*
178 CALL zpotrf( 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 zpotrs( uplo, n, nrhs, a, lda, b, ldb, info )
184*
185 END IF
186 RETURN
187*
188* End of ZPOSV
189*

◆ zposvx()

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

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

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

Purpose:
!>
!> ZPOSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to
!> compute the solution to a complex system of linear equations
!>    A * X = B,
!> where A is an N-by-N Hermitian positive definite matrix and X and B
!> are N-by-NRHS matrices.
!>
!> Error bounds on the solution and a condition estimate are also
!> provided.
!> 
Description:
!>
!> The following steps are performed:
!>
!> 1. If FACT = 'E', real scaling factors are computed to equilibrate
!>    the system:
!>       diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B
!>    Whether or not the system will be equilibrated depends on the
!>    scaling of the matrix A, but if equilibration is used, A is
!>    overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
!>
!> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
!>    factor the matrix A (after equilibration if FACT = 'E') as
!>       A = U**H* U,  if UPLO = 'U', or
!>       A = L * L**H,  if UPLO = 'L',
!>    where U is an upper triangular matrix and L is a lower triangular
!>    matrix.
!>
!> 3. If the leading 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 COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the Hermitian matrix A, except if FACT = 'F' and
!>          EQUED = 'Y', then A must contain the equilibrated matrix
!>          diag(S)*A*diag(S).  If UPLO = 'U', the leading
!>          N-by-N upper triangular part of A contains the upper
!>          triangular part of the matrix A, and the strictly lower
!>          triangular part of A is not referenced.  If UPLO = 'L', the
!>          leading N-by-N lower triangular part of A contains the lower
!>          triangular part of the matrix A, and the strictly upper
!>          triangular part of A is not referenced.  A is not modified if
!>          FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit.
!>
!>          On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
!>          diag(S)*A*diag(S).
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in,out]AF
!>          AF is COMPLEX*16 array, dimension (LDAF,N)
!>          If FACT = 'F', then AF is an input argument and on entry
!>          contains the triangular factor U or L from the Cholesky
!>          factorization A = U**H *U or A = L*L**H, in the same storage
!>          format as A.  If EQUED .ne. 'N', then AF is the factored form
!>          of the equilibrated matrix diag(S)*A*diag(S).
!>
!>          If FACT = 'N', then AF is an output argument and on exit
!>          returns the triangular factor U or L from the Cholesky
!>          factorization A = U**H *U or A = L*L**H of the original
!>          matrix A.
!>
!>          If FACT = 'E', then AF is an output argument and on exit
!>          returns the triangular factor U or L from the Cholesky
!>          factorization A = U**H *U or A = L*L**H of the equilibrated
!>          matrix A (see the description of A for the form of the
!>          equilibrated matrix).
!> 
[in]LDAF
!>          LDAF is INTEGER
!>          The leading dimension of the array AF.  LDAF >= max(1,N).
!> 
[in,out]EQUED
!>          EQUED is CHARACTER*1
!>          Specifies the form of equilibration that was done.
!>          = 'N':  No equilibration (always true if FACT = 'N').
!>          = 'Y':  Equilibration was done, i.e., A has been replaced by
!>                  diag(S) * A * diag(S).
!>          EQUED is an input argument if FACT = 'F'; otherwise, it is an
!>          output argument.
!> 
[in,out]S
!>          S is DOUBLE PRECISION array, dimension (N)
!>          The scale factors for A; not accessed if EQUED = 'N'.  S is
!>          an input argument if FACT = 'F'; otherwise, S is an output
!>          argument.  If FACT = 'F' and EQUED = 'Y', each element of S
!>          must be positive.
!> 
[in,out]B
!>          B is COMPLEX*16 array, dimension (LDB,NRHS)
!>          On entry, the N-by-NRHS righthand side matrix B.
!>          On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
!>          B is overwritten by diag(S) * B.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]X
!>          X is COMPLEX*16 array, dimension (LDX,NRHS)
!>          If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
!>          the original system of equations.  Note that if EQUED = 'Y',
!>          A and B are modified on exit, and the solution to the
!>          equilibrated system is inv(diag(S))*X.
!> 
[in]LDX
!>          LDX is INTEGER
!>          The leading dimension of the array X.  LDX >= max(1,N).
!> 
[out]RCOND
!>          RCOND is DOUBLE PRECISION
!>          The estimate of the reciprocal condition number of the matrix
!>          A after equilibration (if done).  If RCOND is less than the
!>          machine precision (in particular, if RCOND = 0), the matrix
!>          is singular to working precision.  This condition is
!>          indicated by a return code of INFO > 0.
!> 
[out]FERR
!>          FERR is DOUBLE PRECISION array, dimension (NRHS)
!>          The estimated forward error bound for each solution vector
!>          X(j) (the j-th column of the solution matrix X).
!>          If XTRUE is the true solution corresponding to X(j), FERR(j)
!>          is an estimated upper bound for the magnitude of the largest
!>          element in (X(j) - XTRUE) divided by the magnitude of the
!>          largest element in X(j).  The estimate is as reliable as
!>          the estimate for RCOND, and is almost always a slight
!>          overestimate of the true error.
!> 
[out]BERR
!>          BERR is DOUBLE PRECISION array, dimension (NRHS)
!>          The componentwise relative backward error of each solution
!>          vector X(j) (i.e., the smallest relative change in
!>          any element of A or B that makes X(j) an exact solution).
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension (2*N)
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (N)
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value
!>          > 0: if INFO = i, and i is
!>                <= N:  the leading 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 303 of file zposvx.f.

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

◆ zposvxx()

subroutine zposvxx ( character fact,
character uplo,
integer n,
integer nrhs,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldaf, * ) af,
integer ldaf,
character equed,
double precision, dimension( * ) s,
complex*16, dimension( ldb, * ) b,
integer ldb,
complex*16, dimension( ldx, * ) x,
integer ldx,
double precision rcond,
double precision 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 )

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

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

Purpose:
!>
!>    ZPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
!>    to compute the solution to a complex*16 system of linear equations
!>    A * X = B, where A is an N-by-N Hermitian positive definite matrix
!>    and X and B are N-by-NRHS matrices.
!>
!>    If requested, both normwise and maximum componentwise error bounds
!>    are returned. ZPOSVXX 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.
!>
!>    ZPOSVXX 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
!>    ZPOSVXX 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 ZPOSVXX 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 COMPLEX*16 array, dimension (LDA,N)
!>     On entry, the Hermitian matrix A, except if FACT = 'F' and EQUED =
!>     'Y', then A must contain the equilibrated matrix
!>     diag(S)*A*diag(S).  If UPLO = 'U', the leading N-by-N upper
!>     triangular part of A contains the upper triangular part of the
!>     matrix A, and the strictly lower triangular part of A is not
!>     referenced.  If UPLO = 'L', the leading N-by-N lower triangular
!>     part of A contains the lower triangular part of the matrix A, and
!>     the strictly upper triangular part of A is not referenced.  A is
!>     not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
!>     'N' on exit.
!>
!>     On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
!>     diag(S)*A*diag(S).
!> 
[in]LDA
!>          LDA is INTEGER
!>     The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in,out]AF
!>          AF is COMPLEX*16 array, dimension (LDAF,N)
!>     If FACT = 'F', then AF is an input argument and on entry
!>     contains the triangular factor U or L from the Cholesky
!>     factorization A = U**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 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 (2*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 489 of file zposvxx.f.

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