LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine csycon_3 ( character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( * )  E,
integer, dimension( * )  IPIV,
real  ANORM,
real  RCOND,
complex, dimension( * )  WORK,
integer  INFO 
)

CSYCON_3

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

Purpose:
 CSYCON_3 estimates the reciprocal of the condition number (in the
 1-norm) of a complex symmetric matrix A using the factorization
 computed by CSYTRF_RK or CSYTRF_BK:

    A = P*U*D*(U**T)*(P**T) or A = P*L*D*(L**T)*(P**T),

 where U (or L) is unit upper (or lower) triangular matrix,
 U**T (or L**T) is the transpose of U (or L), P is a permutation
 matrix, P**T is the transpose of P, and D is symmetric and block
 diagonal with 1-by-1 and 2-by-2 diagonal blocks.

 An estimate is obtained for norm(inv(A)), and the reciprocal of the
 condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
 This routine uses BLAS3 solver CSYTRS_3.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the details of the factorization are
          stored as an upper or lower triangular matrix:
          = 'U':  Upper triangular, form is A = P*U*D*(U**T)*(P**T);
          = 'L':  Lower triangular, form is A = P*L*D*(L**T)*(P**T).
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]A
          A is COMPLEX array, dimension (LDA,N)
          Diagonal of the block diagonal matrix D and factors U or L
          as computed by CSYTRF_RK and CSYTRF_BK:
            a) ONLY diagonal elements of the symmetric block diagonal
               matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
               (superdiagonal (or subdiagonal) elements of D
                should be provided on entry in array E), and
            b) If UPLO = 'U': factor U in the superdiagonal part of A.
               If UPLO = 'L': factor L in the subdiagonal part of A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in]E
          E is COMPLEX array, dimension (N)
          On entry, contains the superdiagonal (or subdiagonal)
          elements of the symmetric block diagonal matrix D
          with 1-by-1 or 2-by-2 diagonal blocks, where
          If UPLO = 'U': E(i) = D(i-1,i),i=2:N, E(1) not refernced;
          If UPLO = 'L': E(i) = D(i+1,i),i=1:N-1, E(N) not referenced.

          NOTE: For 1-by-1 diagonal block D(k), where
          1 <= k <= N, the element E(k) is not referenced in both
          UPLO = 'U' or UPLO = 'L' cases.
[in]IPIV
          IPIV is INTEGER array, dimension (N)
          Details of the interchanges and the block structure of D
          as determined by CSYTRF_RK or CSYTRF_BK.
[in]ANORM
          ANORM is REAL
          The 1-norm of the original 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]IWORK
          IWORK is INTEGER array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
December 2016
Contributors:

December 2016, Igor Kozachenko, Computer Science Division, University of California, Berkeley

September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas, School of Mathematics, University of Manchester

Definition at line 173 of file csycon_3.f.

173 *
174 * -- LAPACK computational routine (version 3.7.0) --
175 * -- LAPACK is a software package provided by Univ. of Tennessee, --
176 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
177 * December 2016
178 *
179 * .. Scalar Arguments ..
180  CHARACTER uplo
181  INTEGER info, lda, n
182  REAL anorm, rcond
183 * ..
184 * .. Array Arguments ..
185  INTEGER ipiv( * )
186  COMPLEX a( lda, * ), e( * ), work( * )
187 * ..
188 *
189 * =====================================================================
190 *
191 * .. Parameters ..
192  REAL one, zero
193  parameter ( one = 1.0e+0, zero = 0.0e+0 )
194  COMPLEX czero
195  parameter ( czero = ( 0.0e+0, 0.0e+0 ) )
196 * ..
197 * .. Local Scalars ..
198  LOGICAL upper
199  INTEGER i, kase
200  REAL ainvnm
201 * ..
202 * .. Local Arrays ..
203  INTEGER isave( 3 )
204 * ..
205 * .. External Functions ..
206  LOGICAL lsame
207  EXTERNAL lsame
208 * ..
209 * .. External Subroutines ..
210  EXTERNAL clacn2, csytrs_3, xerbla
211 * ..
212 * .. Intrinsic Functions ..
213  INTRINSIC max
214 * ..
215 * .. Executable Statements ..
216 *
217 * Test the input parameters.
218 *
219  info = 0
220  upper = lsame( uplo, 'U' )
221  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
222  info = -1
223  ELSE IF( n.LT.0 ) THEN
224  info = -2
225  ELSE IF( lda.LT.max( 1, n ) ) THEN
226  info = -4
227  ELSE IF( anorm.LT.zero ) THEN
228  info = -7
229  END IF
230  IF( info.NE.0 ) THEN
231  CALL xerbla( 'CSYCON_3', -info )
232  RETURN
233  END IF
234 *
235 * Quick return if possible
236 *
237  rcond = zero
238  IF( n.EQ.0 ) THEN
239  rcond = one
240  RETURN
241  ELSE IF( anorm.LE.zero ) THEN
242  RETURN
243  END IF
244 *
245 * Check that the diagonal matrix D is nonsingular.
246 *
247  IF( upper ) THEN
248 *
249 * Upper triangular storage: examine D from bottom to top
250 *
251  DO i = n, 1, -1
252  IF( ipiv( i ).GT.0 .AND. a( i, i ).EQ.czero )
253  $ RETURN
254  END DO
255  ELSE
256 *
257 * Lower triangular storage: examine D from top to bottom.
258 *
259  DO i = 1, n
260  IF( ipiv( i ).GT.0 .AND. a( i, i ).EQ.czero )
261  $ RETURN
262  END DO
263  END IF
264 *
265 * Estimate the 1-norm of the inverse.
266 *
267  kase = 0
268  30 CONTINUE
269  CALL clacn2( n, work( n+1 ), work, ainvnm, kase, isave )
270  IF( kase.NE.0 ) THEN
271 *
272 * Multiply by inv(L*D*L**T) or inv(U*D*U**T).
273 *
274  CALL csytrs_3( uplo, n, 1, a, lda, e, ipiv, work, n, info )
275  GO TO 30
276  END IF
277 *
278 * Compute the estimate of the reciprocal condition number.
279 *
280  IF( ainvnm.NE.zero )
281  $ rcond = ( one / ainvnm ) / anorm
282 *
283  RETURN
284 *
285 * End of CSYCON_3
286 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine csytrs_3(UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB, INFO)
CSYTRS_3
Definition: csytrs_3.f:167
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: