LAPACK  3.6.0
LAPACK: Linear Algebra PACKage
Collaboration diagram for complex:

Functions

real function cla_porcond_c (UPLO, N, A, LDA, AF, LDAF, C, CAPPLY, INFO, WORK, RWORK)
 CLA_PORCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for Hermitian positive-definite matrices. More...
 
real function cla_porcond_x (UPLO, N, A, LDA, AF, LDAF, X, INFO, WORK, RWORK)
 CLA_PORCOND_X computes the infinity norm condition number of op(A)*diag(x) for Hermitian positive-definite matrices. More...
 
subroutine cla_porfsx_extended (PREC_TYPE, UPLO, N, NRHS, A, LDA, AF, LDAF, COLEQU, C, B, LDB, Y, LDY, BERR_OUT, N_NORMS, ERR_BNDS_NORM, ERR_BNDS_COMP, RES, AYB, DY, Y_TAIL, RCOND, ITHRESH, RTHRESH, DZ_UB, IGNORE_CWISE, INFO)
 CLA_PORFSX_EXTENDED improves the computed solution to a system of linear equations for symmetric or Hermitian positive-definite matrices by performing extra-precise iterative refinement and provides error bounds and backward error estimates for the solution. More...
 
real function cla_porpvgrw (UPLO, NCOLS, A, LDA, AF, LDAF, WORK)
 CLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian positive-definite matrix. More...
 
subroutine cpocon (UPLO, N, A, LDA, ANORM, RCOND, WORK, RWORK, INFO)
 CPOCON More...
 
subroutine cpoequ (N, A, LDA, S, SCOND, AMAX, INFO)
 CPOEQU More...
 
subroutine cpoequb (N, A, LDA, S, SCOND, AMAX, INFO)
 CPOEQUB More...
 
subroutine cporfs (UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
 CPORFS More...
 
subroutine cporfsx (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)
 CPORFSX More...
 
subroutine cpotf2 (UPLO, N, A, LDA, INFO)
 CPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblocked algorithm). More...
 
subroutine cpotrf (UPLO, N, A, LDA, INFO)
 CPOTRF More...
 
recursive subroutine cpotrf2 (UPLO, N, A, LDA, INFO)
 CPOTRF2 More...
 
subroutine cpotri (UPLO, N, A, LDA, INFO)
 CPOTRI More...
 
subroutine cpotrs (UPLO, N, NRHS, A, LDA, B, LDB, INFO)
 CPOTRS More...
 

Detailed Description

This is the group of complex computational functions for PO matrices

Function Documentation

real function cla_porcond_c ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldaf, * )  AF,
integer  LDAF,
real, dimension( * )  C,
logical  CAPPLY,
integer  INFO,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK 
)

CLA_PORCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for Hermitian positive-definite matrices.

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

Purpose:
    CLA_PORCOND_C Computes the infinity norm condition number of
    op(A) * inv(diag(C)) where C is a DOUBLE PRECISION vector
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]A
          A is COMPLEX array, dimension (LDA,N)
     On entry, the N-by-N matrix A
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in]AF
          AF is COMPLEX array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**H*U or A = L*L**H, as computed by CPOTRF.
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in]C
          C is REAL array, dimension (N)
     The vector C in the formula op(A) * inv(diag(C)).
[in]CAPPLY
          CAPPLY is LOGICAL
     If .TRUE. then access the vector C in the formula above.
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
     i > 0:  The ith argument is invalid.
[in]WORK
          WORK is COMPLEX array, dimension (2*N).
     Workspace.
[in]RWORK
          RWORK is REAL array, dimension (N).
     Workspace.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 132 of file cla_porcond_c.f.

132 *
133 * -- LAPACK computational routine (version 3.4.2) --
134 * -- LAPACK is a software package provided by Univ. of Tennessee, --
135 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136 * September 2012
137 *
138 * .. Scalar Arguments ..
139  CHARACTER uplo
140  LOGICAL capply
141  INTEGER n, lda, ldaf, info
142 * ..
143 * .. Array Arguments ..
144  COMPLEX a( lda, * ), af( ldaf, * ), work( * )
145  REAL c( * ), rwork( * )
146 * ..
147 *
148 * =====================================================================
149 *
150 * .. Local Scalars ..
151  INTEGER kase
152  REAL ainvnm, anorm, tmp
153  INTEGER i, j
154  LOGICAL up, upper
155  COMPLEX zdum
156 * ..
157 * .. Local Arrays ..
158  INTEGER isave( 3 )
159 * ..
160 * .. External Functions ..
161  LOGICAL lsame
162  EXTERNAL lsame
163 * ..
164 * .. External Subroutines ..
165  EXTERNAL clacn2, cpotrs, xerbla
166 * ..
167 * .. Intrinsic Functions ..
168  INTRINSIC abs, max, REAL, aimag
169 * ..
170 * .. Statement Functions ..
171  REAL cabs1
172 * ..
173 * .. Statement Function Definitions ..
174  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
175 * ..
176 * .. Executable Statements ..
177 *
178  cla_porcond_c = 0.0e+0
179 *
180  info = 0
181  upper = lsame( uplo, 'U' )
182  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
183  info = -1
184  ELSE IF( n.LT.0 ) THEN
185  info = -2
186  ELSE IF( lda.LT.max( 1, n ) ) THEN
187  info = -4
188  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
189  info = -6
190  END IF
191  IF( info.NE.0 ) THEN
192  CALL xerbla( 'CLA_PORCOND_C', -info )
193  RETURN
194  END IF
195  up = .false.
196  IF ( lsame( uplo, 'U' ) ) up = .true.
197 *
198 * Compute norm of op(A)*op2(C).
199 *
200  anorm = 0.0e+0
201  IF ( up ) THEN
202  DO i = 1, n
203  tmp = 0.0e+0
204  IF ( capply ) THEN
205  DO j = 1, i
206  tmp = tmp + cabs1( a( j, i ) ) / c( j )
207  END DO
208  DO j = i+1, n
209  tmp = tmp + cabs1( a( i, j ) ) / c( j )
210  END DO
211  ELSE
212  DO j = 1, i
213  tmp = tmp + cabs1( a( j, i ) )
214  END DO
215  DO j = i+1, n
216  tmp = tmp + cabs1( a( i, j ) )
217  END DO
218  END IF
219  rwork( i ) = tmp
220  anorm = max( anorm, tmp )
221  END DO
222  ELSE
223  DO i = 1, n
224  tmp = 0.0e+0
225  IF ( capply ) THEN
226  DO j = 1, i
227  tmp = tmp + cabs1( a( i, j ) ) / c( j )
228  END DO
229  DO j = i+1, n
230  tmp = tmp + cabs1( a( j, i ) ) / c( j )
231  END DO
232  ELSE
233  DO j = 1, i
234  tmp = tmp + cabs1( a( i, j ) )
235  END DO
236  DO j = i+1, n
237  tmp = tmp + cabs1( a( j, i ) )
238  END DO
239  END IF
240  rwork( i ) = tmp
241  anorm = max( anorm, tmp )
242  END DO
243  END IF
244 *
245 * Quick return if possible.
246 *
247  IF( n.EQ.0 ) THEN
248  cla_porcond_c = 1.0e+0
249  RETURN
250  ELSE IF( anorm .EQ. 0.0e+0 ) THEN
251  RETURN
252  END IF
253 *
254 * Estimate the norm of inv(op(A)).
255 *
256  ainvnm = 0.0e+0
257 *
258  kase = 0
259  10 CONTINUE
260  CALL clacn2( n, work( n+1 ), work, ainvnm, kase, isave )
261  IF( kase.NE.0 ) THEN
262  IF( kase.EQ.2 ) THEN
263 *
264 * Multiply by R.
265 *
266  DO i = 1, n
267  work( i ) = work( i ) * rwork( i )
268  END DO
269 *
270  IF ( up ) THEN
271  CALL cpotrs( 'U', n, 1, af, ldaf,
272  $ work, n, info )
273  ELSE
274  CALL cpotrs( 'L', n, 1, af, ldaf,
275  $ work, n, info )
276  ENDIF
277 *
278 * Multiply by inv(C).
279 *
280  IF ( capply ) THEN
281  DO i = 1, n
282  work( i ) = work( i ) * c( i )
283  END DO
284  END IF
285  ELSE
286 *
287 * Multiply by inv(C**H).
288 *
289  IF ( capply ) THEN
290  DO i = 1, n
291  work( i ) = work( i ) * c( i )
292  END DO
293  END IF
294 *
295  IF ( up ) THEN
296  CALL cpotrs( 'U', n, 1, af, ldaf,
297  $ work, n, info )
298  ELSE
299  CALL cpotrs( 'L', n, 1, af, ldaf,
300  $ work, n, info )
301  END IF
302 *
303 * Multiply by R.
304 *
305  DO i = 1, n
306  work( i ) = work( i ) * rwork( i )
307  END DO
308  END IF
309  GO TO 10
310  END IF
311 *
312 * Compute the estimate of the reciprocal condition number.
313 *
314  IF( ainvnm .NE. 0.0e+0 )
315  $ cla_porcond_c = 1.0e+0 / ainvnm
316 *
317  RETURN
318 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine cpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
CPOTRS
Definition: cpotrs.f:112
real function cla_porcond_c(UPLO, N, A, LDA, AF, LDAF, C, CAPPLY, INFO, WORK, RWORK)
CLA_PORCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for Hermitian positiv...
subroutine clacn2(N, V, X, EST, KASE, ISAVE)
CLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: clacn2.f:135

Here is the call graph for this function:

Here is the caller graph for this function:

real function cla_porcond_x ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldaf, * )  AF,
integer  LDAF,
complex, dimension( * )  X,
integer  INFO,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK 
)

CLA_PORCOND_X computes the infinity norm condition number of op(A)*diag(x) for Hermitian positive-definite matrices.

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

Purpose:
    CLA_PORCOND_X Computes the infinity norm condition number of
    op(A) * diag(X) where X is a COMPLEX vector.
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]A
          A is COMPLEX array, dimension (LDA,N)
     On entry, the N-by-N matrix A.
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in]AF
          AF is COMPLEX array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**H*U or A = L*L**H, as computed by CPOTRF.
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in]X
          X is COMPLEX array, dimension (N)
     The vector X in the formula op(A) * diag(X).
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
     i > 0:  The ith argument is invalid.
[in]WORK
          WORK is COMPLEX array, dimension (2*N).
     Workspace.
[in]RWORK
          RWORK is REAL array, dimension (N).
     Workspace.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 125 of file cla_porcond_x.f.

125 *
126 * -- LAPACK computational routine (version 3.4.2) --
127 * -- LAPACK is a software package provided by Univ. of Tennessee, --
128 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
129 * September 2012
130 *
131 * .. Scalar Arguments ..
132  CHARACTER uplo
133  INTEGER n, lda, ldaf, info
134 * ..
135 * .. Array Arguments ..
136  COMPLEX a( lda, * ), af( ldaf, * ), work( * ), x( * )
137  REAL rwork( * )
138 * ..
139 *
140 * =====================================================================
141 *
142 * .. Local Scalars ..
143  INTEGER kase, i, j
144  REAL ainvnm, anorm, tmp
145  LOGICAL up, upper
146  COMPLEX zdum
147 * ..
148 * .. Local Arrays ..
149  INTEGER isave( 3 )
150 * ..
151 * .. External Functions ..
152  LOGICAL lsame
153  EXTERNAL lsame
154 * ..
155 * .. External Subroutines ..
156  EXTERNAL clacn2, cpotrs, xerbla
157 * ..
158 * .. Intrinsic Functions ..
159  INTRINSIC abs, max, REAL, aimag
160 * ..
161 * .. Statement Functions ..
162  REAL cabs1
163 * ..
164 * .. Statement Function Definitions ..
165  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
166 * ..
167 * .. Executable Statements ..
168 *
169  cla_porcond_x = 0.0e+0
170 *
171  info = 0
172  upper = lsame( uplo, 'U' )
173  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
174  info = -1
175  ELSE IF ( n.LT.0 ) THEN
176  info = -2
177  ELSE IF( lda.LT.max( 1, n ) ) THEN
178  info = -4
179  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
180  info = -6
181  END IF
182  IF( info.NE.0 ) THEN
183  CALL xerbla( 'CLA_PORCOND_X', -info )
184  RETURN
185  END IF
186  up = .false.
187  IF ( lsame( uplo, 'U' ) ) up = .true.
188 *
189 * Compute norm of op(A)*op2(C).
190 *
191  anorm = 0.0
192  IF ( up ) THEN
193  DO i = 1, n
194  tmp = 0.0e+0
195  DO j = 1, i
196  tmp = tmp + cabs1( a( j, i ) * x( j ) )
197  END DO
198  DO j = i+1, n
199  tmp = tmp + cabs1( a( i, j ) * x( j ) )
200  END DO
201  rwork( i ) = tmp
202  anorm = max( anorm, tmp )
203  END DO
204  ELSE
205  DO i = 1, n
206  tmp = 0.0e+0
207  DO j = 1, i
208  tmp = tmp + cabs1( a( i, j ) * x( j ) )
209  END DO
210  DO j = i+1, n
211  tmp = tmp + cabs1( a( j, i ) * x( j ) )
212  END DO
213  rwork( i ) = tmp
214  anorm = max( anorm, tmp )
215  END DO
216  END IF
217 *
218 * Quick return if possible.
219 *
220  IF( n.EQ.0 ) THEN
221  cla_porcond_x = 1.0e+0
222  RETURN
223  ELSE IF( anorm .EQ. 0.0e+0 ) THEN
224  RETURN
225  END IF
226 *
227 * Estimate the norm of inv(op(A)).
228 *
229  ainvnm = 0.0e+0
230 *
231  kase = 0
232  10 CONTINUE
233  CALL clacn2( n, work( n+1 ), work, ainvnm, kase, isave )
234  IF( kase.NE.0 ) THEN
235  IF( kase.EQ.2 ) THEN
236 *
237 * Multiply by R.
238 *
239  DO i = 1, n
240  work( i ) = work( i ) * rwork( i )
241  END DO
242 *
243  IF ( up ) THEN
244  CALL cpotrs( 'U', n, 1, af, ldaf,
245  $ work, n, info )
246  ELSE
247  CALL cpotrs( 'L', n, 1, af, ldaf,
248  $ work, n, info )
249  ENDIF
250 *
251 * Multiply by inv(X).
252 *
253  DO i = 1, n
254  work( i ) = work( i ) / x( i )
255  END DO
256  ELSE
257 *
258 * Multiply by inv(X**H).
259 *
260  DO i = 1, n
261  work( i ) = work( i ) / x( i )
262  END DO
263 *
264  IF ( up ) THEN
265  CALL cpotrs( 'U', n, 1, af, ldaf,
266  $ work, n, info )
267  ELSE
268  CALL cpotrs( 'L', n, 1, af, ldaf,
269  $ work, n, info )
270  END IF
271 *
272 * Multiply by R.
273 *
274  DO i = 1, n
275  work( i ) = work( i ) * rwork( i )
276  END DO
277  END IF
278  GO TO 10
279  END IF
280 *
281 * Compute the estimate of the reciprocal condition number.
282 *
283  IF( ainvnm .NE. 0.0e+0 )
284  $ cla_porcond_x = 1.0e+0 / ainvnm
285 *
286  RETURN
287 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function cla_porcond_x(UPLO, N, A, LDA, AF, LDAF, X, INFO, WORK, RWORK)
CLA_PORCOND_X computes the infinity norm condition number of op(A)*diag(x) for Hermitian positive-def...
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine cpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
CPOTRS
Definition: cpotrs.f:112
subroutine clacn2(N, V, X, EST, KASE, ISAVE)
CLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: clacn2.f:135

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cla_porfsx_extended ( integer  PREC_TYPE,
character  UPLO,
integer  N,
integer  NRHS,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldaf, * )  AF,
integer  LDAF,
logical  COLEQU,
real, dimension( * )  C,
complex, dimension( ldb, * )  B,
integer  LDB,
complex, dimension( ldy, * )  Y,
integer  LDY,
real, dimension( * )  BERR_OUT,
integer  N_NORMS,
real, dimension( nrhs, * )  ERR_BNDS_NORM,
real, dimension( nrhs, * )  ERR_BNDS_COMP,
complex, dimension( * )  RES,
real, dimension( * )  AYB,
complex, dimension( * )  DY,
complex, dimension( * )  Y_TAIL,
real  RCOND,
integer  ITHRESH,
real  RTHRESH,
real  DZ_UB,
logical  IGNORE_CWISE,
integer  INFO 
)

CLA_PORFSX_EXTENDED improves the computed solution to a system of linear equations for symmetric or Hermitian positive-definite matrices by performing extra-precise iterative refinement and provides error bounds and backward error estimates for the solution.

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

Purpose:
 CLA_PORFSX_EXTENDED improves the computed solution to a system of
 linear equations by performing extra-precise iterative refinement
 and provides error bounds and backward error estimates for the solution.
 This subroutine is called by CPORFSX to perform iterative refinement.
 In addition to normwise error bound, the code provides maximum
 componentwise error bound if possible. See comments for ERR_BNDS_NORM
 and ERR_BNDS_COMP for details of the error bounds. Note that this
 subroutine is only resonsible for setting the second fields of
 ERR_BNDS_NORM and ERR_BNDS_COMP.
Parameters
[in]PREC_TYPE
          PREC_TYPE is INTEGER
     Specifies the intermediate precision to be used in refinement.
     The value is defined by ILAPREC(P) where P is a CHARACTER and
     P    = 'S':  Single
          = 'D':  Double
          = 'I':  Indigenous
          = 'X', 'E':  Extra
[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.
[in]A
          A is COMPLEX array, dimension (LDA,N)
     On entry, the N-by-N matrix A.
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in]AF
          AF is COMPLEX array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by CPOTRF.
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in]COLEQU
          COLEQU is LOGICAL
     If .TRUE. then column equilibration was done to A before calling
     this routine. This is needed to compute the solution and error
     bounds correctly.
[in]C
          C is REAL array, dimension (N)
     The column scale factors for A. If COLEQU = .FALSE., C
     is not accessed. If C is input, each element of C 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]B
          B is COMPLEX array, dimension (LDB,NRHS)
     The right-hand-side matrix B.
[in]LDB
          LDB is INTEGER
     The leading dimension of the array B.  LDB >= max(1,N).
[in,out]Y
          Y is COMPLEX array, dimension
                    (LDY,NRHS)
     On entry, the solution matrix X, as computed by CPOTRS.
     On exit, the improved solution matrix Y.
[in]LDY
          LDY is INTEGER
     The leading dimension of the array Y.  LDY >= max(1,N).
[out]BERR_OUT
          BERR_OUT is REAL array, dimension (NRHS)
     On exit, BERR_OUT(j) contains the componentwise relative backward
     error for right-hand-side j from the formula
         max(i) ( abs(RES(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
     where abs(Z) is the componentwise absolute value of the matrix
     or vector Z. This is computed by CLA_LIN_BERR.
[in]N_NORMS
          N_NORMS is INTEGER
     Determines which error bounds to return (see ERR_BNDS_NORM
     and ERR_BNDS_COMP).
     If N_NORMS >= 1 return normwise error bounds.
     If N_NORMS >= 2 return componentwise error bounds.
[in,out]ERR_BNDS_NORM
          ERR_BNDS_NORM is REAL 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 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" 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) * slamch('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) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". 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.

     This subroutine is only responsible for setting the second field
     above.
     See Lapack Working Note 165 for further details and extra
     cautions.
[in,out]ERR_BNDS_COMP
          ERR_BNDS_COMP is REAL 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 .LT. 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 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" 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) * slamch('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) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". 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.

     This subroutine is only responsible for setting the second field
     above.
     See Lapack Working Note 165 for further details and extra
     cautions.
[in]RES
          RES is COMPLEX array, dimension (N)
     Workspace to hold the intermediate residual.
[in]AYB
          AYB is REAL array, dimension (N)
     Workspace.
[in]DY
          DY is COMPLEX array, dimension (N)
     Workspace to hold the intermediate solution.
[in]Y_TAIL
          Y_TAIL is COMPLEX array, dimension (N)
     Workspace to hold the trailing bits of the intermediate solution.
[in]RCOND
          RCOND is REAL
     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.
[in]ITHRESH
          ITHRESH is INTEGER
     The maximum number of residual computations allowed for
     refinement. The default is 10. For '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.
[in]RTHRESH
          RTHRESH is REAL
     Determines when to stop refinement if the error estimate stops
     decreasing. Refinement will stop when the next solution no longer
     satisfies norm(dx_{i+1}) < RTHRESH * norm(dx_i) where norm(Z) is
     the infinity norm of Z. RTHRESH satisfies 0 < RTHRESH <= 1. The
     default value is 0.5. For 'aggressive' set to 0.9 to permit
     convergence on extremely ill-conditioned matrices. See LAWN 165
     for more details.
[in]DZ_UB
          DZ_UB is REAL
     Determines when to start considering componentwise convergence.
     Componentwise convergence is only considered after each component
     of the solution Y is stable, which we definte as the relative
     change in each component being less than DZ_UB. The default value
     is 0.25, requiring the first bit to be stable. See LAWN 165 for
     more details.
[in]IGNORE_CWISE
          IGNORE_CWISE is LOGICAL
     If .TRUE. then ignore componentwise convergence. Default value
     is .FALSE..
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit.
       < 0:  if INFO = -i, the ith argument to CPOTRS had an illegal
             value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 392 of file cla_porfsx_extended.f.

392 *
393 * -- LAPACK computational routine (version 3.4.2) --
394 * -- LAPACK is a software package provided by Univ. of Tennessee, --
395 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
396 * September 2012
397 *
398 * .. Scalar Arguments ..
399  INTEGER info, lda, ldaf, ldb, ldy, n, nrhs, prec_type,
400  $ n_norms, ithresh
401  CHARACTER uplo
402  LOGICAL colequ, ignore_cwise
403  REAL rthresh, dz_ub
404 * ..
405 * .. Array Arguments ..
406  COMPLEX a( lda, * ), af( ldaf, * ), b( ldb, * ),
407  $ y( ldy, * ), res( * ), dy( * ), y_tail( * )
408  REAL c( * ), ayb( * ), rcond, berr_out( * ),
409  $ err_bnds_norm( nrhs, * ),
410  $ err_bnds_comp( nrhs, * )
411 * ..
412 *
413 * =====================================================================
414 *
415 * .. Local Scalars ..
416  INTEGER uplo2, cnt, i, j, x_state, z_state,
417  $ y_prec_state
418  REAL yk, dyk, ymin, normy, normx, normdx, dxrat,
419  $ dzrat, prevnormdx, prev_dz_z, dxratmax,
420  $ dzratmax, dx_x, dz_z, final_dx_x, final_dz_z,
421  $ eps, hugeval, incr_thresh
422  LOGICAL incr_prec
423  COMPLEX zdum
424 * ..
425 * .. Parameters ..
426  INTEGER unstable_state, working_state, conv_state,
427  $ noprog_state, base_residual, extra_residual,
428  $ extra_y
429  parameter( unstable_state = 0, working_state = 1,
430  $ conv_state = 2, noprog_state = 3 )
431  parameter( base_residual = 0, extra_residual = 1,
432  $ extra_y = 2 )
433  INTEGER final_nrm_err_i, final_cmp_err_i, berr_i
434  INTEGER rcond_i, nrm_rcond_i, nrm_err_i, cmp_rcond_i
435  INTEGER cmp_err_i, piv_growth_i
436  parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
437  $ berr_i = 3 )
438  parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
439  parameter( cmp_rcond_i = 7, cmp_err_i = 8,
440  $ piv_growth_i = 9 )
441  INTEGER la_linrx_itref_i, la_linrx_ithresh_i,
442  $ la_linrx_cwise_i
443  parameter( la_linrx_itref_i = 1,
444  $ la_linrx_ithresh_i = 2 )
445  parameter( la_linrx_cwise_i = 3 )
446  INTEGER la_linrx_trust_i, la_linrx_err_i,
447  $ la_linrx_rcond_i
448  parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
449  parameter( la_linrx_rcond_i = 3 )
450 * ..
451 * .. External Functions ..
452  LOGICAL lsame
453  EXTERNAL ilauplo
454  INTEGER ilauplo
455 * ..
456 * .. External Subroutines ..
457  EXTERNAL caxpy, ccopy, cpotrs, chemv, blas_chemv_x,
458  $ blas_chemv2_x, cla_heamv, cla_wwaddw,
460  REAL slamch
461 * ..
462 * .. Intrinsic Functions ..
463  INTRINSIC abs, REAL, aimag, max, min
464 * ..
465 * .. Statement Functions ..
466  REAL cabs1
467 * ..
468 * .. Statement Function Definitions ..
469  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
470 * ..
471 * .. Executable Statements ..
472 *
473  IF (info.NE.0) RETURN
474  eps = slamch( 'Epsilon' )
475  hugeval = slamch( 'Overflow' )
476 * Force HUGEVAL to Inf
477  hugeval = hugeval * hugeval
478 * Using HUGEVAL may lead to spurious underflows.
479  incr_thresh = REAL(N) * eps
480 
481  IF (lsame(uplo, 'L')) THEN
482  uplo2 = ilauplo( 'L' )
483  ELSE
484  uplo2 = ilauplo( 'U' )
485  ENDIF
486 
487  DO j = 1, nrhs
488  y_prec_state = extra_residual
489  IF (y_prec_state .EQ. extra_y) THEN
490  DO i = 1, n
491  y_tail( i ) = 0.0
492  END DO
493  END IF
494 
495  dxrat = 0.0
496  dxratmax = 0.0
497  dzrat = 0.0
498  dzratmax = 0.0
499  final_dx_x = hugeval
500  final_dz_z = hugeval
501  prevnormdx = hugeval
502  prev_dz_z = hugeval
503  dz_z = hugeval
504  dx_x = hugeval
505 
506  x_state = working_state
507  z_state = unstable_state
508  incr_prec = .false.
509 
510  DO cnt = 1, ithresh
511 *
512 * Compute residual RES = B_s - op(A_s) * Y,
513 * op(A) = A, A**T, or A**H depending on TRANS (and type).
514 *
515  CALL ccopy( n, b( 1, j ), 1, res, 1 )
516  IF (y_prec_state .EQ. base_residual) THEN
517  CALL chemv(uplo, n, cmplx(-1.0), a, lda, y(1,j), 1,
518  $ cmplx(1.0), res, 1)
519  ELSE IF (y_prec_state .EQ. extra_residual) THEN
520  CALL blas_chemv_x(uplo2, n, cmplx(-1.0), a, lda,
521  $ y( 1, j ), 1, cmplx(1.0), res, 1, prec_type)
522  ELSE
523  CALL blas_chemv2_x(uplo2, n, cmplx(-1.0), a, lda,
524  $ y(1, j), y_tail, 1, cmplx(1.0), res, 1, prec_type)
525  END IF
526 
527 ! XXX: RES is no longer needed.
528  CALL ccopy( n, res, 1, dy, 1 )
529  CALL cpotrs( uplo, n, 1, af, ldaf, dy, n, info)
530 *
531 * Calculate relative changes DX_X, DZ_Z and ratios DXRAT, DZRAT.
532 *
533  normx = 0.0
534  normy = 0.0
535  normdx = 0.0
536  dz_z = 0.0
537  ymin = hugeval
538 
539  DO i = 1, n
540  yk = cabs1(y(i, j))
541  dyk = cabs1(dy(i))
542 
543  IF (yk .NE. 0.0) THEN
544  dz_z = max( dz_z, dyk / yk )
545  ELSE IF (dyk .NE. 0.0) THEN
546  dz_z = hugeval
547  END IF
548 
549  ymin = min( ymin, yk )
550 
551  normy = max( normy, yk )
552 
553  IF ( colequ ) THEN
554  normx = max(normx, yk * c(i))
555  normdx = max(normdx, dyk * c(i))
556  ELSE
557  normx = normy
558  normdx = max(normdx, dyk)
559  END IF
560  END DO
561 
562  IF (normx .NE. 0.0) THEN
563  dx_x = normdx / normx
564  ELSE IF (normdx .EQ. 0.0) THEN
565  dx_x = 0.0
566  ELSE
567  dx_x = hugeval
568  END IF
569 
570  dxrat = normdx / prevnormdx
571  dzrat = dz_z / prev_dz_z
572 *
573 * Check termination criteria.
574 *
575  IF (ymin*rcond .LT. incr_thresh*normy
576  $ .AND. y_prec_state .LT. extra_y)
577  $ incr_prec = .true.
578 
579  IF (x_state .EQ. noprog_state .AND. dxrat .LE. rthresh)
580  $ x_state = working_state
581  IF (x_state .EQ. working_state) THEN
582  IF (dx_x .LE. eps) THEN
583  x_state = conv_state
584  ELSE IF (dxrat .GT. rthresh) THEN
585  IF (y_prec_state .NE. extra_y) THEN
586  incr_prec = .true.
587  ELSE
588  x_state = noprog_state
589  END IF
590  ELSE
591  IF (dxrat .GT. dxratmax) dxratmax = dxrat
592  END IF
593  IF (x_state .GT. working_state) final_dx_x = dx_x
594  END IF
595 
596  IF (z_state .EQ. unstable_state .AND. dz_z .LE. dz_ub)
597  $ z_state = working_state
598  IF (z_state .EQ. noprog_state .AND. dzrat .LE. rthresh)
599  $ z_state = working_state
600  IF (z_state .EQ. working_state) THEN
601  IF (dz_z .LE. eps) THEN
602  z_state = conv_state
603  ELSE IF (dz_z .GT. dz_ub) THEN
604  z_state = unstable_state
605  dzratmax = 0.0
606  final_dz_z = hugeval
607  ELSE IF (dzrat .GT. rthresh) THEN
608  IF (y_prec_state .NE. extra_y) THEN
609  incr_prec = .true.
610  ELSE
611  z_state = noprog_state
612  END IF
613  ELSE
614  IF (dzrat .GT. dzratmax) dzratmax = dzrat
615  END IF
616  IF (z_state .GT. working_state) final_dz_z = dz_z
617  END IF
618 
619  IF ( x_state.NE.working_state.AND.
620  $ (ignore_cwise.OR.z_state.NE.working_state) )
621  $ GOTO 666
622 
623  IF (incr_prec) THEN
624  incr_prec = .false.
625  y_prec_state = y_prec_state + 1
626  DO i = 1, n
627  y_tail( i ) = 0.0
628  END DO
629  END IF
630 
631  prevnormdx = normdx
632  prev_dz_z = dz_z
633 *
634 * Update soluton.
635 *
636  IF (y_prec_state .LT. extra_y) THEN
637  CALL caxpy( n, cmplx(1.0), dy, 1, y(1,j), 1 )
638  ELSE
639  CALL cla_wwaddw(n, y(1,j), y_tail, dy)
640  END IF
641 
642  END DO
643 * Target of "IF (Z_STOP .AND. X_STOP)". Sun's f77 won't EXIT.
644  666 CONTINUE
645 *
646 * Set final_* when cnt hits ithresh.
647 *
648  IF (x_state .EQ. working_state) final_dx_x = dx_x
649  IF (z_state .EQ. working_state) final_dz_z = dz_z
650 *
651 * Compute error bounds.
652 *
653  IF (n_norms .GE. 1) THEN
654  err_bnds_norm( j, la_linrx_err_i ) =
655  $ final_dx_x / (1 - dxratmax)
656  END IF
657  IF (n_norms .GE. 2) THEN
658  err_bnds_comp( j, la_linrx_err_i ) =
659  $ final_dz_z / (1 - dzratmax)
660  END IF
661 *
662 * Compute componentwise relative backward error from formula
663 * max(i) ( abs(R(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
664 * where abs(Z) is the componentwise absolute value of the matrix
665 * or vector Z.
666 *
667 * Compute residual RES = B_s - op(A_s) * Y,
668 * op(A) = A, A**T, or A**H depending on TRANS (and type).
669 *
670  CALL ccopy( n, b( 1, j ), 1, res, 1 )
671  CALL chemv(uplo, n, cmplx(-1.0), a, lda, y(1,j), 1, cmplx(1.0),
672  $ res, 1)
673 
674  DO i = 1, n
675  ayb( i ) = cabs1( b( i, j ) )
676  END DO
677 *
678 * Compute abs(op(A_s))*abs(Y) + abs(B_s).
679 *
680  CALL cla_heamv (uplo2, n, 1.0,
681  $ a, lda, y(1, j), 1, 1.0, ayb, 1)
682 
683  CALL cla_lin_berr (n, n, 1, res, ayb, berr_out(j))
684 *
685 * End of loop for each RHS.
686 *
687  END DO
688 *
689  RETURN
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine cla_wwaddw(N, X, Y, W)
CLA_WWADDW adds a vector into a doubled-single vector.
Definition: cla_wwaddw.f:83
integer function ilauplo(UPLO)
ILAUPLO
Definition: ilauplo.f:60
subroutine cla_heamv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CLA_HEAMV computes a matrix-vector product using a Hermitian indefinite matrix to calculate error bou...
Definition: cla_heamv.f:180
subroutine chemv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CHEMV
Definition: chemv.f:156
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine cpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
CPOTRS
Definition: cpotrs.f:112
subroutine caxpy(N, CA, CX, INCX, CY, INCY)
CAXPY
Definition: caxpy.f:53
subroutine cla_lin_berr(N, NZ, NRHS, RES, AYB, BERR)
CLA_LIN_BERR computes a component-wise relative backward error.
Definition: cla_lin_berr.f:103
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:52

Here is the call graph for this function:

Here is the caller graph for this function:

real function cla_porpvgrw ( character*1  UPLO,
integer  NCOLS,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldaf, * )  AF,
integer  LDAF,
real, dimension( * )  WORK 
)

CLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian positive-definite matrix.

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

Purpose:
 CLA_PORPVGRW computes the reciprocal pivot growth factor
 norm(A)/norm(U). The "max absolute element" norm is used. If this is
 much less than 1, 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.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored.
[in]NCOLS
          NCOLS is INTEGER
     The number of columns of the matrix A. NCOLS >= 0.
[in]A
          A is COMPLEX array, dimension (LDA,N)
     On entry, the N-by-N matrix A.
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in]AF
          AF is COMPLEX array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by CPOTRF.
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in]WORK
          WORK is COMPLEX array, dimension (2*N)
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 107 of file cla_porpvgrw.f.

107 *
108 * -- LAPACK computational routine (version 3.4.2) --
109 * -- LAPACK is a software package provided by Univ. of Tennessee, --
110 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
111 * September 2012
112 *
113 * .. Scalar Arguments ..
114  CHARACTER*1 uplo
115  INTEGER ncols, lda, ldaf
116 * ..
117 * .. Array Arguments ..
118  COMPLEX a( lda, * ), af( ldaf, * )
119  REAL work( * )
120 * ..
121 *
122 * =====================================================================
123 *
124 * .. Local Scalars ..
125  INTEGER i, j
126  REAL amax, umax, rpvgrw
127  LOGICAL upper
128  COMPLEX zdum
129 * ..
130 * .. External Functions ..
131  EXTERNAL lsame, claset
132  LOGICAL lsame
133 * ..
134 * .. Intrinsic Functions ..
135  INTRINSIC abs, max, min, REAL, aimag
136 * ..
137 * .. Statement Functions ..
138  REAL cabs1
139 * ..
140 * .. Statement Function Definitions ..
141  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
142 * ..
143 * .. Executable Statements ..
144  upper = lsame( 'Upper', uplo )
145 *
146 * SPOTRF will have factored only the NCOLSxNCOLS leading minor, so
147 * we restrict the growth search to that minor and use only the first
148 * 2*NCOLS workspace entries.
149 *
150  rpvgrw = 1.0
151  DO i = 1, 2*ncols
152  work( i ) = 0.0
153  END DO
154 *
155 * Find the max magnitude entry of each column.
156 *
157  IF ( upper ) THEN
158  DO j = 1, ncols
159  DO i = 1, j
160  work( ncols+j ) =
161  $ max( cabs1( a( i, j ) ), work( ncols+j ) )
162  END DO
163  END DO
164  ELSE
165  DO j = 1, ncols
166  DO i = j, ncols
167  work( ncols+j ) =
168  $ max( cabs1( a( i, j ) ), work( ncols+j ) )
169  END DO
170  END DO
171  END IF
172 *
173 * Now find the max magnitude entry of each column of the factor in
174 * AF. No pivoting, so no permutations.
175 *
176  IF ( lsame( 'Upper', uplo ) ) THEN
177  DO j = 1, ncols
178  DO i = 1, j
179  work( j ) = max( cabs1( af( i, j ) ), work( j ) )
180  END DO
181  END DO
182  ELSE
183  DO j = 1, ncols
184  DO i = j, ncols
185  work( j ) = max( cabs1( af( i, j ) ), work( j ) )
186  END DO
187  END DO
188  END IF
189 *
190 * Compute the *inverse* of the max element growth factor. Dividing
191 * by zero would imply the largest entry of the factor's column is
192 * zero. Than can happen when either the column of A is zero or
193 * massive pivots made the factor underflow to zero. Neither counts
194 * as growth in itself, so simply ignore terms with zero
195 * denominators.
196 *
197  IF ( lsame( 'Upper', uplo ) ) THEN
198  DO i = 1, ncols
199  umax = work( i )
200  amax = work( ncols+i )
201  IF ( umax /= 0.0 ) THEN
202  rpvgrw = min( amax / umax, rpvgrw )
203  END IF
204  END DO
205  ELSE
206  DO i = 1, ncols
207  umax = work( i )
208  amax = work( ncols+i )
209  IF ( umax /= 0.0 ) THEN
210  rpvgrw = min( amax / umax, rpvgrw )
211  END IF
212  END DO
213  END IF
214 
215  cla_porpvgrw = rpvgrw
subroutine claset(UPLO, M, N, ALPHA, BETA, A, LDA)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: claset.f:108
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
real function cla_porpvgrw(UPLO, NCOLS, A, LDA, AF, LDAF, WORK)
CLA_PORPVGRW computes the reciprocal pivot growth factor norm(A)/norm(U) for a symmetric or Hermitian...
Definition: cla_porpvgrw.f:107

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cpocon ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
real  ANORM,
real  RCOND,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CPOCON

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

Purpose:
 CPOCON estimates the reciprocal of the condition number (in the
 1-norm) of a complex Hermitian positive definite matrix using the
 Cholesky factorization A = U**H*U or A = L*L**H computed by CPOTRF.

 An estimate is obtained for norm(inv(A)), and the reciprocal of the
 condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
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]A
          A is COMPLEX array, dimension (LDA,N)
          The triangular factor U or L from the Cholesky factorization
          A = U**H*U or A = L*L**H, as computed by CPOTRF.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in]ANORM
          ANORM is REAL
          The 1-norm (or infinity-norm) of the Hermitian matrix A.
[out]RCOND
          RCOND is REAL
          The reciprocal of the condition number of the matrix A,
          computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is an
          estimate of the 1-norm of inv(A) computed in this routine.
[out]WORK
          WORK is COMPLEX array, dimension (2*N)
[out]RWORK
          RWORK is REAL array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 123 of file cpocon.f.

123 *
124 * -- LAPACK computational routine (version 3.4.0) --
125 * -- LAPACK is a software package provided by Univ. of Tennessee, --
126 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
127 * November 2011
128 *
129 * .. Scalar Arguments ..
130  CHARACTER uplo
131  INTEGER info, lda, n
132  REAL anorm, rcond
133 * ..
134 * .. Array Arguments ..
135  REAL rwork( * )
136  COMPLEX a( lda, * ), work( * )
137 * ..
138 *
139 * =====================================================================
140 *
141 * .. Parameters ..
142  REAL one, zero
143  parameter( one = 1.0e+0, zero = 0.0e+0 )
144 * ..
145 * .. Local Scalars ..
146  LOGICAL upper
147  CHARACTER normin
148  INTEGER ix, kase
149  REAL ainvnm, scale, scalel, scaleu, smlnum
150  COMPLEX zdum
151 * ..
152 * .. Local Arrays ..
153  INTEGER isave( 3 )
154 * ..
155 * .. External Functions ..
156  LOGICAL lsame
157  INTEGER icamax
158  REAL slamch
159  EXTERNAL lsame, icamax, slamch
160 * ..
161 * .. External Subroutines ..
162  EXTERNAL clacn2, clatrs, csrscl, xerbla
163 * ..
164 * .. Intrinsic Functions ..
165  INTRINSIC abs, aimag, max, real
166 * ..
167 * .. Statement Functions ..
168  REAL cabs1
169 * ..
170 * .. Statement Function definitions ..
171  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
172 * ..
173 * .. Executable Statements ..
174 *
175 * Test the input parameters.
176 *
177  info = 0
178  upper = lsame( uplo, 'U' )
179  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
180  info = -1
181  ELSE IF( n.LT.0 ) THEN
182  info = -2
183  ELSE IF( lda.LT.max( 1, n ) ) THEN
184  info = -4
185  ELSE IF( anorm.LT.zero ) THEN
186  info = -5
187  END IF
188  IF( info.NE.0 ) THEN
189  CALL xerbla( 'CPOCON', -info )
190  RETURN
191  END IF
192 *
193 * Quick return if possible
194 *
195  rcond = zero
196  IF( n.EQ.0 ) THEN
197  rcond = one
198  RETURN
199  ELSE IF( anorm.EQ.zero ) THEN
200  RETURN
201  END IF
202 *
203  smlnum = slamch( 'Safe minimum' )
204 *
205 * Estimate the 1-norm of inv(A).
206 *
207  kase = 0
208  normin = 'N'
209  10 CONTINUE
210  CALL clacn2( n, work( n+1 ), work, ainvnm, kase, isave )
211  IF( kase.NE.0 ) THEN
212  IF( upper ) THEN
213 *
214 * Multiply by inv(U**H).
215 *
216  CALL clatrs( 'Upper', 'Conjugate transpose', 'Non-unit',
217  $ normin, n, a, lda, work, scalel, rwork, info )
218  normin = 'Y'
219 *
220 * Multiply by inv(U).
221 *
222  CALL clatrs( 'Upper', 'No transpose', 'Non-unit', normin, n,
223  $ a, lda, work, scaleu, rwork, info )
224  ELSE
225 *
226 * Multiply by inv(L).
227 *
228  CALL clatrs( 'Lower', 'No transpose', 'Non-unit', normin, n,
229  $ a, lda, work, scalel, rwork, info )
230  normin = 'Y'
231 *
232 * Multiply by inv(L**H).
233 *
234  CALL clatrs( 'Lower', 'Conjugate transpose', 'Non-unit',
235  $ normin, n, a, lda, work, scaleu, rwork, info )
236  END IF
237 *
238 * Multiply by 1/SCALE if doing so will not cause overflow.
239 *
240  scale = scalel*scaleu
241  IF( scale.NE.one ) THEN
242  ix = icamax( n, work, 1 )
243  IF( scale.LT.cabs1( work( ix ) )*smlnum .OR. scale.EQ.zero )
244  $ GO TO 20
245  CALL csrscl( n, scale, work, 1 )
246  END IF
247  GO TO 10
248  END IF
249 *
250 * Compute the estimate of the reciprocal condition number.
251 *
252  IF( ainvnm.NE.zero )
253  $ rcond = ( one / ainvnm ) / anorm
254 *
255  20 CONTINUE
256  RETURN
257 *
258 * End of CPOCON
259 *
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
integer function icamax(N, CX, INCX)
ICAMAX
Definition: icamax.f:53
subroutine csrscl(N, SA, SX, INCX)
CSRSCL multiplies a vector by the reciprocal of a real scalar.
Definition: csrscl.f:86
subroutine clatrs(UPLO, TRANS, DIAG, NORMIN, N, A, LDA, X, SCALE, CNORM, INFO)
CLATRS solves a triangular system of equations with the scale factor set to prevent overflow...
Definition: clatrs.f:241
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine clacn2(N, V, X, EST, KASE, ISAVE)
CLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: clacn2.f:135

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cpoequ ( integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
real, dimension( * )  S,
real  SCOND,
real  AMAX,
integer  INFO 
)

CPOEQU

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

Purpose:
 CPOEQU computes row and column scalings intended to equilibrate a
 Hermitian positive definite matrix A and reduce its condition number
 (with respect to the two-norm).  S contains the scale factors,
 S(i) = 1/sqrt(A(i,i)), chosen so that the scaled matrix B with
 elements B(i,j) = S(i)*A(i,j)*S(j) has ones on the diagonal.  This
 choice of S puts the condition number of B within a factor N of the
 smallest possible condition number over all possible diagonal
 scalings.
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]A
          A is COMPLEX array, dimension (LDA,N)
          The N-by-N Hermitian positive definite matrix whose scaling
          factors are to be computed.  Only the diagonal elements of A
          are referenced.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]S
          S is REAL array, dimension (N)
          If INFO = 0, S contains the scale factors for A.
[out]SCOND
          SCOND is REAL
          If INFO = 0, S contains the ratio of the smallest S(i) to
          the largest S(i).  If SCOND >= 0.1 and AMAX is neither too
          large nor too small, it is not worth scaling by S.
[out]AMAX
          AMAX is REAL
          Absolute value of largest matrix element.  If AMAX is very
          close to overflow or very close to underflow, the matrix
          should be scaled.
[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 i-th diagonal element is nonpositive.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 115 of file cpoequ.f.

115 *
116 * -- LAPACK computational routine (version 3.4.0) --
117 * -- LAPACK is a software package provided by Univ. of Tennessee, --
118 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
119 * November 2011
120 *
121 * .. Scalar Arguments ..
122  INTEGER info, lda, n
123  REAL amax, scond
124 * ..
125 * .. Array Arguments ..
126  REAL s( * )
127  COMPLEX a( lda, * )
128 * ..
129 *
130 * =====================================================================
131 *
132 * .. Parameters ..
133  REAL zero, one
134  parameter( zero = 0.0e+0, one = 1.0e+0 )
135 * ..
136 * .. Local Scalars ..
137  INTEGER i
138  REAL smin
139 * ..
140 * .. External Subroutines ..
141  EXTERNAL xerbla
142 * ..
143 * .. Intrinsic Functions ..
144  INTRINSIC max, min, REAL, sqrt
145 * ..
146 * .. Executable Statements ..
147 *
148 * Test the input parameters.
149 *
150  info = 0
151  IF( n.LT.0 ) THEN
152  info = -1
153  ELSE IF( lda.LT.max( 1, n ) ) THEN
154  info = -3
155  END IF
156  IF( info.NE.0 ) THEN
157  CALL xerbla( 'CPOEQU', -info )
158  RETURN
159  END IF
160 *
161 * Quick return if possible
162 *
163  IF( n.EQ.0 ) THEN
164  scond = one
165  amax = zero
166  RETURN
167  END IF
168 *
169 * Find the minimum and maximum diagonal elements.
170 *
171  s( 1 ) = REAL( A( 1, 1 ) )
172  smin = s( 1 )
173  amax = s( 1 )
174  DO 10 i = 2, n
175  s( i ) = REAL( A( I, I ) )
176  smin = min( smin, s( i ) )
177  amax = max( amax, s( i ) )
178  10 CONTINUE
179 *
180  IF( smin.LE.zero ) THEN
181 *
182 * Find the first non-positive diagonal element and return.
183 *
184  DO 20 i = 1, n
185  IF( s( i ).LE.zero ) THEN
186  info = i
187  RETURN
188  END IF
189  20 CONTINUE
190  ELSE
191 *
192 * Set the scale factors to the reciprocals
193 * of the diagonal elements.
194 *
195  DO 30 i = 1, n
196  s( i ) = one / sqrt( s( i ) )
197  30 CONTINUE
198 *
199 * Compute SCOND = min(S(I)) / max(S(I))
200 *
201  scond = sqrt( smin ) / sqrt( amax )
202  END IF
203  RETURN
204 *
205 * End of CPOEQU
206 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cpoequb ( integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
real, dimension( * )  S,
real  SCOND,
real  AMAX,
integer  INFO 
)

CPOEQUB

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

Purpose:
 CPOEQUB computes row and column scalings intended to equilibrate a
 symmetric positive definite matrix A and reduce its condition number
 (with respect to the two-norm).  S contains the scale factors,
 S(i) = 1/sqrt(A(i,i)), chosen so that the scaled matrix B with
 elements B(i,j) = S(i)*A(i,j)*S(j) has ones on the diagonal.  This
 choice of S puts the condition number of B within a factor N of the
 smallest possible condition number over all possible diagonal
 scalings.
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]A
          A is COMPLEX array, dimension (LDA,N)
          The N-by-N symmetric positive definite matrix whose scaling
          factors are to be computed.  Only the diagonal elements of A
          are referenced.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]S
          S is REAL array, dimension (N)
          If INFO = 0, S contains the scale factors for A.
[out]SCOND
          SCOND is REAL
          If INFO = 0, S contains the ratio of the smallest S(i) to
          the largest S(i).  If SCOND >= 0.1 and AMAX is neither too
          large nor too small, it is not worth scaling by S.
[out]AMAX
          AMAX is REAL
          Absolute value of largest matrix element.  If AMAX is very
          close to overflow or very close to underflow, the matrix
          should be scaled.
[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 i-th diagonal element is nonpositive.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 115 of file cpoequb.f.

115 *
116 * -- LAPACK computational routine (version 3.4.0) --
117 * -- LAPACK is a software package provided by Univ. of Tennessee, --
118 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
119 * November 2011
120 *
121 * .. Scalar Arguments ..
122  INTEGER info, lda, n
123  REAL amax, scond
124 * ..
125 * .. Array Arguments ..
126  COMPLEX a( lda, * )
127  REAL s( * )
128 * ..
129 *
130 * =====================================================================
131 *
132 * .. Parameters ..
133  REAL zero, one
134  parameter( zero = 0.0e+0, one = 1.0e+0 )
135 * ..
136 * .. Local Scalars ..
137  INTEGER i
138  REAL smin, base, tmp
139 * ..
140 * .. External Functions ..
141  REAL slamch
142  EXTERNAL slamch
143 * ..
144 * .. External Subroutines ..
145  EXTERNAL xerbla
146 * ..
147 * .. Intrinsic Functions ..
148  INTRINSIC max, min, sqrt, log, int
149 * ..
150 * .. Executable Statements ..
151 *
152 * Test the input parameters.
153 *
154 * Positive definite only performs 1 pass of equilibration.
155 *
156  info = 0
157  IF( n.LT.0 ) THEN
158  info = -1
159  ELSE IF( lda.LT.max( 1, n ) ) THEN
160  info = -3
161  END IF
162  IF( info.NE.0 ) THEN
163  CALL xerbla( 'CPOEQUB', -info )
164  RETURN
165  END IF
166 *
167 * Quick return if possible.
168 *
169  IF( n.EQ.0 ) THEN
170  scond = one
171  amax = zero
172  RETURN
173  END IF
174 
175  base = slamch( 'B' )
176  tmp = -0.5 / log( base )
177 *
178 * Find the minimum and maximum diagonal elements.
179 *
180  s( 1 ) = a( 1, 1 )
181  smin = s( 1 )
182  amax = s( 1 )
183  DO 10 i = 2, n
184  s( i ) = a( i, i )
185  smin = min( smin, s( i ) )
186  amax = max( amax, s( i ) )
187  10 CONTINUE
188 *
189  IF( smin.LE.zero ) THEN
190 *
191 * Find the first non-positive diagonal element and return.
192 *
193  DO 20 i = 1, n
194  IF( s( i ).LE.zero ) THEN
195  info = i
196  RETURN
197  END IF
198  20 CONTINUE
199  ELSE
200 *
201 * Set the scale factors to the reciprocals
202 * of the diagonal elements.
203 *
204  DO 30 i = 1, n
205  s( i ) = base ** int( tmp * log( s( i ) ) )
206  30 CONTINUE
207 *
208 * Compute SCOND = min(S(I)) / max(S(I)).
209 *
210  scond = sqrt( smin ) / sqrt( amax )
211  END IF
212 *
213  RETURN
214 *
215 * End of CPOEQUB
216 *
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cporfs ( character  UPLO,
integer  N,
integer  NRHS,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldaf, * )  AF,
integer  LDAF,
complex, dimension( ldb, * )  B,
integer  LDB,
complex, dimension( ldx, * )  X,
integer  LDX,
real, dimension( * )  FERR,
real, dimension( * )  BERR,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CPORFS

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

Purpose:
 CPORFS improves the computed solution to a system of linear
 equations when the coefficient matrix is Hermitian positive definite,
 and provides error bounds and backward error estimates for the
 solution.
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 matrices B and X.  NRHS >= 0.
[in]A
          A is COMPLEX 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]AF
          AF is COMPLEX array, dimension (LDAF,N)
          The triangular factor U or L from the Cholesky factorization
          A = U**H*U or A = L*L**H, as computed by CPOTRF.
[in]LDAF
          LDAF is INTEGER
          The leading dimension of the array AF.  LDAF >= max(1,N).
[in]B
          B is COMPLEX array, dimension (LDB,NRHS)
          The right hand side matrix B.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in,out]X
          X is COMPLEX array, dimension (LDX,NRHS)
          On entry, the solution matrix X, as computed by CPOTRS.
          On exit, the improved solution matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[out]FERR
          FERR is REAL 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 REAL 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 array, dimension (2*N)
[out]RWORK
          RWORK is REAL array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Internal Parameters:
  ITMAX is the maximum number of steps of iterative refinement.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 185 of file cporfs.f.

185 *
186 * -- LAPACK computational routine (version 3.4.0) --
187 * -- LAPACK is a software package provided by Univ. of Tennessee, --
188 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189 * November 2011
190 *
191 * .. Scalar Arguments ..
192  CHARACTER uplo
193  INTEGER info, lda, ldaf, ldb, ldx, n, nrhs
194 * ..
195 * .. Array Arguments ..
196  REAL berr( * ), ferr( * ), rwork( * )
197  COMPLEX a( lda, * ), af( ldaf, * ), b( ldb, * ),
198  $ work( * ), x( ldx, * )
199 * ..
200 *
201 * ====================================================================
202 *
203 * .. Parameters ..
204  INTEGER itmax
205  parameter( itmax = 5 )
206  REAL zero
207  parameter( zero = 0.0e+0 )
208  COMPLEX one
209  parameter( one = ( 1.0e+0, 0.0e+0 ) )
210  REAL two
211  parameter( two = 2.0e+0 )
212  REAL three
213  parameter( three = 3.0e+0 )
214 * ..
215 * .. Local Scalars ..
216  LOGICAL upper
217  INTEGER count, i, j, k, kase, nz
218  REAL eps, lstres, s, safe1, safe2, safmin, xk
219  COMPLEX zdum
220 * ..
221 * .. Local Arrays ..
222  INTEGER isave( 3 )
223 * ..
224 * .. External Subroutines ..
225  EXTERNAL caxpy, ccopy, chemv, clacn2, cpotrs, xerbla
226 * ..
227 * .. Intrinsic Functions ..
228  INTRINSIC abs, aimag, max, real
229 * ..
230 * .. External Functions ..
231  LOGICAL lsame
232  REAL slamch
233  EXTERNAL lsame, slamch
234 * ..
235 * .. Statement Functions ..
236  REAL cabs1
237 * ..
238 * .. Statement Function definitions ..
239  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
240 * ..
241 * .. Executable Statements ..
242 *
243 * Test the input parameters.
244 *
245  info = 0
246  upper = lsame( uplo, 'U' )
247  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
248  info = -1
249  ELSE IF( n.LT.0 ) THEN
250  info = -2
251  ELSE IF( nrhs.LT.0 ) THEN
252  info = -3
253  ELSE IF( lda.LT.max( 1, n ) ) THEN
254  info = -5
255  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
256  info = -7
257  ELSE IF( ldb.LT.max( 1, n ) ) THEN
258  info = -9
259  ELSE IF( ldx.LT.max( 1, n ) ) THEN
260  info = -11
261  END IF
262  IF( info.NE.0 ) THEN
263  CALL xerbla( 'CPORFS', -info )
264  RETURN
265  END IF
266 *
267 * Quick return if possible
268 *
269  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
270  DO 10 j = 1, nrhs
271  ferr( j ) = zero
272  berr( j ) = zero
273  10 CONTINUE
274  RETURN
275  END IF
276 *
277 * NZ = maximum number of nonzero elements in each row of A, plus 1
278 *
279  nz = n + 1
280  eps = slamch( 'Epsilon' )
281  safmin = slamch( 'Safe minimum' )
282  safe1 = nz*safmin
283  safe2 = safe1 / eps
284 *
285 * Do for each right hand side
286 *
287  DO 140 j = 1, nrhs
288 *
289  count = 1
290  lstres = three
291  20 CONTINUE
292 *
293 * Loop until stopping criterion is satisfied.
294 *
295 * Compute residual R = B - A * X
296 *
297  CALL ccopy( n, b( 1, j ), 1, work, 1 )
298  CALL chemv( uplo, n, -one, a, lda, x( 1, j ), 1, one, work, 1 )
299 *
300 * Compute componentwise relative backward error from formula
301 *
302 * max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
303 *
304 * where abs(Z) is the componentwise absolute value of the matrix
305 * or vector Z. If the i-th component of the denominator is less
306 * than SAFE2, then SAFE1 is added to the i-th components of the
307 * numerator and denominator before dividing.
308 *
309  DO 30 i = 1, n
310  rwork( i ) = cabs1( b( i, j ) )
311  30 CONTINUE
312 *
313 * Compute abs(A)*abs(X) + abs(B).
314 *
315  IF( upper ) THEN
316  DO 50 k = 1, n
317  s = zero
318  xk = cabs1( x( k, j ) )
319  DO 40 i = 1, k - 1
320  rwork( i ) = rwork( i ) + cabs1( a( i, k ) )*xk
321  s = s + cabs1( a( i, k ) )*cabs1( x( i, j ) )
322  40 CONTINUE
323  rwork( k ) = rwork( k ) + abs( REAL( A( K, K ) ) )*xk + s
324  50 CONTINUE
325  ELSE
326  DO 70 k = 1, n
327  s = zero
328  xk = cabs1( x( k, j ) )
329  rwork( k ) = rwork( k ) + abs( REAL( A( K, K ) ) )*xk
330  DO 60 i = k + 1, n
331  rwork( i ) = rwork( i ) + cabs1( a( i, k ) )*xk
332  s = s + cabs1( a( i, k ) )*cabs1( x( i, j ) )
333  60 CONTINUE
334  rwork( k ) = rwork( k ) + s
335  70 CONTINUE
336  END IF
337  s = zero
338  DO 80 i = 1, n
339  IF( rwork( i ).GT.safe2 ) THEN
340  s = max( s, cabs1( work( i ) ) / rwork( i ) )
341  ELSE
342  s = max( s, ( cabs1( work( i ) )+safe1 ) /
343  $ ( rwork( i )+safe1 ) )
344  END IF
345  80 CONTINUE
346  berr( j ) = s
347 *
348 * Test stopping criterion. Continue iterating if
349 * 1) The residual BERR(J) is larger than machine epsilon, and
350 * 2) BERR(J) decreased by at least a factor of 2 during the
351 * last iteration, and
352 * 3) At most ITMAX iterations tried.
353 *
354  IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
355  $ count.LE.itmax ) THEN
356 *
357 * Update solution and try again.
358 *
359  CALL cpotrs( uplo, n, 1, af, ldaf, work, n, info )
360  CALL caxpy( n, one, work, 1, x( 1, j ), 1 )
361  lstres = berr( j )
362  count = count + 1
363  GO TO 20
364  END IF
365 *
366 * Bound error from formula
367 *
368 * norm(X - XTRUE) / norm(X) .le. FERR =
369 * norm( abs(inv(A))*
370 * ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
371 *
372 * where
373 * norm(Z) is the magnitude of the largest component of Z
374 * inv(A) is the inverse of A
375 * abs(Z) is the componentwise absolute value of the matrix or
376 * vector Z
377 * NZ is the maximum number of nonzeros in any row of A, plus 1
378 * EPS is machine epsilon
379 *
380 * The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B))
381 * is incremented by SAFE1 if the i-th component of
382 * abs(A)*abs(X) + abs(B) is less than SAFE2.
383 *
384 * Use CLACN2 to estimate the infinity-norm of the matrix
385 * inv(A) * diag(W),
386 * where W = abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) )))
387 *
388  DO 90 i = 1, n
389  IF( rwork( i ).GT.safe2 ) THEN
390  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i )
391  ELSE
392  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i ) +
393  $ safe1
394  END IF
395  90 CONTINUE
396 *
397  kase = 0
398  100 CONTINUE
399  CALL clacn2( n, work( n+1 ), work, ferr( j ), kase, isave )
400  IF( kase.NE.0 ) THEN
401  IF( kase.EQ.1 ) THEN
402 *
403 * Multiply by diag(W)*inv(A**H).
404 *
405  CALL cpotrs( uplo, n, 1, af, ldaf, work, n, info )
406  DO 110 i = 1, n
407  work( i ) = rwork( i )*work( i )
408  110 CONTINUE
409  ELSE IF( kase.EQ.2 ) THEN
410 *
411 * Multiply by inv(A)*diag(W).
412 *
413  DO 120 i = 1, n
414  work( i ) = rwork( i )*work( i )
415  120 CONTINUE
416  CALL cpotrs( uplo, n, 1, af, ldaf, work, n, info )
417  END IF
418  GO TO 100
419  END IF
420 *
421 * Normalize error.
422 *
423  lstres = zero
424  DO 130 i = 1, n
425  lstres = max( lstres, cabs1( x( i, j ) ) )
426  130 CONTINUE
427  IF( lstres.NE.zero )
428  $ ferr( j ) = ferr( j ) / lstres
429 *
430  140 CONTINUE
431 *
432  RETURN
433 *
434 * End of CPORFS
435 *
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine chemv(UPLO, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CHEMV
Definition: chemv.f:156
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine cpotrs(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
CPOTRS
Definition: cpotrs.f:112
subroutine caxpy(N, CA, CX, INCX, CY, INCY)
CAXPY
Definition: caxpy.f:53
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:52
subroutine clacn2(N, V, X, EST, KASE, ISAVE)
CLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: clacn2.f:135

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cporfsx ( character  UPLO,
character  EQUED,
integer  N,
integer  NRHS,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldaf, * )  AF,
integer  LDAF,
real, dimension( * )  S,
complex, dimension( ldb, * )  B,
integer  LDB,
complex, dimension( ldx, * )  X,
integer  LDX,
real  RCOND,
real, dimension( * )  BERR,
integer  N_ERR_BNDS,
real, dimension( nrhs, * )  ERR_BNDS_NORM,
real, dimension( nrhs, * )  ERR_BNDS_COMP,
integer  NPARAMS,
real, dimension(*)  PARAMS,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CPORFSX

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

Purpose:
    CPORFSX improves the computed solution to a system of linear
    equations when the coefficient matrix is symmetric positive
    definite, and provides error bounds and backward error estimates
    for the solution.  In addition to normwise error bound, the code
    provides maximum componentwise error bound if possible.  See
    comments for ERR_BNDS_NORM and ERR_BNDS_COMP for details of the
    error bounds.

    The original system of linear equations may have been equilibrated
    before calling this routine, as described by arguments EQUED and S
    below. In this case, the solution and error bounds returned are
    for the original unequilibrated system.
     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]UPLO
          UPLO is CHARACTER*1
       = 'U':  Upper triangle of A is stored;
       = 'L':  Lower triangle of A is stored.
[in]EQUED
          EQUED is CHARACTER*1
     Specifies the form of equilibration that was done to A
     before calling this routine. This is needed to compute
     the solution and error bounds correctly.
       = 'N':  No equilibration
       = 'Y':  Both row and column equilibration, i.e., A has been
               replaced by diag(S) * A * diag(S).
               The right hand side B has been changed accordingly.
[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 matrices B and X.  NRHS >= 0.
[in]A
          A is COMPLEX array, dimension (LDA,N)
     The symmetric matrix A.  If UPLO = 'U', the leading N-by-N
     upper triangular part of A contains the upper triangular part
     of the matrix A, and the strictly lower triangular part of A
     is not referenced.  If UPLO = 'L', the leading N-by-N lower
     triangular part of A contains the lower triangular part of
     the matrix A, and the strictly upper triangular part of A is
     not referenced.
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in]AF
          AF is COMPLEX array, dimension (LDAF,N)
     The triangular factor U or L from the Cholesky factorization
     A = U**T*U or A = L*L**T, as computed by SPOTRF.
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in,out]S
          S is REAL 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]B
          B is COMPLEX array, dimension (LDB,NRHS)
     The right hand side matrix B.
[in]LDB
          LDB is INTEGER
     The leading dimension of the array B.  LDB >= max(1,N).
[in,out]X
          X is COMPLEX array, dimension (LDX,NRHS)
     On entry, the solution matrix X, as computed by SGETRS.
     On exit, the improved solution matrix X.
[in]LDX
          LDX is INTEGER
     The leading dimension of the array X.  LDX >= max(1,N).
[out]RCOND
          RCOND is REAL
     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]BERR
          BERR is REAL 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 REAL 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 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" 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) * slamch('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) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". 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 REAL 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 .LT. 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 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * slamch('Epsilon').

     err = 2 "Guaranteed" 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) * slamch('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) * slamch('Epsilon') to determine if the error
              estimate is "guaranteed". 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 .LE. 0, the
     PARAMS array is never referenced and default values are used.
[in,out]PARAMS
          PARAMS is REAL array, dimension NPARAMS
     Specifies algorithm parameters.  If an entry is .LT. 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.0
            = 0.0 : No refinement is performed, and no error bounds are
                    computed.
            = 1.0 : Use the double-precision refinement algorithm,
                    possibly with doubled-single computations if the
                    compilation environment does not support DOUBLE
                    PRECISION.
              (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 array, dimension (2*N)
[out]RWORK
          RWORK is REAL 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.
Date
April 2012

Definition at line 395 of file cporfsx.f.

395 *
396 * -- LAPACK computational routine (version 3.4.1) --
397 * -- LAPACK is a software package provided by Univ. of Tennessee, --
398 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
399 * April 2012
400 *
401 * .. Scalar Arguments ..
402  CHARACTER uplo, equed
403  INTEGER info, lda, ldaf, ldb, ldx, n, nrhs, nparams,
404  $ n_err_bnds
405  REAL rcond
406 * ..
407 * .. Array Arguments ..
408  COMPLEX a( lda, * ), af( ldaf, * ), b( ldb, * ),
409  $ x( ldx, * ), work( * )
410  REAL rwork( * ), s( * ), params(*), berr( * ),
411  $ err_bnds_norm( nrhs, * ),
412  $ err_bnds_comp( nrhs, * )
413 * ..
414 *
415 * ==================================================================
416 *
417 * .. Parameters ..
418  REAL zero, one
419  parameter( zero = 0.0e+0, one = 1.0e+0 )
420  REAL itref_default, ithresh_default,
421  $ componentwise_default
422  REAL rthresh_default, dzthresh_default
423  parameter( itref_default = 1.0 )
424  parameter( ithresh_default = 10.0 )
425  parameter( componentwise_default = 1.0 )
426  parameter( rthresh_default = 0.5 )
427  parameter( dzthresh_default = 0.25 )
428  INTEGER la_linrx_itref_i, la_linrx_ithresh_i,
429  $ la_linrx_cwise_i
430  parameter( la_linrx_itref_i = 1,
431  $ la_linrx_ithresh_i = 2 )
432  parameter( la_linrx_cwise_i = 3 )
433  INTEGER la_linrx_trust_i, la_linrx_err_i,
434  $ la_linrx_rcond_i
435  parameter( la_linrx_trust_i = 1, la_linrx_err_i = 2 )
436  parameter( la_linrx_rcond_i = 3 )
437 * ..
438 * .. Local Scalars ..
439  CHARACTER(1) norm
440  LOGICAL rcequ
441  INTEGER j, prec_type, ref_type
442  INTEGER n_norms
443  REAL anorm, rcond_tmp
444  REAL illrcond_thresh, err_lbnd, cwise_wrong
445  LOGICAL ignore_cwise
446  INTEGER ithresh
447  REAL rthresh, unstable_thresh
448 * ..
449 * .. External Subroutines ..
451 * ..
452 * .. Intrinsic Functions ..
453  INTRINSIC max, sqrt, transfer
454 * ..
455 * .. External Functions ..
456  EXTERNAL lsame, blas_fpinfo_x, ilatrans, ilaprec
459  LOGICAL lsame
460  INTEGER blas_fpinfo_x
461  INTEGER ilatrans, ilaprec
462 * ..
463 * .. Executable Statements ..
464 *
465 * Check the input parameters.
466 *
467  info = 0
468  ref_type = int( itref_default )
469  IF ( nparams .GE. la_linrx_itref_i ) THEN
470  IF ( params( la_linrx_itref_i ) .LT. 0.0 ) THEN
471  params( la_linrx_itref_i ) = itref_default
472  ELSE
473  ref_type = params( la_linrx_itref_i )
474  END IF
475  END IF
476 *
477 * Set default parameters.
478 *
479  illrcond_thresh = REAL( N ) * slamch( 'Epsilon' )
480  ithresh = int( ithresh_default )
481  rthresh = rthresh_default
482  unstable_thresh = dzthresh_default
483  ignore_cwise = componentwise_default .EQ. 0.0
484 *
485  IF ( nparams.GE.la_linrx_ithresh_i ) THEN
486  IF ( params(la_linrx_ithresh_i ).LT.0.0 ) THEN
487  params( la_linrx_ithresh_i ) = ithresh
488  ELSE
489  ithresh = int( params( la_linrx_ithresh_i ) )
490  END IF
491  END IF
492  IF ( nparams.GE.la_linrx_cwise_i ) THEN
493  IF ( params(la_linrx_cwise_i ).LT.0.0 ) THEN
494  IF ( ignore_cwise ) THEN
495  params( la_linrx_cwise_i ) = 0.0
496  ELSE
497  params( la_linrx_cwise_i ) = 1.0
498  END IF
499  ELSE
500  ignore_cwise = params( la_linrx_cwise_i ) .EQ. 0.0
501  END IF
502  END IF
503  IF ( ref_type .EQ. 0 .OR. n_err_bnds .EQ. 0 ) THEN
504  n_norms = 0
505  ELSE IF ( ignore_cwise ) THEN
506  n_norms = 1
507  ELSE
508  n_norms = 2
509  END IF
510 *
511  rcequ = lsame( equed, 'Y' )
512 *
513 * Test input parameters.
514 *
515  IF (.NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
516  info = -1
517  ELSE IF( .NOT.rcequ .AND. .NOT.lsame( equed, 'N' ) ) THEN
518  info = -2
519  ELSE IF( n.LT.0 ) THEN
520  info = -3
521  ELSE IF( nrhs.LT.0 ) THEN
522  info = -4
523  ELSE IF( lda.LT.max( 1, n ) ) THEN
524  info = -6
525  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
526  info = -8
527  ELSE IF( ldb.LT.max( 1, n ) ) THEN
528  info = -11
529  ELSE IF( ldx.LT.max( 1, n ) ) THEN
530  info = -13
531  END IF
532  IF( info.NE.0 ) THEN
533  CALL xerbla( 'CPORFSX', -info )
534  RETURN
535  END IF
536 *
537 * Quick return if possible.
538 *
539  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
540  rcond = 1.0
541  DO j = 1, nrhs
542  berr( j ) = 0.0
543  IF ( n_err_bnds .GE. 1 ) THEN
544  err_bnds_norm( j, la_linrx_trust_i ) = 1.0
545  err_bnds_comp( j, la_linrx_trust_i ) = 1.0
546  END IF
547  IF ( n_err_bnds .GE. 2 ) THEN
548  err_bnds_norm( j, la_linrx_err_i ) = 0.0
549  err_bnds_comp( j, la_linrx_err_i ) = 0.0
550  END IF
551  IF ( n_err_bnds .GE. 3 ) THEN
552  err_bnds_norm( j, la_linrx_rcond_i ) = 1.0
553  err_bnds_comp( j, la_linrx_rcond_i ) = 1.0
554  END IF
555  END DO
556  RETURN
557  END IF
558 *
559 * Default to failure.
560 *
561  rcond = 0.0
562  DO j = 1, nrhs
563  berr( j ) = 1.0
564  IF ( n_err_bnds .GE. 1 ) THEN
565  err_bnds_norm( j, la_linrx_trust_i ) = 1.0
566  err_bnds_comp( j, la_linrx_trust_i ) = 1.0
567  END IF
568  IF ( n_err_bnds .GE. 2 ) THEN
569  err_bnds_norm( j, la_linrx_err_i ) = 1.0
570  err_bnds_comp( j, la_linrx_err_i ) = 1.0
571  END IF
572  IF ( n_err_bnds .GE. 3 ) THEN
573  err_bnds_norm( j, la_linrx_rcond_i ) = 0.0
574  err_bnds_comp( j, la_linrx_rcond_i ) = 0.0
575  END IF
576  END DO
577 *
578 * Compute the norm of A and the reciprocal of the condition
579 * number of A.
580 *
581  norm = 'I'
582  anorm = clanhe( norm, uplo, n, a, lda, rwork )
583  CALL cpocon( uplo, n, af, ldaf, anorm, rcond, work, rwork,
584  $ info )
585 *
586 * Perform refinement on each right-hand side
587 *
588  IF ( ref_type .NE. 0 ) THEN
589 
590  prec_type = ilaprec( 'D' )
591 
592  CALL cla_porfsx_extended( prec_type, uplo, n,
593  $ nrhs, a, lda, af, ldaf, rcequ, s, b,
594  $ ldb, x, ldx, berr, n_norms, err_bnds_norm, err_bnds_comp,
595  $ work, rwork, work(n+1),
596  $ transfer(rwork(1:2*n), (/ (zero, zero) /), n), rcond,
597  $ ithresh, rthresh, unstable_thresh, ignore_cwise,
598  $ info )
599  END IF
600 
601  err_lbnd = max( 10.0, sqrt( REAL( N ) ) ) * slamch( 'Epsilon' )
602  IF ( n_err_bnds .GE. 1 .AND. n_norms .GE. 1 ) THEN
603 *
604 * Compute scaled normwise condition number cond(A*C).
605 *
606  IF ( rcequ ) THEN
607  rcond_tmp = cla_porcond_c( uplo, n, a, lda, af, ldaf,
608  $ s, .true., info, work, rwork )
609  ELSE
610  rcond_tmp = cla_porcond_c( uplo, n, a, lda, af, ldaf,
611  $ s, .false., info, work, rwork )
612  END IF
613  DO j = 1, nrhs
614 *
615 * Cap the error at 1.0.
616 *
617  IF ( n_err_bnds .GE. la_linrx_err_i
618  $ .AND. err_bnds_norm( j, la_linrx_err_i ) .GT. 1.0 )
619  $ err_bnds_norm( j, la_linrx_err_i ) = 1.0
620 *
621 * Threshold the error (see LAWN).
622 *
623  IF ( rcond_tmp .LT. illrcond_thresh ) THEN
624  err_bnds_norm( j, la_linrx_err_i ) = 1.0
625  err_bnds_norm( j, la_linrx_trust_i ) = 0.0
626  IF ( info .LE. n ) info = n + j
627  ELSE IF ( err_bnds_norm( j, la_linrx_err_i ) .LT. err_lbnd )
628  $ THEN
629  err_bnds_norm( j, la_linrx_err_i ) = err_lbnd
630  err_bnds_norm( j, la_linrx_trust_i ) = 1.0
631  END IF
632 *
633 * Save the condition number.
634 *
635  IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
636  err_bnds_norm( j, la_linrx_rcond_i ) = rcond_tmp
637  END IF
638 
639  END DO
640  END IF
641 
642  IF (n_err_bnds .GE. 1 .AND. n_norms .GE. 2) THEN
643 *
644 * Compute componentwise condition number cond(A*diag(Y(:,J))) for
645 * each right-hand side using the current solution as an estimate of
646 * the true solution. If the componentwise error estimate is too
647 * large, then the solution is a lousy estimate of truth and the
648 * estimated RCOND may be too optimistic. To avoid misleading users,
649 * the inverse condition number is set to 0.0 when the estimated
650 * cwise error is at least CWISE_WRONG.
651 *
652  cwise_wrong = sqrt( slamch( 'Epsilon' ) )
653  DO j = 1, nrhs
654  IF (err_bnds_comp( j, la_linrx_err_i ) .LT. cwise_wrong )
655  $ THEN
656  rcond_tmp = cla_porcond_x( uplo, n, a, lda, af, ldaf,
657  $ x(1,j), info, work, rwork )
658  ELSE
659  rcond_tmp = 0.0
660  END IF
661 *
662 * Cap the error at 1.0.
663 *
664  IF ( n_err_bnds .GE. la_linrx_err_i
665  $ .AND. err_bnds_comp( j, la_linrx_err_i ) .GT. 1.0 )
666  $ err_bnds_comp( j, la_linrx_err_i ) = 1.0
667 *
668 * Threshold the error (see LAWN).
669 *
670  IF (rcond_tmp .LT. illrcond_thresh) THEN
671  err_bnds_comp( j, la_linrx_err_i ) = 1.0
672  err_bnds_comp( j, la_linrx_trust_i ) = 0.0
673  IF ( params( la_linrx_cwise_i ) .EQ. 1.0
674  $ .AND. info.LT.n + j ) info = n + j
675  ELSE IF ( err_bnds_comp( j, la_linrx_err_i )
676  $ .LT. err_lbnd ) THEN
677  err_bnds_comp( j, la_linrx_err_i ) = err_lbnd
678  err_bnds_comp( j, la_linrx_trust_i ) = 1.0
679  END IF
680 *
681 * Save the condition number.
682 *
683  IF ( n_err_bnds .GE. la_linrx_rcond_i ) THEN
684  err_bnds_comp( j, la_linrx_rcond_i ) = rcond_tmp
685  END IF
686 
687  END DO
688  END IF
689 *
690  RETURN
691 *
692 * End of CPORFSX
693 *
integer function ilatrans(TRANS)
ILATRANS
Definition: ilatrans.f:60
real function clanhe(NORM, UPLO, N, A, LDA, WORK)
CLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix.
Definition: clanhe.f:126
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function cla_porcond_x(UPLO, N, A, LDA, AF, LDAF, X, INFO, WORK, RWORK)
CLA_PORCOND_X computes the infinity norm condition number of op(A)*diag(x) for Hermitian positive-def...
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
integer function ilaprec(PREC)
ILAPREC
Definition: ilaprec.f:60
subroutine cla_porfsx_extended(PREC_TYPE, UPLO, N, NRHS, A, LDA, AF, LDAF, COLEQU, C, B, LDB, Y, LDY, BERR_OUT, N_NORMS, ERR_BNDS_NORM, ERR_BNDS_COMP, RES, AYB, DY, Y_TAIL, RCOND, ITHRESH, RTHRESH, DZ_UB, IGNORE_CWISE, INFO)
CLA_PORFSX_EXTENDED improves the computed solution to a system of linear equations for symmetric or H...
real function cla_porcond_c(UPLO, N, A, LDA, AF, LDAF, C, CAPPLY, INFO, WORK, RWORK)
CLA_PORCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for Hermitian positiv...
subroutine cpocon(UPLO, N, A, LDA, ANORM, RCOND, WORK, RWORK, INFO)
CPOCON
Definition: cpocon.f:123

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cpotf2 ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
integer  INFO 
)

CPOTF2 computes the Cholesky factorization of a symmetric/Hermitian positive definite matrix (unblocked algorithm).

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

Purpose:
 CPOTF2 computes the Cholesky factorization of a complex Hermitian
 positive definite matrix A.

 The factorization has the form
    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 lower triangular.

 This is the unblocked version of the algorithm, calling Level 2 BLAS.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the upper or lower triangular part of the
          Hermitian matrix A is stored.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is COMPLEX 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).
[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 leading minor of order k is not
               positive definite, and the factorization could not be
               completed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 111 of file cpotf2.f.

111 *
112 * -- LAPACK computational routine (version 3.4.2) --
113 * -- LAPACK is a software package provided by Univ. of Tennessee, --
114 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
115 * September 2012
116 *
117 * .. Scalar Arguments ..
118  CHARACTER uplo
119  INTEGER info, lda, n
120 * ..
121 * .. Array Arguments ..
122  COMPLEX a( lda, * )
123 * ..
124 *
125 * =====================================================================
126 *
127 * .. Parameters ..
128  REAL one, zero
129  parameter( one = 1.0e+0, zero = 0.0e+0 )
130  COMPLEX cone
131  parameter( cone = ( 1.0e+0, 0.0e+0 ) )
132 * ..
133 * .. Local Scalars ..
134  LOGICAL upper
135  INTEGER j
136  REAL ajj
137 * ..
138 * .. External Functions ..
139  LOGICAL lsame, sisnan
140  COMPLEX cdotc
141  EXTERNAL lsame, cdotc, sisnan
142 * ..
143 * .. External Subroutines ..
144  EXTERNAL cgemv, clacgv, csscal, xerbla
145 * ..
146 * .. Intrinsic Functions ..
147  INTRINSIC max, REAL, sqrt
148 * ..
149 * .. Executable Statements ..
150 *
151 * Test the input parameters.
152 *
153  info = 0
154  upper = lsame( uplo, 'U' )
155  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
156  info = -1
157  ELSE IF( n.LT.0 ) THEN
158  info = -2
159  ELSE IF( lda.LT.max( 1, n ) ) THEN
160  info = -4
161  END IF
162  IF( info.NE.0 ) THEN
163  CALL xerbla( 'CPOTF2', -info )
164  RETURN
165  END IF
166 *
167 * Quick return if possible
168 *
169  IF( n.EQ.0 )
170  $ RETURN
171 *
172  IF( upper ) THEN
173 *
174 * Compute the Cholesky factorization A = U**H *U.
175 *
176  DO 10 j = 1, n
177 *
178 * Compute U(J,J) and test for non-positive-definiteness.
179 *
180  ajj = REAL( A( J, J ) ) - cdotc( j-1, a( 1, j ), 1,
181  $ a( 1, j ), 1 )
182  IF( ajj.LE.zero.OR.sisnan( ajj ) ) THEN
183  a( j, j ) = ajj
184  GO TO 30
185  END IF
186  ajj = sqrt( ajj )
187  a( j, j ) = ajj
188 *
189 * Compute elements J+1:N of row J.
190 *
191  IF( j.LT.n ) THEN
192  CALL clacgv( j-1, a( 1, j ), 1 )
193  CALL cgemv( 'Transpose', j-1, n-j, -cone, a( 1, j+1 ),
194  $ lda, a( 1, j ), 1, cone, a( j, j+1 ), lda )
195  CALL clacgv( j-1, a( 1, j ), 1 )
196  CALL csscal( n-j, one / ajj, a( j, j+1 ), lda )
197  END IF
198  10 CONTINUE
199  ELSE
200 *
201 * Compute the Cholesky factorization A = L*L**H.
202 *
203  DO 20 j = 1, n
204 *
205 * Compute L(J,J) and test for non-positive-definiteness.
206 *
207  ajj = REAL( A( J, J ) ) - cdotc( j-1, a( j, 1 ), lda,
208  $ a( j, 1 ), lda )
209  IF( ajj.LE.zero.OR.sisnan( ajj ) ) THEN
210  a( j, j ) = ajj
211  GO TO 30
212  END IF
213  ajj = sqrt( ajj )
214  a( j, j ) = ajj
215 *
216 * Compute elements J+1:N of column J.
217 *
218  IF( j.LT.n ) THEN
219  CALL clacgv( j-1, a( j, 1 ), lda )
220  CALL cgemv( 'No transpose', n-j, j-1, -cone, a( j+1, 1 ),
221  $ lda, a( j, 1 ), lda, cone, a( j+1, j ), 1 )
222  CALL clacgv( j-1, a( j, 1 ), lda )
223  CALL csscal( n-j, one / ajj, a( j+1, j ), 1 )
224  END IF
225  20 CONTINUE
226  END IF
227  GO TO 40
228 *
229  30 CONTINUE
230  info = j
231 *
232  40 CONTINUE
233  RETURN
234 *
235 * End of CPOTF2
236 *
logical function sisnan(SIN)
SISNAN tests input for NaN.
Definition: sisnan.f:61
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
CGEMV
Definition: cgemv.f:160
subroutine csscal(N, SA, CX, INCX)
CSSCAL
Definition: csscal.f:54
subroutine clacgv(N, X, INCX)
CLACGV conjugates a complex vector.
Definition: clacgv.f:76
complex function cdotc(N, CX, INCX, CY, INCY)
CDOTC
Definition: cdotc.f:54
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cpotrf ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
integer  INFO 
)

CPOTRF

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

Purpose:
 CPOTRF computes the Cholesky factorization of a complex Hermitian
 positive definite matrix A.

 The factorization has the form
    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 lower triangular.

 This is the block 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,out]A
          A is COMPLEX 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).
[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 is not
                positive definite, and the factorization could not be
                completed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 109 of file cpotrf.f.

109 *
110 * -- LAPACK computational routine (version 3.6.0) --
111 * -- LAPACK is a software package provided by Univ. of Tennessee, --
112 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
113 * November 2015
114 *
115 * .. Scalar Arguments ..
116  CHARACTER uplo
117  INTEGER info, lda, n
118 * ..
119 * .. Array Arguments ..
120  COMPLEX a( lda, * )
121 * ..
122 *
123 * =====================================================================
124 *
125 * .. Parameters ..
126  REAL one
127  COMPLEX cone
128  parameter( one = 1.0e+0, cone = ( 1.0e+0, 0.0e+0 ) )
129 * ..
130 * .. Local Scalars ..
131  LOGICAL upper
132  INTEGER j, jb, nb
133 * ..
134 * .. External Functions ..
135  LOGICAL lsame
136  INTEGER ilaenv
137  EXTERNAL lsame, ilaenv
138 * ..
139 * .. External Subroutines ..
140  EXTERNAL cgemm, cherk, cpotrf2, ctrsm, xerbla
141 * ..
142 * .. Intrinsic Functions ..
143  INTRINSIC max, min
144 * ..
145 * .. Executable Statements ..
146 *
147 * Test the input parameters.
148 *
149  info = 0
150  upper = lsame( uplo, 'U' )
151  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
152  info = -1
153  ELSE IF( n.LT.0 ) THEN
154  info = -2
155  ELSE IF( lda.LT.max( 1, n ) ) THEN
156  info = -4
157  END IF
158  IF( info.NE.0 ) THEN
159  CALL xerbla( 'CPOTRF', -info )
160  RETURN
161  END IF
162 *
163 * Quick return if possible
164 *
165  IF( n.EQ.0 )
166  $ RETURN
167 *
168 * Determine the block size for this environment.
169 *
170  nb = ilaenv( 1, 'CPOTRF', uplo, n, -1, -1, -1 )
171  IF( nb.LE.1 .OR. nb.GE.n ) THEN
172 *
173 * Use unblocked code.
174 *
175  CALL cpotrf2( uplo, n, a, lda, info )
176  ELSE
177 *
178 * Use blocked code.
179 *
180  IF( upper ) THEN
181 *
182 * Compute the Cholesky factorization A = U**H *U.
183 *
184  DO 10 j = 1, n, nb
185 *
186 * Update and factorize the current diagonal block and test
187 * for non-positive-definiteness.
188 *
189  jb = min( nb, n-j+1 )
190  CALL cherk( 'Upper', 'Conjugate transpose', jb, j-1,
191  $ -one, a( 1, j ), lda, one, a( j, j ), lda )
192  CALL cpotrf2( 'Upper', jb, a( j, j ), lda, info )
193  IF( info.NE.0 )
194  $ GO TO 30
195  IF( j+jb.LE.n ) THEN
196 *
197 * Compute the current block row.
198 *
199  CALL cgemm( 'Conjugate transpose', 'No transpose', jb,
200  $ n-j-jb+1, j-1, -cone, a( 1, j ), lda,
201  $ a( 1, j+jb ), lda, cone, a( j, j+jb ),
202  $ lda )
203  CALL ctrsm( 'Left', 'Upper', 'Conjugate transpose',
204  $ 'Non-unit', jb, n-j-jb+1, cone, a( j, j ),
205  $ lda, a( j, j+jb ), lda )
206  END IF
207  10 CONTINUE
208 *
209  ELSE
210 *
211 * Compute the Cholesky factorization A = L*L**H.
212 *
213  DO 20 j = 1, n, nb
214 *
215 * Update and factorize the current diagonal block and test
216 * for non-positive-definiteness.
217 *
218  jb = min( nb, n-j+1 )
219  CALL cherk( 'Lower', 'No transpose', jb, j-1, -one,
220  $ a( j, 1 ), lda, one, a( j, j ), lda )
221  CALL cpotrf2( 'Lower', jb, a( j, j ), lda, info )
222  IF( info.NE.0 )
223  $ GO TO 30
224  IF( j+jb.LE.n ) THEN
225 *
226 * Compute the current block column.
227 *
228  CALL cgemm( 'No transpose', 'Conjugate transpose',
229  $ n-j-jb+1, jb, j-1, -cone, a( j+jb, 1 ),
230  $ lda, a( j, 1 ), lda, cone, a( j+jb, j ),
231  $ lda )
232  CALL ctrsm( 'Right', 'Lower', 'Conjugate transpose',
233  $ 'Non-unit', n-j-jb+1, jb, cone, a( j, j ),
234  $ lda, a( j+jb, j ), lda )
235  END IF
236  20 CONTINUE
237  END IF
238  END IF
239  GO TO 40
240 *
241  30 CONTINUE
242  info = info + j - 1
243 *
244  40 CONTINUE
245  RETURN
246 *
247 * End of CPOTRF
248 *
subroutine cherk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
CHERK
Definition: cherk.f:175
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:189
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:182
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
recursive subroutine cpotrf2(UPLO, N, A, LDA, INFO)
CPOTRF2
Definition: cpotrf2.f:108

Here is the call graph for this function:

Here is the caller graph for this function:

recursive subroutine cpotrf2 ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
integer  INFO 
)

CPOTRF2

Purpose:
 CPOTRF2 computes the Cholesky factorization of a real symmetric
 positive definite matrix A using the recursive algorithm.

 The factorization has the form
    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 lower triangular.

 This is the recursive version of the algorithm. It divides
 the matrix into four submatrices:

        [  A11 | A12  ]  where A11 is n1 by n1 and A22 is n2 by n2
    A = [ -----|----- ]  with n1 = n/2
        [  A21 | A22  ]       n2 = n-n1

 The subroutine calls itself to factor A11. Update and scale A21
 or A12, update A22 then calls itself to factor A22.
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,out]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the symmetric matrix A.  If UPLO = 'U', the leading
          N-by-N upper triangular part of A contains the upper
          triangular part of the matrix A, and the strictly lower
          triangular part of A is not referenced.  If UPLO = 'L', the
          leading N-by-N lower triangular part of A contains the lower
          triangular part of the matrix A, and the strictly upper
          triangular part of A is not referenced.

          On exit, if INFO = 0, the factor U or L from the Cholesky
          factorization A = U**H*U or A = L*L**H.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= 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 is not
                positive definite, and the factorization could not be
                completed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 108 of file cpotrf2.f.

108 *
109 * -- LAPACK computational routine (version 3.6.0) --
110 * -- LAPACK is a software package provided by Univ. of Tennessee, --
111 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
112 * November 2015
113 *
114 * .. Scalar Arguments ..
115  CHARACTER uplo
116  INTEGER info, lda, n
117 * ..
118 * .. Array Arguments ..
119  COMPLEX a( lda, * )
120 * ..
121 *
122 * =====================================================================
123 *
124 * .. Parameters ..
125  REAL one, zero
126  parameter( one = 1.0e+0, zero = 0.0e+0 )
127  COMPLEX cone
128  parameter( cone = (1.0e+0, 0.0e+0) )
129 * ..
130 * .. Local Scalars ..
131  LOGICAL upper
132  INTEGER n1, n2, iinfo
133  REAL ajj
134 * ..
135 * .. External Functions ..
136  LOGICAL lsame, sisnan
137  EXTERNAL lsame, sisnan
138 * ..
139 * .. External Subroutines ..
140  EXTERNAL cherk, ctrsm, xerbla
141 * ..
142 * .. Intrinsic Functions ..
143  INTRINSIC max, REAL, sqrt
144 * ..
145 * .. Executable Statements ..
146 *
147 * Test the input parameters
148 *
149  info = 0
150  upper = lsame( uplo, 'U' )
151  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
152  info = -1
153  ELSE IF( n.LT.0 ) THEN
154  info = -2
155  ELSE IF( lda.LT.max( 1, n ) ) THEN
156  info = -4
157  END IF
158  IF( info.NE.0 ) THEN
159  CALL xerbla( 'CPOTRF2', -info )
160  RETURN
161  END IF
162 *
163 * Quick return if possible
164 *
165  IF( n.EQ.0 )
166  $ RETURN
167 *
168 * N=1 case
169 *
170  IF( n.EQ.1 ) THEN
171 *
172 * Test for non-positive-definiteness
173 *
174  ajj = REAL( A( 1, 1 ) )
175  IF( ajj.LE.zero.OR.sisnan( ajj ) ) THEN
176  info = 1
177  RETURN
178  END IF
179 *
180 * Factor
181 *
182  a( 1, 1 ) = sqrt( ajj )
183 *
184 * Use recursive code
185 *
186  ELSE
187  n1 = n/2
188  n2 = n-n1
189 *
190 * Factor A11
191 *
192  CALL cpotrf2( uplo, n1, a( 1, 1 ), lda, iinfo )
193  IF ( iinfo.NE.0 ) THEN
194  info = iinfo
195  RETURN
196  END IF
197 *
198 * Compute the Cholesky factorization A = U**H*U
199 *
200  IF( upper ) THEN
201 *
202 * Update and scale A12
203 *
204  CALL ctrsm( 'L', 'U', 'C', 'N', n1, n2, cone,
205  $ a( 1, 1 ), lda, a( 1, n1+1 ), lda )
206 *
207 * Update and factor A22
208 *
209  CALL cherk( uplo, 'C', n2, n1, -one, a( 1, n1+1 ), lda,
210  $ one, a( n1+1, n1+1 ), lda )
211 *
212  CALL cpotrf2( uplo, n2, a( n1+1, n1+1 ), lda, iinfo )
213 *
214  IF ( iinfo.NE.0 ) THEN
215  info = iinfo + n1
216  RETURN
217  END IF
218 *
219 * Compute the Cholesky factorization A = L*L**H
220 *
221  ELSE
222 *
223 * Update and scale A21
224 *
225  CALL ctrsm( 'R', 'L', 'C', 'N', n2, n1, cone,
226  $ a( 1, 1 ), lda, a( n1+1, 1 ), lda )
227 *
228 * Update and factor A22
229 *
230  CALL cherk( uplo, 'N', n2, n1, -one, a( n1+1, 1 ), lda,
231  $ one, a( n1+1, n1+1 ), lda )
232 *
233  CALL cpotrf2( uplo, n2, a( n1+1, n1+1 ), lda, iinfo )
234 *
235  IF ( iinfo.NE.0 ) THEN
236  info = iinfo + n1
237  RETURN
238  END IF
239 *
240  END IF
241  END IF
242  RETURN
243 *
244 * End of CPOTRF2
245 *
subroutine cherk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
CHERK
Definition: cherk.f:175
logical function sisnan(SIN)
SISNAN tests input for NaN.
Definition: sisnan.f:61
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:182
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
recursive subroutine cpotrf2(UPLO, N, A, LDA, INFO)
CPOTRF2
Definition: cpotrf2.f:108

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cpotri ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
integer  INFO 
)

CPOTRI

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

Purpose:
 CPOTRI computes the inverse of a complex Hermitian positive definite
 matrix A using the Cholesky factorization A = U**H*U or A = L*L**H
 computed by CPOTRF.
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,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the triangular factor U or L from the Cholesky
          factorization A = U**H*U or A = L*L**H, as computed by
          CPOTRF.
          On exit, the upper or lower triangle of the (Hermitian)
          inverse of A, overwriting the input factor U or L.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= 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 (i,i) element of the factor U or L is
                zero, and the inverse could not be computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 97 of file cpotri.f.

97 *
98 * -- LAPACK computational routine (version 3.4.0) --
99 * -- LAPACK is a software package provided by Univ. of Tennessee, --
100 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
101 * November 2011
102 *
103 * .. Scalar Arguments ..
104  CHARACTER uplo
105  INTEGER info, lda, n
106 * ..
107 * .. Array Arguments ..
108  COMPLEX a( lda, * )
109 * ..
110 *
111 * =====================================================================
112 *
113 * .. External Functions ..
114  LOGICAL lsame
115  EXTERNAL lsame
116 * ..
117 * .. External Subroutines ..
118  EXTERNAL clauum, ctrtri, xerbla
119 * ..
120 * .. Intrinsic Functions ..
121  INTRINSIC max
122 * ..
123 * .. Executable Statements ..
124 *
125 * Test the input parameters.
126 *
127  info = 0
128  IF( .NOT.lsame( uplo, 'U' ) .AND. .NOT.lsame( uplo, 'L' ) ) THEN
129  info = -1
130  ELSE IF( n.LT.0 ) THEN
131  info = -2
132  ELSE IF( lda.LT.max( 1, n ) ) THEN
133  info = -4
134  END IF
135  IF( info.NE.0 ) THEN
136  CALL xerbla( 'CPOTRI', -info )
137  RETURN
138  END IF
139 *
140 * Quick return if possible
141 *
142  IF( n.EQ.0 )
143  $ RETURN
144 *
145 * Invert the triangular Cholesky factor U or L.
146 *
147  CALL ctrtri( uplo, 'Non-unit', n, a, lda, info )
148  IF( info.GT.0 )
149  $ RETURN
150 *
151 * Form inv(U) * inv(U)**H or inv(L)**H * inv(L).
152 *
153  CALL clauum( uplo, n, a, lda, info )
154 *
155  RETURN
156 *
157 * End of CPOTRI
158 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ctrtri(UPLO, DIAG, N, A, LDA, INFO)
CTRTRI
Definition: ctrtri.f:111
subroutine clauum(UPLO, N, A, LDA, INFO)
CLAUUM computes the product UUH or LHL, where U and L are upper or lower triangular matrices (blocked...
Definition: clauum.f:104
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine cpotrs ( character  UPLO,
integer  N,
integer  NRHS,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldb, * )  B,
integer  LDB,
integer  INFO 
)

CPOTRS

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

Purpose:
 CPOTRS solves a system of linear equations A*X = B with a Hermitian
 positive definite matrix A using the Cholesky factorization 
 A = U**H*U or A = L*L**H computed by CPOTRF.
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]A
          A is COMPLEX array, dimension (LDA,N)
          The triangular factor U or L from the Cholesky factorization
          A = U**H*U or A = L*L**H, as computed by CPOTRF.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is COMPLEX 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]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 112 of file cpotrs.f.

112 *
113 * -- LAPACK computational routine (version 3.4.0) --
114 * -- LAPACK is a software package provided by Univ. of Tennessee, --
115 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
116 * November 2011
117 *
118 * .. Scalar Arguments ..
119  CHARACTER uplo
120  INTEGER info, lda, ldb, n, nrhs
121 * ..
122 * .. Array Arguments ..
123  COMPLEX a( lda, * ), b( ldb, * )
124 * ..
125 *
126 * =====================================================================
127 *
128 * .. Parameters ..
129  COMPLEX one
130  parameter( one = ( 1.0e+0, 0.0e+0 ) )
131 * ..
132 * .. Local Scalars ..
133  LOGICAL upper
134 * ..
135 * .. External Functions ..
136  LOGICAL lsame
137  EXTERNAL lsame
138 * ..
139 * .. External Subroutines ..
140  EXTERNAL ctrsm, xerbla
141 * ..
142 * .. Intrinsic Functions ..
143  INTRINSIC max
144 * ..
145 * .. Executable Statements ..
146 *
147 * Test the input parameters.
148 *
149  info = 0
150  upper = lsame( uplo, 'U' )
151  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
152  info = -1
153  ELSE IF( n.LT.0 ) THEN
154  info = -2
155  ELSE IF( nrhs.LT.0 ) THEN
156  info = -3
157  ELSE IF( lda.LT.max( 1, n ) ) THEN
158  info = -5
159  ELSE IF( ldb.LT.max( 1, n ) ) THEN
160  info = -7
161  END IF
162  IF( info.NE.0 ) THEN
163  CALL xerbla( 'CPOTRS', -info )
164  RETURN
165  END IF
166 *
167 * Quick return if possible
168 *
169  IF( n.EQ.0 .OR. nrhs.EQ.0 )
170  $ RETURN
171 *
172  IF( upper ) THEN
173 *
174 * Solve A*X = B where A = U**H *U.
175 *
176 * Solve U**H *X = B, overwriting B with X.
177 *
178  CALL ctrsm( 'Left', 'Upper', 'Conjugate transpose', 'Non-unit',
179  $ n, nrhs, one, a, lda, b, ldb )
180 *
181 * Solve U*X = B, overwriting B with X.
182 *
183  CALL ctrsm( 'Left', 'Upper', 'No transpose', 'Non-unit', n,
184  $ nrhs, one, a, lda, b, ldb )
185  ELSE
186 *
187 * Solve A*X = B where A = L*L**H.
188 *
189 * Solve L*X = B, overwriting B with X.
190 *
191  CALL ctrsm( 'Left', 'Lower', 'No transpose', 'Non-unit', n,
192  $ nrhs, one, a, lda, b, ldb )
193 *
194 * Solve L**H *X = B, overwriting B with X.
195 *
196  CALL ctrsm( 'Left', 'Lower', 'Conjugate transpose', 'Non-unit',
197  $ n, nrhs, one, a, lda, b, ldb )
198  END IF
199 *
200  RETURN
201 *
202 * End of CPOTRS
203 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:182
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: