LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
cggbal.f
Go to the documentation of this file.
1 *> \brief \b CGGBAL
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download CGGBAL + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cggbal.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cggbal.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cggbal.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE CGGBAL( JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE,
22 * RSCALE, WORK, INFO )
23 *
24 * .. Scalar Arguments ..
25 * CHARACTER JOB
26 * INTEGER IHI, ILO, INFO, LDA, LDB, N
27 * ..
28 * .. Array Arguments ..
29 * REAL LSCALE( * ), RSCALE( * ), WORK( * )
30 * COMPLEX A( LDA, * ), B( LDB, * )
31 * ..
32 *
33 *
34 *> \par Purpose:
35 * =============
36 *>
37 *> \verbatim
38 *>
39 *> CGGBAL balances a pair of general complex matrices (A,B). This
40 *> involves, first, permuting A and B by similarity transformations to
41 *> isolate eigenvalues in the first 1 to ILO$-$1 and last IHI+1 to N
42 *> elements on the diagonal; and second, applying a diagonal similarity
43 *> transformation to rows and columns ILO to IHI to make the rows
44 *> and columns as close in norm as possible. Both steps are optional.
45 *>
46 *> Balancing may reduce the 1-norm of the matrices, and improve the
47 *> accuracy of the computed eigenvalues and/or eigenvectors in the
48 *> generalized eigenvalue problem A*x = lambda*B*x.
49 *> \endverbatim
50 *
51 * Arguments:
52 * ==========
53 *
54 *> \param[in] JOB
55 *> \verbatim
56 *> JOB is CHARACTER*1
57 *> Specifies the operations to be performed on A and B:
58 *> = 'N': none: simply set ILO = 1, IHI = N, LSCALE(I) = 1.0
59 *> and RSCALE(I) = 1.0 for i=1,...,N;
60 *> = 'P': permute only;
61 *> = 'S': scale only;
62 *> = 'B': both permute and scale.
63 *> \endverbatim
64 *>
65 *> \param[in] N
66 *> \verbatim
67 *> N is INTEGER
68 *> The order of the matrices A and B. N >= 0.
69 *> \endverbatim
70 *>
71 *> \param[in,out] A
72 *> \verbatim
73 *> A is COMPLEX array, dimension (LDA,N)
74 *> On entry, the input matrix A.
75 *> On exit, A is overwritten by the balanced matrix.
76 *> If JOB = 'N', A is not referenced.
77 *> \endverbatim
78 *>
79 *> \param[in] LDA
80 *> \verbatim
81 *> LDA is INTEGER
82 *> The leading dimension of the array A. LDA >= max(1,N).
83 *> \endverbatim
84 *>
85 *> \param[in,out] B
86 *> \verbatim
87 *> B is COMPLEX array, dimension (LDB,N)
88 *> On entry, the input matrix B.
89 *> On exit, B is overwritten by the balanced matrix.
90 *> If JOB = 'N', B is not referenced.
91 *> \endverbatim
92 *>
93 *> \param[in] LDB
94 *> \verbatim
95 *> LDB is INTEGER
96 *> The leading dimension of the array B. LDB >= max(1,N).
97 *> \endverbatim
98 *>
99 *> \param[out] ILO
100 *> \verbatim
101 *> ILO is INTEGER
102 *> \endverbatim
103 *>
104 *> \param[out] IHI
105 *> \verbatim
106 *> IHI is INTEGER
107 *> ILO and IHI are set to integers such that on exit
108 *> A(i,j) = 0 and B(i,j) = 0 if i > j and
109 *> j = 1,...,ILO-1 or i = IHI+1,...,N.
110 *> If JOB = 'N' or 'S', ILO = 1 and IHI = N.
111 *> \endverbatim
112 *>
113 *> \param[out] LSCALE
114 *> \verbatim
115 *> LSCALE is REAL array, dimension (N)
116 *> Details of the permutations and scaling factors applied
117 *> to the left side of A and B. If P(j) is the index of the
118 *> row interchanged with row j, and D(j) is the scaling factor
119 *> applied to row j, then
120 *> LSCALE(j) = P(j) for J = 1,...,ILO-1
121 *> = D(j) for J = ILO,...,IHI
122 *> = P(j) for J = IHI+1,...,N.
123 *> The order in which the interchanges are made is N to IHI+1,
124 *> then 1 to ILO-1.
125 *> \endverbatim
126 *>
127 *> \param[out] RSCALE
128 *> \verbatim
129 *> RSCALE is REAL array, dimension (N)
130 *> Details of the permutations and scaling factors applied
131 *> to the right side of A and B. If P(j) is the index of the
132 *> column interchanged with column j, and D(j) is the scaling
133 *> factor applied to column j, then
134 *> RSCALE(j) = P(j) for J = 1,...,ILO-1
135 *> = D(j) for J = ILO,...,IHI
136 *> = P(j) for J = IHI+1,...,N.
137 *> The order in which the interchanges are made is N to IHI+1,
138 *> then 1 to ILO-1.
139 *> \endverbatim
140 *>
141 *> \param[out] WORK
142 *> \verbatim
143 *> WORK is REAL array, dimension (lwork)
144 *> lwork must be at least max(1,6*N) when JOB = 'S' or 'B', and
145 *> at least 1 when JOB = 'N' or 'P'.
146 *> \endverbatim
147 *>
148 *> \param[out] INFO
149 *> \verbatim
150 *> INFO is INTEGER
151 *> = 0: successful exit
152 *> < 0: if INFO = -i, the i-th argument had an illegal value.
153 *> \endverbatim
154 *
155 * Authors:
156 * ========
157 *
158 *> \author Univ. of Tennessee
159 *> \author Univ. of California Berkeley
160 *> \author Univ. of Colorado Denver
161 *> \author NAG Ltd.
162 *
163 *> \date November 2011
164 *
165 *> \ingroup complexGBcomputational
166 *
167 *> \par Further Details:
168 * =====================
169 *>
170 *> \verbatim
171 *>
172 *> See R.C. WARD, Balancing the generalized eigenvalue problem,
173 *> SIAM J. Sci. Stat. Comp. 2 (1981), 141-152.
174 *> \endverbatim
175 *>
176 * =====================================================================
177  SUBROUTINE cggbal( JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE,
178  $ rscale, work, info )
179 *
180 * -- LAPACK computational routine (version 3.4.0) --
181 * -- LAPACK is a software package provided by Univ. of Tennessee, --
182 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 * November 2011
184 *
185 * .. Scalar Arguments ..
186  CHARACTER job
187  INTEGER ihi, ilo, info, lda, ldb, n
188 * ..
189 * .. Array Arguments ..
190  REAL lscale( * ), rscale( * ), work( * )
191  COMPLEX a( lda, * ), b( ldb, * )
192 * ..
193 *
194 * =====================================================================
195 *
196 * .. Parameters ..
197  REAL zero, half, one
198  parameter( zero = 0.0e+0, half = 0.5e+0, one = 1.0e+0 )
199  REAL three, sclfac
200  parameter( three = 3.0e+0, sclfac = 1.0e+1 )
201  COMPLEX czero
202  parameter( czero = ( 0.0e+0, 0.0e+0 ) )
203 * ..
204 * .. Local Scalars ..
205  INTEGER i, icab, iflow, ip1, ir, irab, it, j, jc, jp1,
206  $ k, kount, l, lcab, lm1, lrab, lsfmax, lsfmin,
207  $ m, nr, nrp2
208  REAL alpha, basl, beta, cab, cmax, coef, coef2,
209  $ coef5, cor, ew, ewc, gamma, pgamma, rab, sfmax,
210  $ sfmin, sum, t, ta, tb, tc
211  COMPLEX cdum
212 * ..
213 * .. External Functions ..
214  LOGICAL lsame
215  INTEGER icamax
216  REAL sdot, slamch
217  EXTERNAL lsame, icamax, sdot, slamch
218 * ..
219 * .. External Subroutines ..
220  EXTERNAL csscal, cswap, saxpy, sscal, xerbla
221 * ..
222 * .. Intrinsic Functions ..
223  INTRINSIC abs, aimag, int, log10, max, min, REAL, sign
224 * ..
225 * .. Statement Functions ..
226  REAL cabs1
227 * ..
228 * .. Statement Function definitions ..
229  cabs1( cdum ) = abs( REAL( CDUM ) ) + abs( aimag( cdum ) )
230 * ..
231 * .. Executable Statements ..
232 *
233 * Test the input parameters
234 *
235  info = 0
236  IF( .NOT.lsame( job, 'N' ) .AND. .NOT.lsame( job, 'P' ) .AND.
237  $ .NOT.lsame( job, 'S' ) .AND. .NOT.lsame( job, 'B' ) ) THEN
238  info = -1
239  ELSE IF( n.LT.0 ) THEN
240  info = -2
241  ELSE IF( lda.LT.max( 1, n ) ) THEN
242  info = -4
243  ELSE IF( ldb.LT.max( 1, n ) ) THEN
244  info = -6
245  END IF
246  IF( info.NE.0 ) THEN
247  CALL xerbla( 'CGGBAL', -info )
248  return
249  END IF
250 *
251 * Quick return if possible
252 *
253  IF( n.EQ.0 ) THEN
254  ilo = 1
255  ihi = n
256  return
257  END IF
258 *
259  IF( n.EQ.1 ) THEN
260  ilo = 1
261  ihi = n
262  lscale( 1 ) = one
263  rscale( 1 ) = one
264  return
265  END IF
266 *
267  IF( lsame( job, 'N' ) ) THEN
268  ilo = 1
269  ihi = n
270  DO 10 i = 1, n
271  lscale( i ) = one
272  rscale( i ) = one
273  10 continue
274  return
275  END IF
276 *
277  k = 1
278  l = n
279  IF( lsame( job, 'S' ) )
280  $ go to 190
281 *
282  go to 30
283 *
284 * Permute the matrices A and B to isolate the eigenvalues.
285 *
286 * Find row with one nonzero in columns 1 through L
287 *
288  20 continue
289  l = lm1
290  IF( l.NE.1 )
291  $ go to 30
292 *
293  rscale( 1 ) = one
294  lscale( 1 ) = one
295  go to 190
296 *
297  30 continue
298  lm1 = l - 1
299  DO 80 i = l, 1, -1
300  DO 40 j = 1, lm1
301  jp1 = j + 1
302  IF( a( i, j ).NE.czero .OR. b( i, j ).NE.czero )
303  $ go to 50
304  40 continue
305  j = l
306  go to 70
307 *
308  50 continue
309  DO 60 j = jp1, l
310  IF( a( i, j ).NE.czero .OR. b( i, j ).NE.czero )
311  $ go to 80
312  60 continue
313  j = jp1 - 1
314 *
315  70 continue
316  m = l
317  iflow = 1
318  go to 160
319  80 continue
320  go to 100
321 *
322 * Find column with one nonzero in rows K through N
323 *
324  90 continue
325  k = k + 1
326 *
327  100 continue
328  DO 150 j = k, l
329  DO 110 i = k, lm1
330  ip1 = i + 1
331  IF( a( i, j ).NE.czero .OR. b( i, j ).NE.czero )
332  $ go to 120
333  110 continue
334  i = l
335  go to 140
336  120 continue
337  DO 130 i = ip1, l
338  IF( a( i, j ).NE.czero .OR. b( i, j ).NE.czero )
339  $ go to 150
340  130 continue
341  i = ip1 - 1
342  140 continue
343  m = k
344  iflow = 2
345  go to 160
346  150 continue
347  go to 190
348 *
349 * Permute rows M and I
350 *
351  160 continue
352  lscale( m ) = i
353  IF( i.EQ.m )
354  $ go to 170
355  CALL cswap( n-k+1, a( i, k ), lda, a( m, k ), lda )
356  CALL cswap( n-k+1, b( i, k ), ldb, b( m, k ), ldb )
357 *
358 * Permute columns M and J
359 *
360  170 continue
361  rscale( m ) = j
362  IF( j.EQ.m )
363  $ go to 180
364  CALL cswap( l, a( 1, j ), 1, a( 1, m ), 1 )
365  CALL cswap( l, b( 1, j ), 1, b( 1, m ), 1 )
366 *
367  180 continue
368  go to( 20, 90 )iflow
369 *
370  190 continue
371  ilo = k
372  ihi = l
373 *
374  IF( lsame( job, 'P' ) ) THEN
375  DO 195 i = ilo, ihi
376  lscale( i ) = one
377  rscale( i ) = one
378  195 continue
379  return
380  END IF
381 *
382  IF( ilo.EQ.ihi )
383  $ return
384 *
385 * Balance the submatrix in rows ILO to IHI.
386 *
387  nr = ihi - ilo + 1
388  DO 200 i = ilo, ihi
389  rscale( i ) = zero
390  lscale( i ) = zero
391 *
392  work( i ) = zero
393  work( i+n ) = zero
394  work( i+2*n ) = zero
395  work( i+3*n ) = zero
396  work( i+4*n ) = zero
397  work( i+5*n ) = zero
398  200 continue
399 *
400 * Compute right side vector in resulting linear equations
401 *
402  basl = log10( sclfac )
403  DO 240 i = ilo, ihi
404  DO 230 j = ilo, ihi
405  IF( a( i, j ).EQ.czero ) THEN
406  ta = zero
407  go to 210
408  END IF
409  ta = log10( cabs1( a( i, j ) ) ) / basl
410 *
411  210 continue
412  IF( b( i, j ).EQ.czero ) THEN
413  tb = zero
414  go to 220
415  END IF
416  tb = log10( cabs1( b( i, j ) ) ) / basl
417 *
418  220 continue
419  work( i+4*n ) = work( i+4*n ) - ta - tb
420  work( j+5*n ) = work( j+5*n ) - ta - tb
421  230 continue
422  240 continue
423 *
424  coef = one / REAL( 2*nr )
425  coef2 = coef*coef
426  coef5 = half*coef2
427  nrp2 = nr + 2
428  beta = zero
429  it = 1
430 *
431 * Start generalized conjugate gradient iteration
432 *
433  250 continue
434 *
435  gamma = sdot( nr, work( ilo+4*n ), 1, work( ilo+4*n ), 1 ) +
436  $ sdot( nr, work( ilo+5*n ), 1, work( ilo+5*n ), 1 )
437 *
438  ew = zero
439  ewc = zero
440  DO 260 i = ilo, ihi
441  ew = ew + work( i+4*n )
442  ewc = ewc + work( i+5*n )
443  260 continue
444 *
445  gamma = coef*gamma - coef2*( ew**2+ewc**2 ) - coef5*( ew-ewc )**2
446  IF( gamma.EQ.zero )
447  $ go to 350
448  IF( it.NE.1 )
449  $ beta = gamma / pgamma
450  t = coef5*( ewc-three*ew )
451  tc = coef5*( ew-three*ewc )
452 *
453  CALL sscal( nr, beta, work( ilo ), 1 )
454  CALL sscal( nr, beta, work( ilo+n ), 1 )
455 *
456  CALL saxpy( nr, coef, work( ilo+4*n ), 1, work( ilo+n ), 1 )
457  CALL saxpy( nr, coef, work( ilo+5*n ), 1, work( ilo ), 1 )
458 *
459  DO 270 i = ilo, ihi
460  work( i ) = work( i ) + tc
461  work( i+n ) = work( i+n ) + t
462  270 continue
463 *
464 * Apply matrix to vector
465 *
466  DO 300 i = ilo, ihi
467  kount = 0
468  sum = zero
469  DO 290 j = ilo, ihi
470  IF( a( i, j ).EQ.czero )
471  $ go to 280
472  kount = kount + 1
473  sum = sum + work( j )
474  280 continue
475  IF( b( i, j ).EQ.czero )
476  $ go to 290
477  kount = kount + 1
478  sum = sum + work( j )
479  290 continue
480  work( i+2*n ) = REAL( kount )*work( i+n ) + sum
481  300 continue
482 *
483  DO 330 j = ilo, ihi
484  kount = 0
485  sum = zero
486  DO 320 i = ilo, ihi
487  IF( a( i, j ).EQ.czero )
488  $ go to 310
489  kount = kount + 1
490  sum = sum + work( i+n )
491  310 continue
492  IF( b( i, j ).EQ.czero )
493  $ go to 320
494  kount = kount + 1
495  sum = sum + work( i+n )
496  320 continue
497  work( j+3*n ) = REAL( kount )*work( j ) + sum
498  330 continue
499 *
500  sum = sdot( nr, work( ilo+n ), 1, work( ilo+2*n ), 1 ) +
501  $ sdot( nr, work( ilo ), 1, work( ilo+3*n ), 1 )
502  alpha = gamma / sum
503 *
504 * Determine correction to current iteration
505 *
506  cmax = zero
507  DO 340 i = ilo, ihi
508  cor = alpha*work( i+n )
509  IF( abs( cor ).GT.cmax )
510  $ cmax = abs( cor )
511  lscale( i ) = lscale( i ) + cor
512  cor = alpha*work( i )
513  IF( abs( cor ).GT.cmax )
514  $ cmax = abs( cor )
515  rscale( i ) = rscale( i ) + cor
516  340 continue
517  IF( cmax.LT.half )
518  $ go to 350
519 *
520  CALL saxpy( nr, -alpha, work( ilo+2*n ), 1, work( ilo+4*n ), 1 )
521  CALL saxpy( nr, -alpha, work( ilo+3*n ), 1, work( ilo+5*n ), 1 )
522 *
523  pgamma = gamma
524  it = it + 1
525  IF( it.LE.nrp2 )
526  $ go to 250
527 *
528 * End generalized conjugate gradient iteration
529 *
530  350 continue
531  sfmin = slamch( 'S' )
532  sfmax = one / sfmin
533  lsfmin = int( log10( sfmin ) / basl+one )
534  lsfmax = int( log10( sfmax ) / basl )
535  DO 360 i = ilo, ihi
536  irab = icamax( n-ilo+1, a( i, ilo ), lda )
537  rab = abs( a( i, irab+ilo-1 ) )
538  irab = icamax( n-ilo+1, b( i, ilo ), ldb )
539  rab = max( rab, abs( b( i, irab+ilo-1 ) ) )
540  lrab = int( log10( rab+sfmin ) / basl+one )
541  ir = lscale( i ) + sign( half, lscale( i ) )
542  ir = min( max( ir, lsfmin ), lsfmax, lsfmax-lrab )
543  lscale( i ) = sclfac**ir
544  icab = icamax( ihi, a( 1, i ), 1 )
545  cab = abs( a( icab, i ) )
546  icab = icamax( ihi, b( 1, i ), 1 )
547  cab = max( cab, abs( b( icab, i ) ) )
548  lcab = int( log10( cab+sfmin ) / basl+one )
549  jc = rscale( i ) + sign( half, rscale( i ) )
550  jc = min( max( jc, lsfmin ), lsfmax, lsfmax-lcab )
551  rscale( i ) = sclfac**jc
552  360 continue
553 *
554 * Row scaling of matrices A and B
555 *
556  DO 370 i = ilo, ihi
557  CALL csscal( n-ilo+1, lscale( i ), a( i, ilo ), lda )
558  CALL csscal( n-ilo+1, lscale( i ), b( i, ilo ), ldb )
559  370 continue
560 *
561 * Column scaling of matrices A and B
562 *
563  DO 380 j = ilo, ihi
564  CALL csscal( ihi, rscale( j ), a( 1, j ), 1 )
565  CALL csscal( ihi, rscale( j ), b( 1, j ), 1 )
566  380 continue
567 *
568  return
569 *
570 * End of CGGBAL
571 *
572  END