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