LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zgebal.f
Go to the documentation of this file.
1*> \brief \b ZGEBAL
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZGEBAL + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgebal.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgebal.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgebal.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZGEBAL( JOB, N, A, LDA, ILO, IHI, SCALE, INFO )
20*
21* .. Scalar Arguments ..
22* CHARACTER JOB
23* INTEGER IHI, ILO, INFO, LDA, N
24* ..
25* .. Array Arguments ..
26* DOUBLE PRECISION SCALE( * )
27* COMPLEX*16 A( LDA, * )
28* ..
29*
30*
31*> \par Purpose:
32* =============
33*>
34*> \verbatim
35*>
36*> ZGEBAL balances a general complex matrix A. This involves, first,
37*> permuting A by a similarity transformation to isolate eigenvalues
38*> in the first 1 to ILO-1 and last IHI+1 to N elements on the
39*> diagonal; and second, applying a diagonal similarity transformation
40*> to rows and columns ILO to IHI to make the rows and columns as
41*> close in norm as possible. Both steps are optional.
42*>
43*> Balancing may reduce the 1-norm of the matrix, and improve the
44*> accuracy of the computed eigenvalues and/or eigenvectors.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] JOB
51*> \verbatim
52*> JOB is CHARACTER*1
53*> Specifies the operations to be performed on A:
54*> = 'N': none: simply set ILO = 1, IHI = N, SCALE(I) = 1.0
55*> for i = 1,...,N;
56*> = 'P': permute only;
57*> = 'S': scale only;
58*> = 'B': both permute and scale.
59*> \endverbatim
60*>
61*> \param[in] N
62*> \verbatim
63*> N is INTEGER
64*> The order of the matrix A. N >= 0.
65*> \endverbatim
66*>
67*> \param[in,out] A
68*> \verbatim
69*> A is COMPLEX*16 array, dimension (LDA,N)
70*> On entry, the input matrix A.
71*> On exit, A is overwritten by the balanced matrix.
72*> If JOB = 'N', A is not referenced.
73*> See Further Details.
74*> \endverbatim
75*>
76*> \param[in] LDA
77*> \verbatim
78*> LDA is INTEGER
79*> The leading dimension of the array A. LDA >= max(1,N).
80*> \endverbatim
81*>
82*> \param[out] ILO
83*> \verbatim
84*> ILO is INTEGER
85*> \endverbatim
86*>
87*> \param[out] IHI
88*> \verbatim
89*> IHI is INTEGER
90*> ILO and IHI are set to integers such that on exit
91*> A(i,j) = 0 if i > j and j = 1,...,ILO-1 or I = IHI+1,...,N.
92*> If JOB = 'N' or 'S', ILO = 1 and IHI = N.
93*> \endverbatim
94*>
95*> \param[out] SCALE
96*> \verbatim
97*> SCALE is DOUBLE PRECISION array, dimension (N)
98*> Details of the permutations and scaling factors applied to
99*> A. If P(j) is the index of the row and column interchanged
100*> with row and column j and D(j) is the scaling factor
101*> applied to row and column j, then
102*> SCALE(j) = P(j) for j = 1,...,ILO-1
103*> = D(j) for j = ILO,...,IHI
104*> = P(j) for j = IHI+1,...,N.
105*> The order in which the interchanges are made is N to IHI+1,
106*> then 1 to ILO-1.
107*> \endverbatim
108*>
109*> \param[out] INFO
110*> \verbatim
111*> INFO is INTEGER
112*> = 0: successful exit.
113*> < 0: if INFO = -i, the i-th argument had an illegal value.
114*> \endverbatim
115*
116* Authors:
117* ========
118*
119*> \author Univ. of Tennessee
120*> \author Univ. of California Berkeley
121*> \author Univ. of Colorado Denver
122*> \author NAG Ltd.
123*
124*> \ingroup gebal
125*
126*> \par Further Details:
127* =====================
128*>
129*> \verbatim
130*>
131*> The permutations consist of row and column interchanges which put
132*> the matrix in the form
133*>
134*> ( T1 X Y )
135*> P A P = ( 0 B Z )
136*> ( 0 0 T2 )
137*>
138*> where T1 and T2 are upper triangular matrices whose eigenvalues lie
139*> along the diagonal. The column indices ILO and IHI mark the starting
140*> and ending columns of the submatrix B. Balancing consists of applying
141*> a diagonal similarity transformation inv(D) * B * D to make the
142*> 1-norms of each row of B and its corresponding column nearly equal.
143*> The output matrix is
144*>
145*> ( T1 X*D Y )
146*> ( 0 inv(D)*B*D inv(D)*Z ).
147*> ( 0 0 T2 )
148*>
149*> Information about the permutations P and the diagonal matrix D is
150*> returned in the vector SCALE.
151*>
152*> This subroutine is based on the EISPACK routine CBAL.
153*>
154*> Modified by Tzu-Yi Chen, Computer Science Division, University of
155*> California at Berkeley, USA
156*>
157*> Refactored by Evert Provoost, Department of Computer Science,
158*> KU Leuven, Belgium
159*> \endverbatim
160*>
161* =====================================================================
162 SUBROUTINE zgebal( JOB, N, A, LDA, ILO, IHI, SCALE, INFO )
163*
164* -- LAPACK computational routine --
165* -- LAPACK is a software package provided by Univ. of Tennessee, --
166* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
167*
168* .. Scalar Arguments ..
169 CHARACTER JOB
170 INTEGER IHI, ILO, INFO, LDA, N
171* ..
172* .. Array Arguments ..
173 DOUBLE PRECISION SCALE( * )
174 COMPLEX*16 A( LDA, * )
175* ..
176*
177* =====================================================================
178*
179* .. Parameters ..
180 DOUBLE PRECISION ZERO, ONE
181 parameter( zero = 0.0d+0, one = 1.0d+0 )
182 DOUBLE PRECISION SCLFAC
183 parameter( sclfac = 2.0d+0 )
184 DOUBLE PRECISION FACTOR
185 parameter( factor = 0.95d+0 )
186* ..
187* .. Local Scalars ..
188 LOGICAL NOCONV, CANSWAP
189 INTEGER I, ICA, IRA, J, K, L
190 DOUBLE PRECISION C, CA, F, G, R, RA, S, SFMAX1, SFMAX2, SFMIN1,
191 $ SFMIN2
192* ..
193* .. External Functions ..
194 LOGICAL DISNAN, LSAME
195 INTEGER IZAMAX
196 DOUBLE PRECISION DLAMCH, DZNRM2
197 EXTERNAL disnan, lsame, izamax, dlamch,
198 $ dznrm2
199* ..
200* .. External Subroutines ..
201 EXTERNAL xerbla, zdscal, zswap
202* ..
203* .. Intrinsic Functions ..
204 INTRINSIC abs, dble, dimag, max, min
205*
206* Test the input parameters
207*
208 info = 0
209 IF( .NOT.lsame( job, 'N' ) .AND.
210 $ .NOT.lsame( job, 'P' ) .AND.
211 $ .NOT.lsame( job, 'S' ) .AND.
212 $ .NOT.lsame( job, 'B' ) ) THEN
213 info = -1
214 ELSE IF( n.LT.0 ) THEN
215 info = -2
216 ELSE IF( lda.LT.max( 1, n ) ) THEN
217 info = -4
218 END IF
219 IF( info.NE.0 ) THEN
220 CALL xerbla( 'ZGEBAL', -info )
221 RETURN
222 END IF
223*
224* Quick returns.
225*
226 IF( n.EQ.0 ) THEN
227 ilo = 1
228 ihi = 0
229 RETURN
230 END IF
231*
232 IF( lsame( job, 'N' ) ) THEN
233 DO i = 1, n
234 scale( i ) = one
235 END DO
236 ilo = 1
237 ihi = n
238 RETURN
239 END IF
240*
241* Permutation to isolate eigenvalues if possible.
242*
243 k = 1
244 l = n
245*
246 IF( .NOT.lsame( job, 'S' ) ) THEN
247*
248* Row and column exchange.
249*
250 noconv = .true.
251 DO WHILE( noconv )
252*
253* Search for rows isolating an eigenvalue and push them down.
254*
255 noconv = .false.
256 DO i = l, 1, -1
257 canswap = .true.
258 DO j = 1, l
259 IF( i.NE.j .AND. ( dble( a( i, j ) ).NE.zero .OR.
260 $ dimag( a( i, j ) ).NE.zero ) ) THEN
261 canswap = .false.
262 EXIT
263 END IF
264 END DO
265*
266 IF( canswap ) THEN
267 scale( l ) = i
268 IF( i.NE.l ) THEN
269 CALL zswap( l, a( 1, i ), 1, a( 1, l ), 1 )
270 CALL zswap( n-k+1, a( i, k ), lda, a( l, k ),
271 $ lda )
272 END IF
273 noconv = .true.
274*
275 IF( l.EQ.1 ) THEN
276 ilo = 1
277 ihi = 1
278 RETURN
279 END IF
280*
281 l = l - 1
282 END IF
283 END DO
284*
285 END DO
286
287 noconv = .true.
288 DO WHILE( noconv )
289*
290* Search for columns isolating an eigenvalue and push them left.
291*
292 noconv = .false.
293 DO j = k, l
294 canswap = .true.
295 DO i = k, l
296 IF( i.NE.j .AND. ( dble( a( i, j ) ).NE.zero .OR.
297 $ dimag( a( i, j ) ).NE.zero ) ) THEN
298 canswap = .false.
299 EXIT
300 END IF
301 END DO
302*
303 IF( canswap ) THEN
304 scale( k ) = j
305 IF( j.NE.k ) THEN
306 CALL zswap( l, a( 1, j ), 1, a( 1, k ), 1 )
307 CALL zswap( n-k+1, a( j, k ), lda, a( k, k ),
308 $ lda )
309 END IF
310 noconv = .true.
311*
312 k = k + 1
313 END IF
314 END DO
315*
316 END DO
317*
318 END IF
319*
320* Initialize SCALE for non-permuted submatrix.
321*
322 DO i = k, l
323 scale( i ) = one
324 END DO
325*
326* If we only had to permute, we are done.
327*
328 IF( lsame( job, 'P' ) ) THEN
329 ilo = k
330 ihi = l
331 RETURN
332 END IF
333*
334* Balance the submatrix in rows K to L.
335*
336* Iterative loop for norm reduction.
337*
338 sfmin1 = dlamch( 'S' ) / dlamch( 'P' )
339 sfmax1 = one / sfmin1
340 sfmin2 = sfmin1*sclfac
341 sfmax2 = one / sfmin2
342*
343 noconv = .true.
344 DO WHILE( noconv )
345 noconv = .false.
346*
347 DO i = k, l
348*
349 c = dznrm2( l-k+1, a( k, i ), 1 )
350 r = dznrm2( l-k+1, a( i, k ), lda )
351 ica = izamax( l, a( 1, i ), 1 )
352 ca = abs( a( ica, i ) )
353 ira = izamax( n-k+1, a( i, k ), lda )
354 ra = abs( a( i, ira+k-1 ) )
355*
356* Guard against zero C or R due to underflow.
357*
358 IF( c.EQ.zero .OR. r.EQ.zero ) cycle
359*
360* Exit if NaN to avoid infinite loop
361*
362 IF( disnan( c+ca+r+ra ) ) THEN
363 info = -3
364 CALL xerbla( 'ZGEBAL', -info )
365 RETURN
366 END IF
367*
368 g = r / sclfac
369 f = one
370 s = c + r
371*
372 DO WHILE( c.LT.g .AND. max( f, c, ca ).LT.sfmax2 .AND.
373 $ min( r, g, ra ).GT.sfmin2 )
374 f = f*sclfac
375 c = c*sclfac
376 ca = ca*sclfac
377 r = r / sclfac
378 g = g / sclfac
379 ra = ra / sclfac
380 END DO
381*
382 g = c / sclfac
383*
384 DO WHILE( g.GE.r .AND. max( r, ra ).LT.sfmax2 .AND.
385 $ min( f, c, g, ca ).GT.sfmin2 )
386 f = f / sclfac
387 c = c / sclfac
388 g = g / sclfac
389 ca = ca / sclfac
390 r = r*sclfac
391 ra = ra*sclfac
392 END DO
393*
394* Now balance.
395*
396 IF( ( c+r ).GE.factor*s ) cycle
397 IF( f.LT.one .AND. scale( i ).LT.one ) THEN
398 IF( f*scale( i ).LE.sfmin1 ) cycle
399 END IF
400 IF( f.GT.one .AND. scale( i ).GT.one ) THEN
401 IF( scale( i ).GE.sfmax1 / f ) cycle
402 END IF
403 g = one / f
404 scale( i ) = scale( i )*f
405 noconv = .true.
406*
407 CALL zdscal( n-k+1, g, a( i, k ), lda )
408 CALL zdscal( l, f, a( 1, i ), 1 )
409*
410 END DO
411*
412 END DO
413*
414 ilo = k
415 ihi = l
416*
417 RETURN
418*
419* End of ZGEBAL
420*
421 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgebal(job, n, a, lda, ilo, ihi, scale, info)
ZGEBAL
Definition zgebal.f:163
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81