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