LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
sgges3.f
Go to the documentation of this file.
1 *> \brief <b> SGGES3 computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices (blocked algorithm)</b>
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download SGGES3 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgges3.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgges3.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgges3.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE SGGES3( JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B,
22 * $ LDB, SDIM, ALPHAR, ALPHAI, BETA, VSL, LDVSL,
23 * $ VSR, LDVSR, WORK, LWORK, BWORK, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER JOBVSL, JOBVSR, SORT
27 * INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM
28 * ..
29 * .. Array Arguments ..
30 * LOGICAL BWORK( * )
31 * REAL A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
32 * $ B( LDB, * ), BETA( * ), VSL( LDVSL, * ),
33 * $ VSR( LDVSR, * ), WORK( * )
34 * ..
35 * .. Function Arguments ..
36 * LOGICAL SELCTG
37 * EXTERNAL SELCTG
38 * ..
39 *
40 *
41 *> \par Purpose:
42 * =============
43 *>
44 *> \verbatim
45 *>
46 *> SGGES3 computes for a pair of N-by-N real nonsymmetric matrices (A,B),
47 *> the generalized eigenvalues, the generalized real Schur form (S,T),
48 *> optionally, the left and/or right matrices of Schur vectors (VSL and
49 *> VSR). This gives the generalized Schur factorization
50 *>
51 *> (A,B) = ( (VSL)*S*(VSR)**T, (VSL)*T*(VSR)**T )
52 *>
53 *> Optionally, it also orders the eigenvalues so that a selected cluster
54 *> of eigenvalues appears in the leading diagonal blocks of the upper
55 *> quasi-triangular matrix S and the upper triangular matrix T.The
56 *> leading columns of VSL and VSR then form an orthonormal basis for the
57 *> corresponding left and right eigenspaces (deflating subspaces).
58 *>
59 *> (If only the generalized eigenvalues are needed, use the driver
60 *> SGGEV instead, which is faster.)
61 *>
62 *> A generalized eigenvalue for a pair of matrices (A,B) is a scalar w
63 *> or a ratio alpha/beta = w, such that A - w*B is singular. It is
64 *> usually represented as the pair (alpha,beta), as there is a
65 *> reasonable interpretation for beta=0 or both being zero.
66 *>
67 *> A pair of matrices (S,T) is in generalized real Schur form if T is
68 *> upper triangular with non-negative diagonal and S is block upper
69 *> triangular with 1-by-1 and 2-by-2 blocks. 1-by-1 blocks correspond
70 *> to real generalized eigenvalues, while 2-by-2 blocks of S will be
71 *> "standardized" by making the corresponding elements of T have the
72 *> form:
73 *> [ a 0 ]
74 *> [ 0 b ]
75 *>
76 *> and the pair of corresponding 2-by-2 blocks in S and T will have a
77 *> complex conjugate pair of generalized eigenvalues.
78 *>
79 *> \endverbatim
80 *
81 * Arguments:
82 * ==========
83 *
84 *> \param[in] JOBVSL
85 *> \verbatim
86 *> JOBVSL is CHARACTER*1
87 *> = 'N': do not compute the left Schur vectors;
88 *> = 'V': compute the left Schur vectors.
89 *> \endverbatim
90 *>
91 *> \param[in] JOBVSR
92 *> \verbatim
93 *> JOBVSR is CHARACTER*1
94 *> = 'N': do not compute the right Schur vectors;
95 *> = 'V': compute the right Schur vectors.
96 *> \endverbatim
97 *>
98 *> \param[in] SORT
99 *> \verbatim
100 *> SORT is CHARACTER*1
101 *> Specifies whether or not to order the eigenvalues on the
102 *> diagonal of the generalized Schur form.
103 *> = 'N': Eigenvalues are not ordered;
104 *> = 'S': Eigenvalues are ordered (see SELCTG);
105 *> \endverbatim
106 *>
107 *> \param[in] SELCTG
108 *> \verbatim
109 *> SELCTG is a LOGICAL FUNCTION of three REAL arguments
110 *> SELCTG must be declared EXTERNAL in the calling subroutine.
111 *> If SORT = 'N', SELCTG is not referenced.
112 *> If SORT = 'S', SELCTG is used to select eigenvalues to sort
113 *> to the top left of the Schur form.
114 *> An eigenvalue (ALPHAR(j)+ALPHAI(j))/BETA(j) is selected if
115 *> SELCTG(ALPHAR(j),ALPHAI(j),BETA(j)) is true; i.e. if either
116 *> one of a complex conjugate pair of eigenvalues is selected,
117 *> then both complex eigenvalues are selected.
118 *>
119 *> Note that in the ill-conditioned case, a selected complex
120 *> eigenvalue may no longer satisfy SELCTG(ALPHAR(j),ALPHAI(j),
121 *> BETA(j)) = .TRUE. after ordering. INFO is to be set to N+2
122 *> in this case.
123 *> \endverbatim
124 *>
125 *> \param[in] N
126 *> \verbatim
127 *> N is INTEGER
128 *> The order of the matrices A, B, VSL, and VSR. N >= 0.
129 *> \endverbatim
130 *>
131 *> \param[in,out] A
132 *> \verbatim
133 *> A is REAL array, dimension (LDA, N)
134 *> On entry, the first of the pair of matrices.
135 *> On exit, A has been overwritten by its generalized Schur
136 *> form S.
137 *> \endverbatim
138 *>
139 *> \param[in] LDA
140 *> \verbatim
141 *> LDA is INTEGER
142 *> The leading dimension of A. LDA >= max(1,N).
143 *> \endverbatim
144 *>
145 *> \param[in,out] B
146 *> \verbatim
147 *> B is REAL array, dimension (LDB, N)
148 *> On entry, the second of the pair of matrices.
149 *> On exit, B has been overwritten by its generalized Schur
150 *> form T.
151 *> \endverbatim
152 *>
153 *> \param[in] LDB
154 *> \verbatim
155 *> LDB is INTEGER
156 *> The leading dimension of B. LDB >= max(1,N).
157 *> \endverbatim
158 *>
159 *> \param[out] SDIM
160 *> \verbatim
161 *> SDIM is INTEGER
162 *> If SORT = 'N', SDIM = 0.
163 *> If SORT = 'S', SDIM = number of eigenvalues (after sorting)
164 *> for which SELCTG is true. (Complex conjugate pairs for which
165 *> SELCTG is true for either eigenvalue count as 2.)
166 *> \endverbatim
167 *>
168 *> \param[out] ALPHAR
169 *> \verbatim
170 *> ALPHAR is REAL array, dimension (N)
171 *> \endverbatim
172 *>
173 *> \param[out] ALPHAI
174 *> \verbatim
175 *> ALPHAI is REAL array, dimension (N)
176 *> \endverbatim
177 *>
178 *> \param[out] BETA
179 *> \verbatim
180 *> BETA is REAL array, dimension (N)
181 *> On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will
182 *> be the generalized eigenvalues. ALPHAR(j) + ALPHAI(j)*i,
183 *> and BETA(j),j=1,...,N are the diagonals of the complex Schur
184 *> form (S,T) that would result if the 2-by-2 diagonal blocks of
185 *> the real Schur form of (A,B) were further reduced to
186 *> triangular form using 2-by-2 complex unitary transformations.
187 *> If ALPHAI(j) is zero, then the j-th eigenvalue is real; if
188 *> positive, then the j-th and (j+1)-st eigenvalues are a
189 *> complex conjugate pair, with ALPHAI(j+1) negative.
190 *>
191 *> Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j)
192 *> may easily over- or underflow, and BETA(j) may even be zero.
193 *> Thus, the user should avoid naively computing the ratio.
194 *> However, ALPHAR and ALPHAI will be always less than and
195 *> usually comparable with norm(A) in magnitude, and BETA always
196 *> less than and usually comparable with norm(B).
197 *> \endverbatim
198 *>
199 *> \param[out] VSL
200 *> \verbatim
201 *> VSL is REAL array, dimension (LDVSL,N)
202 *> If JOBVSL = 'V', VSL will contain the left Schur vectors.
203 *> Not referenced if JOBVSL = 'N'.
204 *> \endverbatim
205 *>
206 *> \param[in] LDVSL
207 *> \verbatim
208 *> LDVSL is INTEGER
209 *> The leading dimension of the matrix VSL. LDVSL >=1, and
210 *> if JOBVSL = 'V', LDVSL >= N.
211 *> \endverbatim
212 *>
213 *> \param[out] VSR
214 *> \verbatim
215 *> VSR is REAL array, dimension (LDVSR,N)
216 *> If JOBVSR = 'V', VSR will contain the right Schur vectors.
217 *> Not referenced if JOBVSR = 'N'.
218 *> \endverbatim
219 *>
220 *> \param[in] LDVSR
221 *> \verbatim
222 *> LDVSR is INTEGER
223 *> The leading dimension of the matrix VSR. LDVSR >= 1, and
224 *> if JOBVSR = 'V', LDVSR >= N.
225 *> \endverbatim
226 *>
227 *> \param[out] WORK
228 *> \verbatim
229 *> WORK is REAL array, dimension (MAX(1,LWORK))
230 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
231 *> \endverbatim
232 *>
233 *> \param[in] LWORK
234 *> \verbatim
235 *> LWORK is INTEGER
236 *> The dimension of the array WORK.
237 *>
238 *> If LWORK = -1, then a workspace query is assumed; the routine
239 *> only calculates the optimal size of the WORK array, returns
240 *> this value as the first entry of the WORK array, and no error
241 *> message related to LWORK is issued by XERBLA.
242 *> \endverbatim
243 *>
244 *> \param[out] BWORK
245 *> \verbatim
246 *> BWORK is LOGICAL array, dimension (N)
247 *> Not referenced if SORT = 'N'.
248 *> \endverbatim
249 *>
250 *> \param[out] INFO
251 *> \verbatim
252 *> INFO is INTEGER
253 *> = 0: successful exit
254 *> < 0: if INFO = -i, the i-th argument had an illegal value.
255 *> = 1,...,N:
256 *> The QZ iteration failed. (A,B) are not in Schur
257 *> form, but ALPHAR(j), ALPHAI(j), and BETA(j) should
258 *> be correct for j=INFO+1,...,N.
259 *> > N: =N+1: other than QZ iteration failed in SHGEQZ.
260 *> =N+2: after reordering, roundoff changed values of
261 *> some complex eigenvalues so that leading
262 *> eigenvalues in the Generalized Schur form no
263 *> longer satisfy SELCTG=.TRUE. This could also
264 *> be caused due to scaling.
265 *> =N+3: reordering failed in STGSEN.
266 *> \endverbatim
267 *
268 * Authors:
269 * ========
270 *
271 *> \author Univ. of Tennessee
272 *> \author Univ. of California Berkeley
273 *> \author Univ. of Colorado Denver
274 *> \author NAG Ltd.
275 *
276 *> \date January 2015
277 *
278 *> \ingroup realGEeigen
279 *
280 * =====================================================================
281  SUBROUTINE sgges3( JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B,
282  $ ldb, sdim, alphar, alphai, beta, vsl, ldvsl,
283  $ vsr, ldvsr, work, lwork, bwork, info )
284 *
285 * -- LAPACK driver routine (version 3.6.0) --
286 * -- LAPACK is a software package provided by Univ. of Tennessee, --
287 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
288 * January 2015
289 *
290 * .. Scalar Arguments ..
291  CHARACTER JOBVSL, JOBVSR, SORT
292  INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM
293 * ..
294 * .. Array Arguments ..
295  LOGICAL BWORK( * )
296  REAL A( lda, * ), ALPHAI( * ), ALPHAR( * ),
297  $ b( ldb, * ), beta( * ), vsl( ldvsl, * ),
298  $ vsr( ldvsr, * ), work( * )
299 * ..
300 * .. Function Arguments ..
301  LOGICAL SELCTG
302  EXTERNAL selctg
303 * ..
304 *
305 * =====================================================================
306 *
307 * .. Parameters ..
308  REAL ZERO, ONE
309  parameter ( zero = 0.0e+0, one = 1.0e+0 )
310 * ..
311 * .. Local Scalars ..
312  LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
313  $ lquery, lst2sl, wantst
314  INTEGER I, ICOLS, IERR, IHI, IJOBVL, IJOBVR, ILEFT,
315  $ ilo, ip, iright, irows, itau, iwrk, lwkopt
316  REAL ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PVSL,
317  $ pvsr, safmax, safmin, smlnum
318 * ..
319 * .. Local Arrays ..
320  INTEGER IDUM( 1 )
321  REAL DIF( 2 )
322 * ..
323 * .. External Subroutines ..
324  EXTERNAL sgeqrf, sggbak, sggbal, sgghd3, shgeqz, slabad,
326  $ xerbla
327 * ..
328 * .. External Functions ..
329  LOGICAL LSAME
330  REAL SLAMCH, SLANGE
331  EXTERNAL lsame, slamch, slange
332 * ..
333 * .. Intrinsic Functions ..
334  INTRINSIC abs, max, sqrt
335 * ..
336 * .. Executable Statements ..
337 *
338 * Decode the input arguments
339 *
340  IF( lsame( jobvsl, 'N' ) ) THEN
341  ijobvl = 1
342  ilvsl = .false.
343  ELSE IF( lsame( jobvsl, 'V' ) ) THEN
344  ijobvl = 2
345  ilvsl = .true.
346  ELSE
347  ijobvl = -1
348  ilvsl = .false.
349  END IF
350 *
351  IF( lsame( jobvsr, 'N' ) ) THEN
352  ijobvr = 1
353  ilvsr = .false.
354  ELSE IF( lsame( jobvsr, 'V' ) ) THEN
355  ijobvr = 2
356  ilvsr = .true.
357  ELSE
358  ijobvr = -1
359  ilvsr = .false.
360  END IF
361 *
362  wantst = lsame( sort, 'S' )
363 *
364 * Test the input arguments
365 *
366  info = 0
367  lquery = ( lwork.EQ.-1 )
368  IF( ijobvl.LE.0 ) THEN
369  info = -1
370  ELSE IF( ijobvr.LE.0 ) THEN
371  info = -2
372  ELSE IF( ( .NOT.wantst ) .AND. ( .NOT.lsame( sort, 'N' ) ) ) THEN
373  info = -3
374  ELSE IF( n.LT.0 ) THEN
375  info = -5
376  ELSE IF( lda.LT.max( 1, n ) ) THEN
377  info = -7
378  ELSE IF( ldb.LT.max( 1, n ) ) THEN
379  info = -9
380  ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
381  info = -15
382  ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
383  info = -17
384  ELSE IF( lwork.LT.6*n+16 .AND. .NOT.lquery ) THEN
385  info = -19
386  END IF
387 *
388 * Compute workspace
389 *
390  IF( info.EQ.0 ) THEN
391  CALL sgeqrf( n, n, b, ldb, work, work, -1, ierr )
392  lwkopt = max( 6*n+16, 3*n+int( work( 1 ) ) )
393  CALL sormqr( 'L', 'T', n, n, n, b, ldb, work, a, lda, work,
394  $ -1, ierr )
395  lwkopt = max( lwkopt, 3*n+int( work( 1 ) ) )
396  IF( ilvsl ) THEN
397  CALL sorgqr( n, n, n, vsl, ldvsl, work, work, -1, ierr )
398  lwkopt = max( lwkopt, 3*n+int( work( 1 ) ) )
399  END IF
400  CALL sgghd3( jobvsl, jobvsr, n, 1, n, a, lda, b, ldb, vsl,
401  $ ldvsl, vsr, ldvsr, work, -1, ierr )
402  lwkopt = max( lwkopt, 3*n+int( work( 1 ) ) )
403  CALL shgeqz( 'S', jobvsl, jobvsr, n, 1, n, a, lda, b, ldb,
404  $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
405  $ work, -1, ierr )
406  lwkopt = max( lwkopt, 2*n+int( work( 1 ) ) )
407  IF( wantst ) THEN
408  CALL stgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb,
409  $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
410  $ sdim, pvsl, pvsr, dif, work, -1, idum, 1,
411  $ ierr )
412  lwkopt = max( lwkopt, 2*n+int( work( 1 ) ) )
413  END IF
414  work( 1 ) = lwkopt
415  END IF
416 *
417  IF( info.NE.0 ) THEN
418  CALL xerbla( 'SGGES3 ', -info )
419  RETURN
420  ELSE IF( lquery ) THEN
421  RETURN
422  END IF
423 *
424 * Quick return if possible
425 *
426  IF( n.EQ.0 ) THEN
427  sdim = 0
428  RETURN
429  END IF
430 *
431 * Get machine constants
432 *
433  eps = slamch( 'P' )
434  safmin = slamch( 'S' )
435  safmax = one / safmin
436  CALL slabad( safmin, safmax )
437  smlnum = sqrt( safmin ) / eps
438  bignum = one / smlnum
439 *
440 * Scale A if max element outside range [SMLNUM,BIGNUM]
441 *
442  anrm = slange( 'M', n, n, a, lda, work )
443  ilascl = .false.
444  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
445  anrmto = smlnum
446  ilascl = .true.
447  ELSE IF( anrm.GT.bignum ) THEN
448  anrmto = bignum
449  ilascl = .true.
450  END IF
451  IF( ilascl )
452  $ CALL slascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
453 *
454 * Scale B if max element outside range [SMLNUM,BIGNUM]
455 *
456  bnrm = slange( 'M', n, n, b, ldb, work )
457  ilbscl = .false.
458  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
459  bnrmto = smlnum
460  ilbscl = .true.
461  ELSE IF( bnrm.GT.bignum ) THEN
462  bnrmto = bignum
463  ilbscl = .true.
464  END IF
465  IF( ilbscl )
466  $ CALL slascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
467 *
468 * Permute the matrix to make it more nearly triangular
469 *
470  ileft = 1
471  iright = n + 1
472  iwrk = iright + n
473  CALL sggbal( 'P', n, a, lda, b, ldb, ilo, ihi, work( ileft ),
474  $ work( iright ), work( iwrk ), ierr )
475 *
476 * Reduce B to triangular form (QR decomposition of B)
477 *
478  irows = ihi + 1 - ilo
479  icols = n + 1 - ilo
480  itau = iwrk
481  iwrk = itau + irows
482  CALL sgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
483  $ work( iwrk ), lwork+1-iwrk, ierr )
484 *
485 * Apply the orthogonal transformation to matrix A
486 *
487  CALL sormqr( 'L', 'T', irows, icols, irows, b( ilo, ilo ), ldb,
488  $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
489  $ lwork+1-iwrk, ierr )
490 *
491 * Initialize VSL
492 *
493  IF( ilvsl ) THEN
494  CALL slaset( 'Full', n, n, zero, one, vsl, ldvsl )
495  IF( irows.GT.1 ) THEN
496  CALL slacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
497  $ vsl( ilo+1, ilo ), ldvsl )
498  END IF
499  CALL sorgqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
500  $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
501  END IF
502 *
503 * Initialize VSR
504 *
505  IF( ilvsr )
506  $ CALL slaset( 'Full', n, n, zero, one, vsr, ldvsr )
507 *
508 * Reduce to generalized Hessenberg form
509 *
510  CALL sgghd3( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
511  $ ldvsl, vsr, ldvsr, work( iwrk ), lwork+1-iwrk, ierr )
512 *
513 * Perform QZ algorithm, computing Schur vectors if desired
514 *
515  iwrk = itau
516  CALL shgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
517  $ alphar, alphai, beta, vsl, ldvsl, vsr, ldvsr,
518  $ work( iwrk ), lwork+1-iwrk, ierr )
519  IF( ierr.NE.0 ) THEN
520  IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
521  info = ierr
522  ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
523  info = ierr - n
524  ELSE
525  info = n + 1
526  END IF
527  GO TO 40
528  END IF
529 *
530 * Sort eigenvalues ALPHA/BETA if desired
531 *
532  sdim = 0
533  IF( wantst ) THEN
534 *
535 * Undo scaling on eigenvalues before SELCTGing
536 *
537  IF( ilascl ) THEN
538  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n,
539  $ ierr )
540  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n,
541  $ ierr )
542  END IF
543  IF( ilbscl )
544  $ CALL slascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
545 *
546 * Select eigenvalues
547 *
548  DO 10 i = 1, n
549  bwork( i ) = selctg( alphar( i ), alphai( i ), beta( i ) )
550  10 CONTINUE
551 *
552  CALL stgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb, alphar,
553  $ alphai, beta, vsl, ldvsl, vsr, ldvsr, sdim, pvsl,
554  $ pvsr, dif, work( iwrk ), lwork-iwrk+1, idum, 1,
555  $ ierr )
556  IF( ierr.EQ.1 )
557  $ info = n + 3
558 *
559  END IF
560 *
561 * Apply back-permutation to VSL and VSR
562 *
563  IF( ilvsl )
564  $ CALL sggbak( 'P', 'L', n, ilo, ihi, work( ileft ),
565  $ work( iright ), n, vsl, ldvsl, ierr )
566 *
567  IF( ilvsr )
568  $ CALL sggbak( 'P', 'R', n, ilo, ihi, work( ileft ),
569  $ work( iright ), n, vsr, ldvsr, ierr )
570 *
571 * Check if unscaling would cause over/underflow, if so, rescale
572 * (ALPHAR(I),ALPHAI(I),BETA(I)) so BETA(I) is on the order of
573 * B(I,I) and ALPHAR(I) and ALPHAI(I) are on the order of A(I,I)
574 *
575  IF( ilascl )THEN
576  DO 50 i = 1, n
577  IF( alphai( i ).NE.zero ) THEN
578  IF( ( alphar( i )/safmax ).GT.( anrmto/anrm ) .OR.
579  $ ( safmin/alphar( i ) ).GT.( anrm/anrmto ) ) THEN
580  work( 1 ) = abs( a( i, i )/alphar( i ) )
581  beta( i ) = beta( i )*work( 1 )
582  alphar( i ) = alphar( i )*work( 1 )
583  alphai( i ) = alphai( i )*work( 1 )
584  ELSE IF( ( alphai( i )/safmax ).GT.( anrmto/anrm ) .OR.
585  $ ( safmin/alphai( i ) ).GT.( anrm/anrmto ) ) THEN
586  work( 1 ) = abs( a( i, i+1 )/alphai( i ) )
587  beta( i ) = beta( i )*work( 1 )
588  alphar( i ) = alphar( i )*work( 1 )
589  alphai( i ) = alphai( i )*work( 1 )
590  END IF
591  END IF
592  50 CONTINUE
593  END IF
594 *
595  IF( ilbscl )THEN
596  DO 60 i = 1, n
597  IF( alphai( i ).NE.zero ) THEN
598  IF( ( beta( i )/safmax ).GT.( bnrmto/bnrm ) .OR.
599  $ ( safmin/beta( i ) ).GT.( bnrm/bnrmto ) ) THEN
600  work( 1 ) = abs(b( i, i )/beta( i ))
601  beta( i ) = beta( i )*work( 1 )
602  alphar( i ) = alphar( i )*work( 1 )
603  alphai( i ) = alphai( i )*work( 1 )
604  END IF
605  END IF
606  60 CONTINUE
607  END IF
608 *
609 * Undo scaling
610 *
611  IF( ilascl ) THEN
612  CALL slascl( 'H', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
613  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n, ierr )
614  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n, ierr )
615  END IF
616 *
617  IF( ilbscl ) THEN
618  CALL slascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
619  CALL slascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
620  END IF
621 *
622  IF( wantst ) THEN
623 *
624 * Check if reordering is correct
625 *
626  lastsl = .true.
627  lst2sl = .true.
628  sdim = 0
629  ip = 0
630  DO 30 i = 1, n
631  cursl = selctg( alphar( i ), alphai( i ), beta( i ) )
632  IF( alphai( i ).EQ.zero ) THEN
633  IF( cursl )
634  $ sdim = sdim + 1
635  ip = 0
636  IF( cursl .AND. .NOT.lastsl )
637  $ info = n + 2
638  ELSE
639  IF( ip.EQ.1 ) THEN
640 *
641 * Last eigenvalue of conjugate pair
642 *
643  cursl = cursl .OR. lastsl
644  lastsl = cursl
645  IF( cursl )
646  $ sdim = sdim + 2
647  ip = -1
648  IF( cursl .AND. .NOT.lst2sl )
649  $ info = n + 2
650  ELSE
651 *
652 * First eigenvalue of conjugate pair
653 *
654  ip = 1
655  END IF
656  END IF
657  lst2sl = lastsl
658  lastsl = cursl
659  30 CONTINUE
660 *
661  END IF
662 *
663  40 CONTINUE
664 *
665  work( 1 ) = lwkopt
666 *
667  RETURN
668 *
669 * End of SGGES3
670 *
671  END
subroutine sggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
SGGBAL
Definition: sggbal.f:179
subroutine stgsen(IJOB, WANTQ, WANTZ, SELECT, N, A, LDA, B, LDB, ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, M, PL, PR, DIF, WORK, LWORK, IWORK, LIWORK, INFO)
STGSEN
Definition: stgsen.f:453
subroutine sormqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
SORMQR
Definition: sormqr.f:170
subroutine slascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
SLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: slascl.f:145
subroutine slabad(SMALL, LARGE)
SLABAD
Definition: slabad.f:76
subroutine sgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
SGEQRF
Definition: sgeqrf.f:138
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine sgghd3(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, WORK, LWORK, INFO)
SGGHD3
Definition: sgghd3.f:232
subroutine shgeqz(JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, INFO)
SHGEQZ
Definition: shgeqz.f:306
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
subroutine sgges3(JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B, LDB, SDIM, ALPHAR, ALPHAI, BETA, VSL, LDVSL, VSR, LDVSR, WORK, LWORK, BWORK, INFO)
SGGES3 computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE...
Definition: sgges3.f:284
subroutine sggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
SGGBAK
Definition: sggbak.f:149
subroutine sorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
SORGQR
Definition: sorgqr.f:130