LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dggev.f
Go to the documentation of this file.
1*> \brief <b> DGGEV 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 DGGEV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dggev.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dggev.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dggev.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DGGEV( 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* DOUBLE PRECISION 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*> DGGEV 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
116*> \endverbatim
117*>
118*> \param[out] ALPHAI
119*> \verbatim
120*> ALPHAI is DOUBLE PRECISION array, dimension (N)
121*> \endverbatim
122*>
123*> \param[out] BETA
124*> \verbatim
125*> BETA is DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DHGEQZ.
210*> =N+2: error return from DTGEVC.
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*> \ingroup ggev
222*
223* =====================================================================
224 SUBROUTINE dggev( JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI,
225 $ BETA, VL, LDVL, VR, LDVR, WORK, LWORK, 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 DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
237 $ b( ldb, * ), beta( * ), vl( ldvl, * ),
238 $ vr( ldvr, * ), work( * )
239* ..
240*
241* =====================================================================
242*
243* .. Parameters ..
244 DOUBLE PRECISION ZERO, ONE
245 parameter( zero = 0.0d+0, one = 1.0d+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, maxwrk,
252 $ minwrk
253 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS,
254 $ smlnum, temp
255* ..
256* .. Local Arrays ..
257 LOGICAL LDUMMA( 1 )
258* ..
259* .. External Subroutines ..
260 EXTERNAL dgeqrf, dggbak, dggbal, dgghrd, dhgeqz, dlacpy,
262* ..
263* .. External Functions ..
264 LOGICAL LSAME
265 INTEGER ILAENV
266 DOUBLE PRECISION DLAMCH, DLANGE
267 EXTERNAL lsame, ilaenv, dlamch, dlange
268* ..
269* .. Intrinsic Functions ..
270 INTRINSIC abs, max, sqrt
271* ..
272* .. Executable Statements ..
273*
274* Decode the input arguments
275*
276 IF( lsame( jobvl, 'N' ) ) THEN
277 ijobvl = 1
278 ilvl = .false.
279 ELSE IF( lsame( jobvl, 'V' ) ) THEN
280 ijobvl = 2
281 ilvl = .true.
282 ELSE
283 ijobvl = -1
284 ilvl = .false.
285 END IF
286*
287 IF( lsame( jobvr, 'N' ) ) THEN
288 ijobvr = 1
289 ilvr = .false.
290 ELSE IF( lsame( jobvr, 'V' ) ) THEN
291 ijobvr = 2
292 ilvr = .true.
293 ELSE
294 ijobvr = -1
295 ilvr = .false.
296 END IF
297 ilv = ilvl .OR. ilvr
298*
299* Test the input arguments
300*
301 info = 0
302 lquery = ( lwork.EQ.-1 )
303 IF( ijobvl.LE.0 ) THEN
304 info = -1
305 ELSE IF( ijobvr.LE.0 ) THEN
306 info = -2
307 ELSE IF( n.LT.0 ) THEN
308 info = -3
309 ELSE IF( lda.LT.max( 1, n ) ) THEN
310 info = -5
311 ELSE IF( ldb.LT.max( 1, n ) ) THEN
312 info = -7
313 ELSE IF( ldvl.LT.1 .OR. ( ilvl .AND. ldvl.LT.n ) ) THEN
314 info = -12
315 ELSE IF( ldvr.LT.1 .OR. ( ilvr .AND. ldvr.LT.n ) ) THEN
316 info = -14
317 END IF
318*
319* Compute workspace
320* (Note: Comments in the code beginning "Workspace:" describe the
321* minimal amount of workspace needed at that point in the code,
322* as well as the preferred amount for good performance.
323* NB refers to the optimal block size for the immediately
324* following subroutine, as returned by ILAENV. The workspace is
325* computed assuming ILO = 1 and IHI = N, the worst case.)
326*
327 IF( info.EQ.0 ) THEN
328 minwrk = max( 1, 8*n )
329 maxwrk = max( 1, n*( 7 +
330 $ ilaenv( 1, 'DGEQRF', ' ', n, 1, n, 0 ) ) )
331 maxwrk = max( maxwrk, n*( 7 +
332 $ ilaenv( 1, 'DORMQR', ' ', n, 1, n, 0 ) ) )
333 IF( ilvl ) THEN
334 maxwrk = max( maxwrk, n*( 7 +
335 $ ilaenv( 1, 'DORGQR', ' ', n, 1, n, -1 ) ) )
336 END IF
337 work( 1 ) = maxwrk
338*
339 IF( lwork.LT.minwrk .AND. .NOT.lquery )
340 $ info = -16
341 END IF
342*
343 IF( info.NE.0 ) THEN
344 CALL xerbla( 'DGGEV ', -info )
345 RETURN
346 ELSE IF( lquery ) THEN
347 RETURN
348 END IF
349*
350* Quick return if possible
351*
352 IF( n.EQ.0 )
353 $ RETURN
354*
355* Get machine constants
356*
357 eps = dlamch( 'P' )
358 smlnum = dlamch( 'S' )
359 bignum = one / smlnum
360 smlnum = sqrt( smlnum ) / eps
361 bignum = one / smlnum
362*
363* Scale A if max element outside range [SMLNUM,BIGNUM]
364*
365 anrm = dlange( 'M', n, n, a, lda, work )
366 ilascl = .false.
367 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
368 anrmto = smlnum
369 ilascl = .true.
370 ELSE IF( anrm.GT.bignum ) THEN
371 anrmto = bignum
372 ilascl = .true.
373 END IF
374 IF( ilascl )
375 $ CALL dlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
376*
377* Scale B if max element outside range [SMLNUM,BIGNUM]
378*
379 bnrm = dlange( 'M', n, n, b, ldb, work )
380 ilbscl = .false.
381 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
382 bnrmto = smlnum
383 ilbscl = .true.
384 ELSE IF( bnrm.GT.bignum ) THEN
385 bnrmto = bignum
386 ilbscl = .true.
387 END IF
388 IF( ilbscl )
389 $ CALL dlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
390*
391* Permute the matrices A, B to isolate eigenvalues if possible
392* (Workspace: need 6*N)
393*
394 ileft = 1
395 iright = n + 1
396 iwrk = iright + n
397 CALL dggbal( 'P', n, a, lda, b, ldb, ilo, ihi, work( ileft ),
398 $ work( iright ), work( iwrk ), ierr )
399*
400* Reduce B to triangular form (QR decomposition of B)
401* (Workspace: need N, prefer N*NB)
402*
403 irows = ihi + 1 - ilo
404 IF( ilv ) THEN
405 icols = n + 1 - ilo
406 ELSE
407 icols = irows
408 END IF
409 itau = iwrk
410 iwrk = itau + irows
411 CALL dgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
412 $ work( iwrk ), lwork+1-iwrk, ierr )
413*
414* Apply the orthogonal transformation to matrix A
415* (Workspace: need N, prefer N*NB)
416*
417 CALL dormqr( 'L', 'T', irows, icols, irows, b( ilo, ilo ), ldb,
418 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
419 $ lwork+1-iwrk, ierr )
420*
421* Initialize VL
422* (Workspace: need N, prefer N*NB)
423*
424 IF( ilvl ) THEN
425 CALL dlaset( 'Full', n, n, zero, one, vl, ldvl )
426 IF( irows.GT.1 ) THEN
427 CALL dlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
428 $ vl( ilo+1, ilo ), ldvl )
429 END IF
430 CALL dorgqr( irows, irows, irows, vl( ilo, ilo ), ldvl,
431 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
432 END IF
433*
434* Initialize VR
435*
436 IF( ilvr )
437 $ CALL dlaset( 'Full', n, n, zero, one, vr, ldvr )
438*
439* Reduce to generalized Hessenberg form
440* (Workspace: none needed)
441*
442 IF( ilv ) THEN
443*
444* Eigenvectors requested -- work on whole matrix.
445*
446 CALL dgghrd( jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb, vl,
447 $ ldvl, vr, ldvr, ierr )
448 ELSE
449 CALL dgghrd( 'N', 'N', irows, 1, irows, a( ilo, ilo ), lda,
450 $ b( ilo, ilo ), ldb, vl, ldvl, vr, ldvr, ierr )
451 END IF
452*
453* Perform QZ algorithm (Compute eigenvalues, and optionally, the
454* Schur forms and Schur vectors)
455* (Workspace: need N)
456*
457 iwrk = itau
458 IF( ilv ) THEN
459 chtemp = 'S'
460 ELSE
461 chtemp = 'E'
462 END IF
463 CALL dhgeqz( chtemp, jobvl, jobvr, n, ilo, ihi, a, lda, b, ldb,
464 $ alphar, alphai, beta, vl, ldvl, vr, ldvr,
465 $ work( iwrk ), lwork+1-iwrk, ierr )
466 IF( ierr.NE.0 ) THEN
467 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
468 info = ierr
469 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
470 info = ierr - n
471 ELSE
472 info = n + 1
473 END IF
474 GO TO 110
475 END IF
476*
477* Compute Eigenvectors
478* (Workspace: need 6*N)
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 dtgevc( chtemp, 'B', ldumma, n, a, lda, b, ldb, vl, ldvl,
491 $ vr, ldvr, n, in, work( iwrk ), ierr )
492 IF( ierr.NE.0 ) THEN
493 info = n + 2
494 GO TO 110
495 END IF
496*
497* Undo balancing on VL and VR and normalization
498* (Workspace: none needed)
499*
500 IF( ilvl ) THEN
501 CALL dggbak( 'P', 'L', n, ilo, ihi, work( ileft ),
502 $ work( iright ), n, vl, ldvl, ierr )
503 DO 50 jc = 1, n
504 IF( alphai( jc ).LT.zero )
505 $ GO TO 50
506 temp = zero
507 IF( alphai( jc ).EQ.zero ) THEN
508 DO 10 jr = 1, n
509 temp = max( temp, abs( vl( jr, jc ) ) )
510 10 CONTINUE
511 ELSE
512 DO 20 jr = 1, n
513 temp = max( temp, abs( vl( jr, jc ) )+
514 $ abs( vl( jr, jc+1 ) ) )
515 20 CONTINUE
516 END IF
517 IF( temp.LT.smlnum )
518 $ GO TO 50
519 temp = one / temp
520 IF( alphai( jc ).EQ.zero ) THEN
521 DO 30 jr = 1, n
522 vl( jr, jc ) = vl( jr, jc )*temp
523 30 CONTINUE
524 ELSE
525 DO 40 jr = 1, n
526 vl( jr, jc ) = vl( jr, jc )*temp
527 vl( jr, jc+1 ) = vl( jr, jc+1 )*temp
528 40 CONTINUE
529 END IF
530 50 CONTINUE
531 END IF
532 IF( ilvr ) THEN
533 CALL dggbak( 'P', 'R', n, ilo, ihi, work( ileft ),
534 $ work( iright ), n, vr, ldvr, ierr )
535 DO 100 jc = 1, n
536 IF( alphai( jc ).LT.zero )
537 $ GO TO 100
538 temp = zero
539 IF( alphai( jc ).EQ.zero ) THEN
540 DO 60 jr = 1, n
541 temp = max( temp, abs( vr( jr, jc ) ) )
542 60 CONTINUE
543 ELSE
544 DO 70 jr = 1, n
545 temp = max( temp, abs( vr( jr, jc ) )+
546 $ abs( vr( jr, jc+1 ) ) )
547 70 CONTINUE
548 END IF
549 IF( temp.LT.smlnum )
550 $ GO TO 100
551 temp = one / temp
552 IF( alphai( jc ).EQ.zero ) THEN
553 DO 80 jr = 1, n
554 vr( jr, jc ) = vr( jr, jc )*temp
555 80 CONTINUE
556 ELSE
557 DO 90 jr = 1, n
558 vr( jr, jc ) = vr( jr, jc )*temp
559 vr( jr, jc+1 ) = vr( jr, jc+1 )*temp
560 90 CONTINUE
561 END IF
562 100 CONTINUE
563 END IF
564*
565* End of eigenvector calculation
566*
567 END IF
568*
569* Undo scaling if necessary
570*
571 110 CONTINUE
572*
573 IF( ilascl ) THEN
574 CALL dlascl( 'G', 0, 0, anrmto, anrm, n, 1, alphar, n, ierr )
575 CALL dlascl( 'G', 0, 0, anrmto, anrm, n, 1, alphai, n, ierr )
576 END IF
577*
578 IF( ilbscl ) THEN
579 CALL dlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
580 END IF
581*
582 work( 1 ) = maxwrk
583 RETURN
584*
585* End of DGGEV
586*
587 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgeqrf(m, n, a, lda, tau, work, lwork, info)
DGEQRF
Definition dgeqrf.f:146
subroutine dggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
DGGBAK
Definition dggbak.f:147
subroutine dggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
DGGBAL
Definition dggbal.f:177
subroutine dggev(jobvl, jobvr, n, a, lda, b, ldb, alphar, alphai, beta, vl, ldvl, vr, ldvr, work, lwork, info)
DGGEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
Definition dggev.f:226
subroutine dgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
DGGHRD
Definition dgghrd.f:207
subroutine dhgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alphar, alphai, beta, q, ldq, z, ldz, work, lwork, info)
DHGEQZ
Definition dhgeqz.f:304
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:143
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:110
subroutine dtgevc(side, howmny, select, n, s, lds, p, ldp, vl, ldvl, vr, ldvr, mm, m, work, info)
DTGEVC
Definition dtgevc.f:295
subroutine dorgqr(m, n, k, a, lda, tau, work, lwork, info)
DORGQR
Definition dorgqr.f:128
subroutine dormqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
DORMQR
Definition dormqr.f:167