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

◆ zgebal()

subroutine zgebal ( character  JOB,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer  ILO,
integer  IHI,
double precision, dimension( * )  SCALE,
integer  INFO 
)

ZGEBAL

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

Purpose:
 ZGEBAL balances a general complex matrix A.  This involves, first,
 permuting A by a similarity transformation to isolate eigenvalues
 in the first 1 to ILO-1 and last IHI+1 to N elements on the
 diagonal; and second, applying a diagonal similarity transformation
 to rows and columns ILO to IHI to make the rows and columns as
 close in norm as possible.  Both steps are optional.

 Balancing may reduce the 1-norm of the matrix, and improve the
 accuracy of the computed eigenvalues and/or eigenvectors.
Parameters
[in]JOB
          JOB is CHARACTER*1
          Specifies the operations to be performed on A:
          = 'N':  none:  simply set ILO = 1, IHI = N, SCALE(I) = 1.0
                  for i = 1,...,N;
          = 'P':  permute only;
          = 'S':  scale only;
          = 'B':  both permute and scale.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the input matrix A.
          On exit,  A is overwritten by the balanced matrix.
          If JOB = 'N', A is not referenced.
          See Further Details.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]ILO
          ILO is INTEGER
[out]IHI
          IHI is INTEGER
          ILO and IHI are set to INTEGER such that on exit
          A(i,j) = 0 if i > j and j = 1,...,ILO-1 or I = IHI+1,...,N.
          If JOB = 'N' or 'S', ILO = 1 and IHI = N.
[out]SCALE
          SCALE is DOUBLE PRECISION array, dimension (N)
          Details of the permutations and scaling factors applied to
          A.  If P(j) is the index of the row and column interchanged
          with row and column j and D(j) is the scaling factor
          applied to row and column j, then
          SCALE(j) = P(j)    for j = 1,...,ILO-1
                   = D(j)    for j = ILO,...,IHI
                   = P(j)    for j = IHI+1,...,N.
          The order in which the interchanges are made is N to IHI+1,
          then 1 to ILO-1.
[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.
Further Details:
  The permutations consist of row and column interchanges which put
  the matrix in the form

             ( T1   X   Y  )
     P A P = (  0   B   Z  )
             (  0   0   T2 )

  where T1 and T2 are upper triangular matrices whose eigenvalues lie
  along the diagonal.  The column indices ILO and IHI mark the starting
  and ending columns of the submatrix B. Balancing consists of applying
  a diagonal similarity transformation inv(D) * B * D to make the
  1-norms of each row of B and its corresponding column nearly equal.
  The output matrix is

     ( T1     X*D          Y    )
     (  0  inv(D)*B*D  inv(D)*Z ).
     (  0      0           T2   )

  Information about the permutations P and the diagonal matrix D is
  returned in the vector SCALE.

  This subroutine is based on the EISPACK routine CBAL.

  Modified by Tzu-Yi Chen, Computer Science Division, University of
    California at Berkeley, USA

Definition at line 161 of file zgebal.f.

162*
163* -- LAPACK computational routine --
164* -- LAPACK is a software package provided by Univ. of Tennessee, --
165* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
166*
167* .. Scalar Arguments ..
168 CHARACTER JOB
169 INTEGER IHI, ILO, INFO, LDA, N
170* ..
171* .. Array Arguments ..
172 DOUBLE PRECISION SCALE( * )
173 COMPLEX*16 A( LDA, * )
174* ..
175*
176* =====================================================================
177*
178* .. Parameters ..
179 DOUBLE PRECISION ZERO, ONE
180 parameter( zero = 0.0d+0, one = 1.0d+0 )
181 DOUBLE PRECISION SCLFAC
182 parameter( sclfac = 2.0d+0 )
183 DOUBLE PRECISION FACTOR
184 parameter( factor = 0.95d+0 )
185* ..
186* .. Local Scalars ..
187 LOGICAL NOCONV
188 INTEGER I, ICA, IEXC, IRA, J, K, L, M
189 DOUBLE PRECISION C, CA, F, G, R, RA, S, SFMAX1, SFMAX2, SFMIN1,
190 $ SFMIN2
191* ..
192* .. External Functions ..
193 LOGICAL DISNAN, LSAME
194 INTEGER IZAMAX
195 DOUBLE PRECISION DLAMCH, DZNRM2
196 EXTERNAL disnan, lsame, izamax, dlamch, dznrm2
197* ..
198* .. External Subroutines ..
199 EXTERNAL xerbla, zdscal, zswap
200* ..
201* .. Intrinsic Functions ..
202 INTRINSIC abs, dble, dimag, max, min
203*
204* Test the input parameters
205*
206 info = 0
207 IF( .NOT.lsame( job, 'N' ) .AND. .NOT.lsame( job, 'P' ) .AND.
208 $ .NOT.lsame( job, 'S' ) .AND. .NOT.lsame( job, 'B' ) ) THEN
209 info = -1
210 ELSE IF( n.LT.0 ) THEN
211 info = -2
212 ELSE IF( lda.LT.max( 1, n ) ) THEN
213 info = -4
214 END IF
215 IF( info.NE.0 ) THEN
216 CALL xerbla( 'ZGEBAL', -info )
217 RETURN
218 END IF
219*
220 k = 1
221 l = n
222*
223 IF( n.EQ.0 )
224 $ GO TO 210
225*
226 IF( lsame( job, 'N' ) ) THEN
227 DO 10 i = 1, n
228 scale( i ) = one
229 10 CONTINUE
230 GO TO 210
231 END IF
232*
233 IF( lsame( job, 'S' ) )
234 $ GO TO 120
235*
236* Permutation to isolate eigenvalues if possible
237*
238 GO TO 50
239*
240* Row and column exchange.
241*
242 20 CONTINUE
243 scale( m ) = j
244 IF( j.EQ.m )
245 $ GO TO 30
246*
247 CALL zswap( l, a( 1, j ), 1, a( 1, m ), 1 )
248 CALL zswap( n-k+1, a( j, k ), lda, a( m, k ), lda )
249*
250 30 CONTINUE
251 GO TO ( 40, 80 )iexc
252*
253* Search for rows isolating an eigenvalue and push them down.
254*
255 40 CONTINUE
256 IF( l.EQ.1 )
257 $ GO TO 210
258 l = l - 1
259*
260 50 CONTINUE
261 DO 70 j = l, 1, -1
262*
263 DO 60 i = 1, l
264 IF( i.EQ.j )
265 $ GO TO 60
266 IF( dble( a( j, i ) ).NE.zero .OR. dimag( a( j, i ) ).NE.
267 $ zero )GO TO 70
268 60 CONTINUE
269*
270 m = l
271 iexc = 1
272 GO TO 20
273 70 CONTINUE
274*
275 GO TO 90
276*
277* Search for columns isolating an eigenvalue and push them left.
278*
279 80 CONTINUE
280 k = k + 1
281*
282 90 CONTINUE
283 DO 110 j = k, l
284*
285 DO 100 i = k, l
286 IF( i.EQ.j )
287 $ GO TO 100
288 IF( dble( a( i, j ) ).NE.zero .OR. dimag( a( i, j ) ).NE.
289 $ zero )GO TO 110
290 100 CONTINUE
291*
292 m = k
293 iexc = 2
294 GO TO 20
295 110 CONTINUE
296*
297 120 CONTINUE
298 DO 130 i = k, l
299 scale( i ) = one
300 130 CONTINUE
301*
302 IF( lsame( job, 'P' ) )
303 $ GO TO 210
304*
305* Balance the submatrix in rows K to L.
306*
307* Iterative loop for norm reduction
308*
309 sfmin1 = dlamch( 'S' ) / dlamch( 'P' )
310 sfmax1 = one / sfmin1
311 sfmin2 = sfmin1*sclfac
312 sfmax2 = one / sfmin2
313 140 CONTINUE
314 noconv = .false.
315*
316 DO 200 i = k, l
317*
318 c = dznrm2( l-k+1, a( k, i ), 1 )
319 r = dznrm2( l-k+1, a( i, k ), lda )
320 ica = izamax( l, a( 1, i ), 1 )
321 ca = abs( a( ica, i ) )
322 ira = izamax( n-k+1, a( i, k ), lda )
323 ra = abs( a( i, ira+k-1 ) )
324*
325* Guard against zero C or R due to underflow.
326*
327 IF( c.EQ.zero .OR. r.EQ.zero )
328 $ GO TO 200
329 g = r / sclfac
330 f = one
331 s = c + r
332 160 CONTINUE
333 IF( c.GE.g .OR. max( f, c, ca ).GE.sfmax2 .OR.
334 $ min( r, g, ra ).LE.sfmin2 )GO TO 170
335 IF( disnan( c+f+ca+r+g+ra ) ) THEN
336*
337* Exit if NaN to avoid infinite loop
338*
339 info = -3
340 CALL xerbla( 'ZGEBAL', -info )
341 RETURN
342 END IF
343 f = f*sclfac
344 c = c*sclfac
345 ca = ca*sclfac
346 r = r / sclfac
347 g = g / sclfac
348 ra = ra / sclfac
349 GO TO 160
350*
351 170 CONTINUE
352 g = c / sclfac
353 180 CONTINUE
354 IF( g.LT.r .OR. max( r, ra ).GE.sfmax2 .OR.
355 $ min( f, c, g, ca ).LE.sfmin2 )GO TO 190
356 f = f / sclfac
357 c = c / sclfac
358 g = g / sclfac
359 ca = ca / sclfac
360 r = r*sclfac
361 ra = ra*sclfac
362 GO TO 180
363*
364* Now balance.
365*
366 190 CONTINUE
367 IF( ( c+r ).GE.factor*s )
368 $ GO TO 200
369 IF( f.LT.one .AND. scale( i ).LT.one ) THEN
370 IF( f*scale( i ).LE.sfmin1 )
371 $ GO TO 200
372 END IF
373 IF( f.GT.one .AND. scale( i ).GT.one ) THEN
374 IF( scale( i ).GE.sfmax1 / f )
375 $ GO TO 200
376 END IF
377 g = one / f
378 scale( i ) = scale( i )*f
379 noconv = .true.
380*
381 CALL zdscal( n-k+1, g, a( i, k ), lda )
382 CALL zdscal( l, f, a( 1, i ), 1 )
383*
384 200 CONTINUE
385*
386 IF( noconv )
387 $ GO TO 140
388*
389 210 CONTINUE
390 ilo = k
391 ihi = l
392*
393 RETURN
394*
395* End of ZGEBAL
396*
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
logical function disnan(DIN)
DISNAN tests input for NaN.
Definition: disnan.f:59
integer function izamax(N, ZX, INCX)
IZAMAX
Definition: izamax.f:71
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:81
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:78
real(wp) function dznrm2(n, x, incx)
DZNRM2
Definition: dznrm2.f90:90
Here is the call graph for this function:
Here is the caller graph for this function: