LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zheequb ( character  UPLO,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  S,
double precision  SCOND,
double precision  AMAX,
complex*16, dimension( * )  WORK,
integer  INFO 
)

ZHEEQUB

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

Purpose:
 ZHEEQUB computes row and column scalings intended to equilibrate a
 Hermitian 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]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The N-by-N Hermitian 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 DOUBLE PRECISION array, dimension (N)
          If INFO = 0, S contains the scale factors for A.
[out]SCOND
          SCOND is DOUBLE PRECISION
          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 DOUBLE PRECISION
          Absolute value of largest matrix element.  If AMAX is very
          close to overflow or very close to underflow, the matrix
          should be scaled.
[out]WORK
          WORK is COMPLEX*16 array, dimension (3*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-th diagonal element is nonpositive.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
April 2012

Definition at line 128 of file zheequb.f.

128 *
129 * -- LAPACK computational routine (version 3.4.1) --
130 * -- LAPACK is a software package provided by Univ. of Tennessee, --
131 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
132 * April 2012
133 *
134 * .. Scalar Arguments ..
135  INTEGER info, lda, n
136  DOUBLE PRECISION amax, scond
137  CHARACTER uplo
138 * ..
139 * .. Array Arguments ..
140  COMPLEX*16 a( lda, * ), work( * )
141  DOUBLE PRECISION s( * )
142 * ..
143 *
144 * =====================================================================
145 *
146 * .. Parameters ..
147  DOUBLE PRECISION one, zero
148  parameter ( one = 1.0d+0, zero = 0.0d+0 )
149  INTEGER max_iter
150  parameter ( max_iter = 100 )
151 * ..
152 * .. Local Scalars ..
153  INTEGER i, j, iter
154  DOUBLE PRECISION avg, std, tol, c0, c1, c2, t, u, si, d,
155  $ base, smin, smax, smlnum, bignum, scale, sumsq
156  LOGICAL up
157  COMPLEX*16 zdum
158 * ..
159 * .. External Functions ..
160  DOUBLE PRECISION dlamch
161  LOGICAL lsame
162  EXTERNAL dlamch, lsame
163 * ..
164 * .. External Subroutines ..
165  EXTERNAL zlassq
166 * ..
167 * .. Intrinsic Functions ..
168  INTRINSIC abs, dble, dimag, int, log, max, min, sqrt
169 * ..
170 * .. Statement Functions ..
171  DOUBLE PRECISION cabs1
172 * ..
173 * .. Statement Function Definitions ..
174  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
175 *
176 * Test input parameters.
177 *
178  info = 0
179  IF (.NOT. ( lsame( uplo, 'U' ) .OR. 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  END IF
186  IF ( info .NE. 0 ) THEN
187  CALL xerbla( 'ZHEEQUB', -info )
188  RETURN
189  END IF
190 
191  up = lsame( uplo, 'U' )
192  amax = zero
193 *
194 * Quick return if possible.
195 *
196  IF ( n .EQ. 0 ) THEN
197  scond = one
198  RETURN
199  END IF
200 
201  DO i = 1, n
202  s( i ) = zero
203  END DO
204 
205  amax = zero
206  IF ( up ) THEN
207  DO j = 1, n
208  DO i = 1, j-1
209  s( i ) = max( s( i ), cabs1( a( i, j ) ) )
210  s( j ) = max( s( j ), cabs1( a( i, j ) ) )
211  amax = max( amax, cabs1( a( i, j ) ) )
212  END DO
213  s( j ) = max( s( j ), cabs1( a( j, j ) ) )
214  amax = max( amax, cabs1( a( j, j ) ) )
215  END DO
216  ELSE
217  DO j = 1, n
218  s( j ) = max( s( j ), cabs1( a( j, j ) ) )
219  amax = max( amax, cabs1( a( j, j ) ) )
220  DO i = j+1, n
221  s( i ) = max( s( i ), cabs1( a( i, j ) ) )
222  s( j ) = max( s( j ), cabs1( a( i, j ) ) )
223  amax = max( amax, cabs1( a(i, j ) ) )
224  END DO
225  END DO
226  END IF
227  DO j = 1, n
228  s( j ) = 1.0d+0 / s( j )
229  END DO
230 
231  tol = one / sqrt( 2.0d0 * n )
232 
233  DO iter = 1, max_iter
234  scale = 0.0d+0
235  sumsq = 0.0d+0
236 * beta = |A|s
237  DO i = 1, n
238  work( i ) = zero
239  END DO
240  IF ( up ) THEN
241  DO j = 1, n
242  DO i = 1, j-1
243  t = cabs1( a( i, j ) )
244  work( i ) = work( i ) + cabs1( a( i, j ) ) * s( j )
245  work( j ) = work( j ) + cabs1( a( i, j ) ) * s( i )
246  END DO
247  work( j ) = work( j ) + cabs1( a( j, j ) ) * s( j )
248  END DO
249  ELSE
250  DO j = 1, n
251  work( j ) = work( j ) + cabs1( a( j, j ) ) * s( j )
252  DO i = j+1, n
253  t = cabs1( a( i, j ) )
254  work( i ) = work( i ) + cabs1( a( i, j ) ) * s( j )
255  work( j ) = work( j ) + cabs1( a( i, j ) ) * s( i )
256  END DO
257  END DO
258  END IF
259 
260 * avg = s^T beta / n
261  avg = 0.0d+0
262  DO i = 1, n
263  avg = avg + s( i )*work( i )
264  END DO
265  avg = avg / n
266 
267  std = 0.0d+0
268  DO i = 2*n+1, 3*n
269  work( i ) = s( i-2*n ) * work( i-2*n ) - avg
270  END DO
271  CALL zlassq( n, work( 2*n+1 ), 1, scale, sumsq )
272  std = scale * sqrt( sumsq / n )
273 
274  IF ( std .LT. tol * avg ) GOTO 999
275 
276  DO i = 1, n
277  t = cabs1( a( i, i ) )
278  si = s( i )
279  c2 = ( n-1 ) * t
280  c1 = ( n-2 ) * ( work( i ) - t*si )
281  c0 = -(t*si)*si + 2*work( i )*si - n*avg
282 
283  d = c1*c1 - 4*c0*c2
284  IF ( d .LE. 0 ) THEN
285  info = -1
286  RETURN
287  END IF
288  si = -2*c0 / ( c1 + sqrt( d ) )
289 
290  d = si - s(i)
291  u = zero
292  IF ( up ) THEN
293  DO j = 1, i
294  t = cabs1( a( j, i ) )
295  u = u + s( j )*t
296  work( j ) = work( j ) + d*t
297  END DO
298  DO j = i+1,n
299  t = cabs1( a( i, j ) )
300  u = u + s( j )*t
301  work( j ) = work( j ) + d*t
302  END DO
303  ELSE
304  DO j = 1, i
305  t = cabs1( a( i, j ) )
306  u = u + s( j )*t
307  work( j ) = work( j ) + d*t
308  END DO
309  DO j = i+1,n
310  t = cabs1( a( j, i ) )
311  u = u + s( j )*t
312  work( j ) = work( j ) + d*t
313  END DO
314  END IF
315  avg = avg + ( u + work( i ) ) * d / n
316  s( i ) = si
317  END DO
318 
319  END DO
320 
321  999 CONTINUE
322 
323  smlnum = dlamch( 'SAFEMIN' )
324  bignum = one / smlnum
325  smin = bignum
326  smax = zero
327  t = one / sqrt( avg )
328  base = dlamch( 'B' )
329  u = one / log( base )
330  DO i = 1, n
331  s( i ) = base ** int( u * log( s( i ) * t ) )
332  smin = min( smin, s( i ) )
333  smax = max( smax, s( i ) )
334  END DO
335  scond = max( smin, smlnum ) / min( smax, bignum )
336 
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlassq(N, X, INCX, SCALE, SUMSQ)
ZLASSQ updates a sum of squares represented in scaled form.
Definition: zlassq.f:108
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: