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

SSYCON_3

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

Purpose:
 SSYCON_3 estimates the reciprocal of the condition number (in the
 1-norm) of a real symmetric matrix A using the factorization
 computed by DSYTRF_RK or DSYTRF_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 SSYTRS_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 REAL array, dimension (LDA,N)
          Diagonal of the block diagonal matrix D and factors U or L
          as computed by SSYTRF_RK and SSYTRF_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 REAL 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 SSYTRF_RK or SSYTRF_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 REAL 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 ssycon_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( * ), iwork( * )
186  REAL 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 * ..
195 * .. Local Scalars ..
196  LOGICAL upper
197  INTEGER i, kase
198  REAL ainvnm
199 * ..
200 * .. Local Arrays ..
201  INTEGER isave( 3 )
202 * ..
203 * .. External Functions ..
204  LOGICAL lsame
205  EXTERNAL lsame
206 * ..
207 * .. External Subroutines ..
208  EXTERNAL slacn2, ssytrs_3, xerbla
209 * ..
210 * .. Intrinsic Functions ..
211  INTRINSIC max
212 * ..
213 * .. Executable Statements ..
214 *
215 * Test the input parameters.
216 *
217  info = 0
218  upper = lsame( uplo, 'U' )
219  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
220  info = -1
221  ELSE IF( n.LT.0 ) THEN
222  info = -2
223  ELSE IF( lda.LT.max( 1, n ) ) THEN
224  info = -4
225  ELSE IF( anorm.LT.zero ) THEN
226  info = -7
227  END IF
228  IF( info.NE.0 ) THEN
229  CALL xerbla( 'SSYCON_3', -info )
230  RETURN
231  END IF
232 *
233 * Quick return if possible
234 *
235  rcond = zero
236  IF( n.EQ.0 ) THEN
237  rcond = one
238  RETURN
239  ELSE IF( anorm.LE.zero ) THEN
240  RETURN
241  END IF
242 *
243 * Check that the diagonal matrix D is nonsingular.
244 *
245  IF( upper ) THEN
246 *
247 * Upper triangular storage: examine D from bottom to top
248 *
249  DO i = n, 1, -1
250  IF( ipiv( i ).GT.0 .AND. a( i, i ).EQ.zero )
251  $ RETURN
252  END DO
253  ELSE
254 *
255 * Lower triangular storage: examine D from top to bottom.
256 *
257  DO i = 1, n
258  IF( ipiv( i ).GT.0 .AND. a( i, i ).EQ.zero )
259  $ RETURN
260  END DO
261  END IF
262 *
263 * Estimate the 1-norm of the inverse.
264 *
265  kase = 0
266  30 CONTINUE
267  CALL slacn2( n, work( n+1 ), work, iwork, ainvnm, kase, isave )
268  IF( kase.NE.0 ) THEN
269 *
270 * Multiply by inv(L*D*L**T) or inv(U*D*U**T).
271 *
272  CALL ssytrs_3( uplo, n, 1, a, lda, e, ipiv, work, n, info )
273  GO TO 30
274  END IF
275 *
276 * Compute the estimate of the reciprocal condition number.
277 *
278  IF( ainvnm.NE.zero )
279  $ rcond = ( one / ainvnm ) / anorm
280 *
281  RETURN
282 *
283 * End of DSYCON_3
284 *
subroutine ssytrs_3(UPLO, N, NRHS, A, LDA, E, IPIV, B, LDB, INFO)
SSYTRS_3
Definition: ssytrs_3.f:167
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine slacn2(N, V, X, ISGN, EST, KASE, ISAVE)
SLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: slacn2.f:138

Here is the call graph for this function:

Here is the caller graph for this function: