LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dsycon_3()

subroutine dsycon_3 ( character  uplo,
integer  n,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( * )  e,
integer, dimension( * )  ipiv,
double precision  anorm,
double precision  rcond,
double precision, dimension( * )  work,
integer, dimension( * )  iwork,
integer  info 
)

DSYCON_3

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

Purpose:
 DSYCON_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 DSYTRS_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 DOUBLE PRECISION array, dimension (LDA,N)
          Diagonal of the block diagonal matrix D and factors U or L
          as computed by DSYTRF_RK and DSYTRF_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 DOUBLE PRECISION 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 referenced;
          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 DSYTRF_RK or DSYTRF_BK.
[in]ANORM
          ANORM is DOUBLE PRECISION
          The 1-norm of the original matrix A.
[out]RCOND
          RCOND is DOUBLE PRECISION
          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 DOUBLE PRECISION 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.
Contributors:
  June 2017,  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 169 of file dsycon_3.f.

171*
172* -- LAPACK computational routine --
173* -- LAPACK is a software package provided by Univ. of Tennessee, --
174* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176* .. Scalar Arguments ..
177 CHARACTER UPLO
178 INTEGER INFO, LDA, N
179 DOUBLE PRECISION ANORM, RCOND
180* ..
181* .. Array Arguments ..
182 INTEGER IPIV( * ), IWORK( * )
183 DOUBLE PRECISION A( LDA, * ), E( * ), WORK( * )
184* ..
185*
186* =====================================================================
187*
188* .. Parameters ..
189 DOUBLE PRECISION ONE, ZERO
190 parameter( one = 1.0d+0, zero = 0.0d+0 )
191* ..
192* .. Local Scalars ..
193 LOGICAL UPPER
194 INTEGER I, KASE
195 DOUBLE PRECISION AINVNM
196* ..
197* .. Local Arrays ..
198 INTEGER ISAVE( 3 )
199* ..
200* .. External Functions ..
201 LOGICAL LSAME
202 EXTERNAL lsame
203* ..
204* .. External Subroutines ..
205 EXTERNAL dlacn2, dsytrs_3, xerbla
206* ..
207* .. Intrinsic Functions ..
208 INTRINSIC max
209* ..
210* .. Executable Statements ..
211*
212* Test the input parameters.
213*
214 info = 0
215 upper = lsame( uplo, 'U' )
216 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
217 info = -1
218 ELSE IF( n.LT.0 ) THEN
219 info = -2
220 ELSE IF( lda.LT.max( 1, n ) ) THEN
221 info = -4
222 ELSE IF( anorm.LT.zero ) THEN
223 info = -7
224 END IF
225 IF( info.NE.0 ) THEN
226 CALL xerbla( 'DSYCON_3', -info )
227 RETURN
228 END IF
229*
230* Quick return if possible
231*
232 rcond = zero
233 IF( n.EQ.0 ) THEN
234 rcond = one
235 RETURN
236 ELSE IF( anorm.LE.zero ) THEN
237 RETURN
238 END IF
239*
240* Check that the diagonal matrix D is nonsingular.
241*
242 IF( upper ) THEN
243*
244* Upper triangular storage: examine D from bottom to top
245*
246 DO i = n, 1, -1
247 IF( ipiv( i ).GT.0 .AND. a( i, i ).EQ.zero )
248 $ RETURN
249 END DO
250 ELSE
251*
252* Lower triangular storage: examine D from top to bottom.
253*
254 DO i = 1, n
255 IF( ipiv( i ).GT.0 .AND. a( i, i ).EQ.zero )
256 $ RETURN
257 END DO
258 END IF
259*
260* Estimate the 1-norm of the inverse.
261*
262 kase = 0
263 30 CONTINUE
264 CALL dlacn2( n, work( n+1 ), work, iwork, ainvnm, kase, isave )
265 IF( kase.NE.0 ) THEN
266*
267* Multiply by inv(L*D*L**T) or inv(U*D*U**T).
268*
269 CALL dsytrs_3( uplo, n, 1, a, lda, e, ipiv, work, n, info )
270 GO TO 30
271 END IF
272*
273* Compute the estimate of the reciprocal condition number.
274*
275 IF( ainvnm.NE.zero )
276 $ rcond = ( one / ainvnm ) / anorm
277*
278 RETURN
279*
280* End of DSYCON_3
281*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dsytrs_3(uplo, n, nrhs, a, lda, e, ipiv, b, ldb, info)
DSYTRS_3
Definition dsytrs_3.f:165
subroutine dlacn2(n, v, x, isgn, est, kase, isave)
DLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition dlacn2.f:136
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: