LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
sggev.f
Go to the documentation of this file.
1 *> \brief <b> SGGEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices</b>
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download SGGEV + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sggev.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sggev.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sggev.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE SGGEV( JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI,
22 * BETA, VL, LDVL, VR, LDVR, WORK, LWORK, INFO )
23 *
24 * .. Scalar Arguments ..
25 * CHARACTER JOBVL, JOBVR
26 * INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, N
27 * ..
28 * .. Array Arguments ..
29 * REAL A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
30 * $ B( LDB, * ), BETA( * ), VL( LDVL, * ),
31 * $ VR( LDVR, * ), WORK( * )
32 * ..
33 *
34 *
35 *> \par Purpose:
36 * =============
37 *>
38 *> \verbatim
39 *>
40 *> SGGEV computes for a pair of N-by-N real nonsymmetric matrices (A,B)
41 *> the generalized eigenvalues, and optionally, the left and/or right
42 *> generalized eigenvectors.
43 *>
44 *> A generalized eigenvalue for a pair of matrices (A,B) is a scalar
45 *> lambda or a ratio alpha/beta = lambda, such that A - lambda*B is
46 *> singular. It is usually represented as the pair (alpha,beta), as
47 *> there is a reasonable interpretation for beta=0, and even for both
48 *> being zero.
49 *>
50 *> The right eigenvector v(j) corresponding to the eigenvalue lambda(j)
51 *> of (A,B) satisfies
52 *>
53 *> A * v(j) = lambda(j) * B * v(j).
54 *>
55 *> The left eigenvector u(j) corresponding to the eigenvalue lambda(j)
56 *> of (A,B) satisfies
57 *>
58 *> u(j)**H * A = lambda(j) * u(j)**H * B .
59 *>
60 *> where u(j)**H is the conjugate-transpose of u(j).
61 *>
62 *> \endverbatim
63 *
64 * Arguments:
65 * ==========
66 *
67 *> \param[in] JOBVL
68 *> \verbatim
69 *> JOBVL is CHARACTER*1
70 *> = 'N': do not compute the left generalized eigenvectors;
71 *> = 'V': compute the left generalized eigenvectors.
72 *> \endverbatim
73 *>
74 *> \param[in] JOBVR
75 *> \verbatim
76 *> JOBVR is CHARACTER*1
77 *> = 'N': do not compute the right generalized eigenvectors;
78 *> = 'V': compute the right generalized eigenvectors.
79 *> \endverbatim
80 *>
81 *> \param[in] N
82 *> \verbatim
83 *> N is INTEGER
84 *> The order of the matrices A, B, VL, and VR. N >= 0.
85 *> \endverbatim
86 *>
87 *> \param[in,out] A
88 *> \verbatim
89 *> A is REAL array, dimension (LDA, N)
90 *> On entry, the matrix A in the pair (A,B).
91 *> On exit, A has been overwritten.
92 *> \endverbatim
93 *>
94 *> \param[in] LDA
95 *> \verbatim
96 *> LDA is INTEGER
97 *> The leading dimension of A. LDA >= max(1,N).
98 *> \endverbatim
99 *>
100 *> \param[in,out] B
101 *> \verbatim
102 *> B is REAL array, dimension (LDB, N)
103 *> On entry, the matrix B in the pair (A,B).
104 *> On exit, B has been overwritten.
105 *> \endverbatim
106 *>
107 *> \param[in] LDB
108 *> \verbatim
109 *> LDB is INTEGER
110 *> The leading dimension of B. LDB >= max(1,N).
111 *> \endverbatim
112 *>
113 *> \param[out] ALPHAR
114 *> \verbatim
115 *> ALPHAR is REAL array, dimension (N)
116 *> \endverbatim
117 *>
118 *> \param[out] ALPHAI
119 *> \verbatim
120 *> ALPHAI is REAL array, dimension (N)
121 *> \endverbatim
122 *>
123 *> \param[out] BETA
124 *> \verbatim
125 *> BETA is REAL array, dimension (N)
126 *> On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will
127 *> be the generalized eigenvalues. If ALPHAI(j) is zero, then
128 *> the j-th eigenvalue is real; if positive, then the j-th and
129 *> (j+1)-st eigenvalues are a complex conjugate pair, with
130 *> ALPHAI(j+1) negative.
131 *>
132 *> Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j)
133 *> may easily over- or underflow, and BETA(j) may even be zero.
134 *> Thus, the user should avoid naively computing the ratio
135 *> alpha/beta. However, ALPHAR and ALPHAI will be always less
136 *> than and usually comparable with norm(A) in magnitude, and
137 *> BETA always less than and usually comparable with norm(B).
138 *> \endverbatim
139 *>
140 *> \param[out] VL
141 *> \verbatim
142 *> VL is REAL array, dimension (LDVL,N)
143 *> If JOBVL = 'V', the left eigenvectors u(j) are stored one
144 *> after another in the columns of VL, in the same order as
145 *> their eigenvalues. If the j-th eigenvalue is real, then
146 *> u(j) = VL(:,j), the j-th column of VL. If the j-th and
147 *> (j+1)-th eigenvalues form a complex conjugate pair, then
148 *> u(j) = VL(:,j)+i*VL(:,j+1) and u(j+1) = VL(:,j)-i*VL(:,j+1).
149 *> Each eigenvector is scaled so the largest component has
150 *> abs(real part)+abs(imag. part)=1.
151 *> Not referenced if JOBVL = 'N'.
152 *> \endverbatim
153 *>
154 *> \param[in] LDVL
155 *> \verbatim
156 *> LDVL is INTEGER
157 *> The leading dimension of the matrix VL. LDVL >= 1, and
158 *> if JOBVL = 'V', LDVL >= N.
159 *> \endverbatim
160 *>
161 *> \param[out] VR
162 *> \verbatim
163 *> VR is REAL array, dimension (LDVR,N)
164 *> If JOBVR = 'V', the right eigenvectors v(j) are stored one
165 *> after another in the columns of VR, in the same order as
166 *> their eigenvalues. If the j-th eigenvalue is real, then
167 *> v(j) = VR(:,j), the j-th column of VR. If the j-th and
168 *> (j+1)-th eigenvalues form a complex conjugate pair, then
169 *> v(j) = VR(:,j)+i*VR(:,j+1) and v(j+1) = VR(:,j)-i*VR(:,j+1).
170 *> Each eigenvector is scaled so the largest component has
171 *> abs(real part)+abs(imag. part)=1.
172 *> Not referenced if JOBVR = 'N'.
173 *> \endverbatim
174 *>
175 *> \param[in] LDVR
176 *> \verbatim
177 *> LDVR is INTEGER
178 *> The leading dimension of the matrix VR. LDVR >= 1, and
179 *> if JOBVR = 'V', LDVR >= N.
180 *> \endverbatim
181 *>
182 *> \param[out] WORK
183 *> \verbatim
184 *> WORK is REAL array, dimension (MAX(1,LWORK))
185 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
186 *> \endverbatim
187 *>
188 *> \param[in] LWORK
189 *> \verbatim
190 *> LWORK is INTEGER
191 *> The dimension of the array WORK. LWORK >= max(1,8*N).
192 *> For good performance, LWORK must generally be larger.
193 *>
194 *> If LWORK = -1, then a workspace query is assumed; the routine
195 *> only calculates the optimal size of the WORK array, returns
196 *> this value as the first entry of the WORK array, and no error
197 *> message related to LWORK is issued by XERBLA.
198 *> \endverbatim
199 *>
200 *> \param[out] INFO
201 *> \verbatim
202 *> INFO is INTEGER
203 *> = 0: successful exit
204 *> < 0: if INFO = -i, the i-th argument had an illegal value.
205 *> = 1,...,N:
206 *> The QZ iteration failed. No eigenvectors have been
207 *> calculated, but ALPHAR(j), ALPHAI(j), and BETA(j)
208 *> should be correct for j=INFO+1,...,N.
209 *> > N: =N+1: other than QZ iteration failed in SHGEQZ.
210 *> =N+2: error return from STGEVC.
211 *> \endverbatim
212 *
213 * Authors:
214 * ========
215 *
216 *> \author Univ. of Tennessee
217 *> \author Univ. of California Berkeley
218 *> \author Univ. of Colorado Denver
219 *> \author NAG Ltd.
220 *
221 *> \date April 2012
222 *
223 *> \ingroup realGEeigen
224 *
225 * =====================================================================
226  SUBROUTINE sggev( JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI,
227  $ beta, vl, ldvl, vr, ldvr, work, lwork, info )
228 *
229 * -- LAPACK driver routine (version 3.4.1) --
230 * -- LAPACK is a software package provided by Univ. of Tennessee, --
231 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
232 * April 2012
233 *
234 * .. Scalar Arguments ..
235  CHARACTER JOBVL, JOBVR
236  INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, N
237 * ..
238 * .. Array Arguments ..
239  REAL A( lda, * ), ALPHAI( * ), ALPHAR( * ),
240  $ b( ldb, * ), beta( * ), vl( ldvl, * ),
241  $ vr( ldvr, * ), work( * )
242 * ..
243 *
244 * =====================================================================
245 *
246 * .. Parameters ..
247  REAL ZERO, ONE
248  parameter ( zero = 0.0e+0, one = 1.0e+0 )
249 * ..
250 * .. Local Scalars ..
251  LOGICAL ILASCL, ILBSCL, ILV, ILVL, ILVR, LQUERY
252  CHARACTER CHTEMP
253  INTEGER ICOLS, IERR, IHI, IJOBVL, IJOBVR, ILEFT, ILO,
254  $ in, iright, irows, itau, iwrk, jc, jr, maxwrk,
255  $ minwrk
256  REAL ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS,
257  $ smlnum, temp
258 * ..
259 * .. Local Arrays ..
260  LOGICAL LDUMMA( 1 )
261 * ..
262 * .. External Subroutines ..
263  EXTERNAL sgeqrf, sggbak, sggbal, sgghrd, shgeqz, slabad,
265  $ xerbla
266 * ..
267 * .. External Functions ..
268  LOGICAL LSAME
269  INTEGER ILAENV
270  REAL SLAMCH, SLANGE
271  EXTERNAL lsame, ilaenv, slamch, slange
272 * ..
273 * .. Intrinsic Functions ..
274  INTRINSIC abs, max, sqrt
275 * ..
276 * .. Executable Statements ..
277 *
278 * Decode the input arguments
279 *
280  IF( lsame( jobvl, 'N' ) ) THEN
281  ijobvl = 1
282  ilvl = .false.
283  ELSE IF( lsame( jobvl, 'V' ) ) THEN
284  ijobvl = 2
285  ilvl = .true.
286  ELSE
287  ijobvl = -1
288  ilvl = .false.
289  END IF
290 *
291  IF( lsame( jobvr, 'N' ) ) THEN
292  ijobvr = 1
293  ilvr = .false.
294  ELSE IF( lsame( jobvr, 'V' ) ) THEN
295  ijobvr = 2
296  ilvr = .true.
297  ELSE
298  ijobvr = -1
299  ilvr = .false.
300  END IF
301  ilv = ilvl .OR. ilvr
302 *
303 * Test the input arguments
304 *
305  info = 0
306  lquery = ( lwork.EQ.-1 )
307  IF( ijobvl.LE.0 ) THEN
308  info = -1
309  ELSE IF( ijobvr.LE.0 ) THEN
310  info = -2
311  ELSE IF( n.LT.0 ) THEN
312  info = -3
313  ELSE IF( lda.LT.max( 1, n ) ) THEN
314  info = -5
315  ELSE IF( ldb.LT.max( 1, n ) ) THEN
316  info = -7
317  ELSE IF( ldvl.LT.1 .OR. ( ilvl .AND. ldvl.LT.n ) ) THEN
318  info = -12
319  ELSE IF( ldvr.LT.1 .OR. ( ilvr .AND. ldvr.LT.n ) ) THEN
320  info = -14
321  END IF
322 *
323 * Compute workspace
324 * (Note: Comments in the code beginning "Workspace:" describe the
325 * minimal amount of workspace needed at that point in the code,
326 * as well as the preferred amount for good performance.
327 * NB refers to the optimal block size for the immediately
328 * following subroutine, as returned by ILAENV. The workspace is
329 * computed assuming ILO = 1 and IHI = N, the worst case.)
330 *
331  IF( info.EQ.0 ) THEN
332  minwrk = max( 1, 8*n )
333  maxwrk = max( 1, n*( 7 +
334  $ ilaenv( 1, 'SGEQRF', ' ', n, 1, n, 0 ) ) )
335  maxwrk = max( maxwrk, n*( 7 +
336  $ ilaenv( 1, 'SORMQR', ' ', n, 1, n, 0 ) ) )
337  IF( ilvl ) THEN
338  maxwrk = max( maxwrk, n*( 7 +
339  $ ilaenv( 1, 'SORGQR', ' ', n, 1, n, -1 ) ) )
340  END IF
341  work( 1 ) = maxwrk
342 *
343  IF( lwork.LT.minwrk .AND. .NOT.lquery )
344  $ info = -16
345  END IF
346 *
347  IF( info.NE.0 ) THEN
348  CALL xerbla( 'SGGEV ', -info )
349  RETURN
350  ELSE IF( lquery ) THEN
351  RETURN
352  END IF
353 *
354 * Quick return if possible
355 *
356  IF( n.EQ.0 )
357  $ RETURN
358 *
359 * Get machine constants
360 *
361  eps = slamch( 'P' )
362  smlnum = slamch( 'S' )
363  bignum = one / smlnum
364  CALL slabad( smlnum, bignum )
365  smlnum = sqrt( smlnum ) / eps
366  bignum = one / smlnum
367 *
368 * Scale A if max element outside range [SMLNUM,BIGNUM]
369 *
370  anrm = slange( 'M', n, n, a, lda, work )
371  ilascl = .false.
372  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
373  anrmto = smlnum
374  ilascl = .true.
375  ELSE IF( anrm.GT.bignum ) THEN
376  anrmto = bignum
377  ilascl = .true.
378  END IF
379  IF( ilascl )
380  $ CALL slascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
381 *
382 * Scale B if max element outside range [SMLNUM,BIGNUM]
383 *
384  bnrm = slange( 'M', n, n, b, ldb, work )
385  ilbscl = .false.
386  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
387  bnrmto = smlnum
388  ilbscl = .true.
389  ELSE IF( bnrm.GT.bignum ) THEN
390  bnrmto = bignum
391  ilbscl = .true.
392  END IF
393  IF( ilbscl )
394  $ CALL slascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
395 *
396 * Permute the matrices A, B to isolate eigenvalues if possible
397 * (Workspace: need 6*N)
398 *
399  ileft = 1
400  iright = n + 1
401  iwrk = iright + n
402  CALL sggbal( 'P', n, a, lda, b, ldb, ilo, ihi, work( ileft ),
403  $ work( iright ), work( iwrk ), ierr )
404 *
405 * Reduce B to triangular form (QR decomposition of B)
406 * (Workspace: need N, prefer N*NB)
407 *
408  irows = ihi + 1 - ilo
409  IF( ilv ) THEN
410  icols = n + 1 - ilo
411  ELSE
412  icols = irows
413  END IF
414  itau = iwrk
415  iwrk = itau + irows
416  CALL sgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
417  $ work( iwrk ), lwork+1-iwrk, ierr )
418 *
419 * Apply the orthogonal transformation to matrix A
420 * (Workspace: need N, prefer N*NB)
421 *
422  CALL sormqr( 'L', 'T', irows, icols, irows, b( ilo, ilo ), ldb,
423  $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
424  $ lwork+1-iwrk, ierr )
425 *
426 * Initialize VL
427 * (Workspace: need N, prefer N*NB)
428 *
429  IF( ilvl ) THEN
430  CALL slaset( 'Full', n, n, zero, one, vl, ldvl )
431  IF( irows.GT.1 ) THEN
432  CALL slacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
433  $ vl( ilo+1, ilo ), ldvl )
434  END IF
435  CALL sorgqr( irows, irows, irows, vl( ilo, ilo ), ldvl,
436  $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
437  END IF
438 *
439 * Initialize VR
440 *
441  IF( ilvr )
442  $ CALL slaset( 'Full', n, n, zero, one, vr, ldvr )
443 *
444 * Reduce to generalized Hessenberg form
445 * (Workspace: none needed)
446 *
447  IF( ilv ) THEN
448 *
449 * Eigenvectors requested -- work on whole matrix.
450 *
451  CALL sgghrd( jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb, vl,
452  $ ldvl, vr, ldvr, ierr )
453  ELSE
454  CALL sgghrd( 'N', 'N', irows, 1, irows, a( ilo, ilo ), lda,
455  $ b( ilo, ilo ), ldb, vl, ldvl, vr, ldvr, ierr )
456  END IF
457 *
458 * Perform QZ algorithm (Compute eigenvalues, and optionally, the
459 * Schur forms and Schur vectors)
460 * (Workspace: need N)
461 *
462  iwrk = itau
463  IF( ilv ) THEN
464  chtemp = 'S'
465  ELSE
466  chtemp = 'E'
467  END IF
468  CALL shgeqz( chtemp, jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb,
469  $ alphar, alphai, beta, vl, ldvl, vr, ldvr,
470  $ work( iwrk ), lwork+1-iwrk, ierr )
471  IF( ierr.NE.0 ) THEN
472  IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
473  info = ierr
474  ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
475  info = ierr - n
476  ELSE
477  info = n + 1
478  END IF
479  GO TO 110
480  END IF
481 *
482 * Compute Eigenvectors
483 * (Workspace: need 6*N)
484 *
485  IF( ilv ) THEN
486  IF( ilvl ) THEN
487  IF( ilvr ) THEN
488  chtemp = 'B'
489  ELSE
490  chtemp = 'L'
491  END IF
492  ELSE
493  chtemp = 'R'
494  END IF
495  CALL stgevc( chtemp, 'B', ldumma, n, a, lda, b, ldb, vl, ldvl,
496  $ vr, ldvr, n, in, work( iwrk ), ierr )
497  IF( ierr.NE.0 ) THEN
498  info = n + 2
499  GO TO 110
500  END IF
501 *
502 * Undo balancing on VL and VR and normalization
503 * (Workspace: none needed)
504 *
505  IF( ilvl ) THEN
506  CALL sggbak( 'P', 'L', n, ilo, ihi, work( ileft ),
507  $ work( iright ), n, vl, ldvl, ierr )
508  DO 50 jc = 1, n
509  IF( alphai( jc ).LT.zero )
510  $ GO TO 50
511  temp = zero
512  IF( alphai( jc ).EQ.zero ) THEN
513  DO 10 jr = 1, n
514  temp = max( temp, abs( vl( jr, jc ) ) )
515  10 CONTINUE
516  ELSE
517  DO 20 jr = 1, n
518  temp = max( temp, abs( vl( jr, jc ) )+
519  $ abs( vl( jr, jc+1 ) ) )
520  20 CONTINUE
521  END IF
522  IF( temp.LT.smlnum )
523  $ GO TO 50
524  temp = one / temp
525  IF( alphai( jc ).EQ.zero ) THEN
526  DO 30 jr = 1, n
527  vl( jr, jc ) = vl( jr, jc )*temp
528  30 CONTINUE
529  ELSE
530  DO 40 jr = 1, n
531  vl( jr, jc ) = vl( jr, jc )*temp
532  vl( jr, jc+1 ) = vl( jr, jc+1 )*temp
533  40 CONTINUE
534  END IF
535  50 CONTINUE
536  END IF
537  IF( ilvr ) THEN
538  CALL sggbak( 'P', 'R', n, ilo, ihi, work( ileft ),
539  $ work( iright ), n, vr, ldvr, ierr )
540  DO 100 jc = 1, n
541  IF( alphai( jc ).LT.zero )
542  $ GO TO 100
543  temp = zero
544  IF( alphai( jc ).EQ.zero ) THEN
545  DO 60 jr = 1, n
546  temp = max( temp, abs( vr( jr, jc ) ) )
547  60 CONTINUE
548  ELSE
549  DO 70 jr = 1, n
550  temp = max( temp, abs( vr( jr, jc ) )+
551  $ abs( vr( jr, jc+1 ) ) )
552  70 CONTINUE
553  END IF
554  IF( temp.LT.smlnum )
555  $ GO TO 100
556  temp = one / temp
557  IF( alphai( jc ).EQ.zero ) THEN
558  DO 80 jr = 1, n
559  vr( jr, jc ) = vr( jr, jc )*temp
560  80 CONTINUE
561  ELSE
562  DO 90 jr = 1, n
563  vr( jr, jc ) = vr( jr, jc )*temp
564  vr( jr, jc+1 ) = vr( jr, jc+1 )*temp
565  90 CONTINUE
566  END IF
567  100 CONTINUE
568  END IF
569 *
570 * End of eigenvector calculation
571 *
572  END IF
573 *
574 * Undo scaling if necessary
575 *
576  110 CONTINUE
577 *
578  IF( ilascl ) THEN
579  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n, ierr )
580  CALL slascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n, ierr )
581  END IF
582 *
583  IF( ilbscl ) THEN
584  CALL slascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
585  END IF
586 *
587  work( 1 ) = maxwrk
588  RETURN
589 *
590 * End of SGGEV
591 *
592  END
subroutine sggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
SGGBAL
Definition: sggbal.f:179
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 stgevc(SIDE, HOWMNY, SELECT, N, S, LDS, P, LDP, VL, LDVL, VR, LDVR, MM, M, WORK, INFO)
STGEVC
Definition: stgevc.f:297
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 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 sggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
SGGBAK
Definition: sggbak.f:149
subroutine sggev(JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, WORK, LWORK, INFO)
SGGEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices ...
Definition: sggev.f:228
subroutine sgghrd(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, INFO)
SGGHRD
Definition: sgghrd.f:209
subroutine sorgqr(M, N, K, A, LDA, TAU, WORK, LWORK, INFO)
SORGQR
Definition: sorgqr.f:130