LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
slagv2.f
Go to the documentation of this file.
1*> \brief \b SLAGV2 computes the Generalized Schur factorization of a real 2-by-2 matrix pencil (A,B) where B is upper triangular.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SLAGV2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slagv2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slagv2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slagv2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SLAGV2( A, LDA, B, LDB, ALPHAR, ALPHAI, BETA, CSL, SNL,
20* CSR, SNR )
21*
22* .. Scalar Arguments ..
23* INTEGER LDA, LDB
24* REAL CSL, CSR, SNL, SNR
25* ..
26* .. Array Arguments ..
27* REAL A( LDA, * ), ALPHAI( 2 ), ALPHAR( 2 ),
28* $ B( LDB, * ), BETA( 2 )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> SLAGV2 computes the Generalized Schur factorization of a real 2-by-2
38*> matrix pencil (A,B) where B is upper triangular. This routine
39*> computes orthogonal (rotation) matrices given by CSL, SNL and CSR,
40*> SNR such that
41*>
42*> 1) if the pencil (A,B) has two real eigenvalues (include 0/0 or 1/0
43*> types), then
44*>
45*> [ a11 a12 ] := [ CSL SNL ] [ a11 a12 ] [ CSR -SNR ]
46*> [ 0 a22 ] [ -SNL CSL ] [ a21 a22 ] [ SNR CSR ]
47*>
48*> [ b11 b12 ] := [ CSL SNL ] [ b11 b12 ] [ CSR -SNR ]
49*> [ 0 b22 ] [ -SNL CSL ] [ 0 b22 ] [ SNR CSR ],
50*>
51*> 2) if the pencil (A,B) has a pair of complex conjugate eigenvalues,
52*> then
53*>
54*> [ a11 a12 ] := [ CSL SNL ] [ a11 a12 ] [ CSR -SNR ]
55*> [ a21 a22 ] [ -SNL CSL ] [ a21 a22 ] [ SNR CSR ]
56*>
57*> [ b11 0 ] := [ CSL SNL ] [ b11 b12 ] [ CSR -SNR ]
58*> [ 0 b22 ] [ -SNL CSL ] [ 0 b22 ] [ SNR CSR ]
59*>
60*> where b11 >= b22 > 0.
61*>
62*> \endverbatim
63*
64* Arguments:
65* ==========
66*
67*> \param[in,out] A
68*> \verbatim
69*> A is REAL array, dimension (LDA, 2)
70*> On entry, the 2 x 2 matrix A.
71*> On exit, A is overwritten by the ``A-part'' of the
72*> generalized Schur form.
73*> \endverbatim
74*>
75*> \param[in] LDA
76*> \verbatim
77*> LDA is INTEGER
78*> THe leading dimension of the array A. LDA >= 2.
79*> \endverbatim
80*>
81*> \param[in,out] B
82*> \verbatim
83*> B is REAL array, dimension (LDB, 2)
84*> On entry, the upper triangular 2 x 2 matrix B.
85*> On exit, B is overwritten by the ``B-part'' of the
86*> generalized Schur form.
87*> \endverbatim
88*>
89*> \param[in] LDB
90*> \verbatim
91*> LDB is INTEGER
92*> THe leading dimension of the array B. LDB >= 2.
93*> \endverbatim
94*>
95*> \param[out] ALPHAR
96*> \verbatim
97*> ALPHAR is REAL array, dimension (2)
98*> \endverbatim
99*>
100*> \param[out] ALPHAI
101*> \verbatim
102*> ALPHAI is REAL array, dimension (2)
103*> \endverbatim
104*>
105*> \param[out] BETA
106*> \verbatim
107*> BETA is REAL array, dimension (2)
108*> (ALPHAR(k)+i*ALPHAI(k))/BETA(k) are the eigenvalues of the
109*> pencil (A,B), k=1,2, i = sqrt(-1). Note that BETA(k) may
110*> be zero.
111*> \endverbatim
112*>
113*> \param[out] CSL
114*> \verbatim
115*> CSL is REAL
116*> The cosine of the left rotation matrix.
117*> \endverbatim
118*>
119*> \param[out] SNL
120*> \verbatim
121*> SNL is REAL
122*> The sine of the left rotation matrix.
123*> \endverbatim
124*>
125*> \param[out] CSR
126*> \verbatim
127*> CSR is REAL
128*> The cosine of the right rotation matrix.
129*> \endverbatim
130*>
131*> \param[out] SNR
132*> \verbatim
133*> SNR is REAL
134*> The sine of the right rotation matrix.
135*> \endverbatim
136*
137* Authors:
138* ========
139*
140*> \author Univ. of Tennessee
141*> \author Univ. of California Berkeley
142*> \author Univ. of Colorado Denver
143*> \author NAG Ltd.
144*
145*> \ingroup lagv2
146*
147*> \par Contributors:
148* ==================
149*>
150*> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
151*
152* =====================================================================
153 SUBROUTINE slagv2( A, LDA, B, LDB, ALPHAR, ALPHAI, BETA, CSL,
154 $ SNL,
155 $ CSR, SNR )
156*
157* -- LAPACK auxiliary routine --
158* -- LAPACK is a software package provided by Univ. of Tennessee, --
159* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
160*
161* .. Scalar Arguments ..
162 INTEGER LDA, LDB
163 REAL CSL, CSR, SNL, SNR
164* ..
165* .. Array Arguments ..
166 REAL A( LDA, * ), ALPHAI( 2 ), ALPHAR( 2 ),
167 $ B( LDB, * ), BETA( 2 )
168* ..
169*
170* =====================================================================
171*
172* .. Parameters ..
173 REAL ZERO, ONE
174 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
175* ..
176* .. Local Scalars ..
177 REAL ANORM, ASCALE, BNORM, BSCALE, H1, H2, H3, QQ,
178 $ R, RR, SAFMIN, SCALE1, SCALE2, T, ULP, WI, WR1,
179 $ wr2
180* ..
181* .. External Subroutines ..
182 EXTERNAL slag2, slartg, slasv2, srot
183* ..
184* .. External Functions ..
185 REAL SLAMCH, SLAPY2
186 EXTERNAL SLAMCH, SLAPY2
187* ..
188* .. Intrinsic Functions ..
189 INTRINSIC abs, max
190* ..
191* .. Executable Statements ..
192*
193 safmin = slamch( 'S' )
194 ulp = slamch( 'P' )
195*
196* Scale A
197*
198 anorm = max( abs( a( 1, 1 ) )+abs( a( 2, 1 ) ),
199 $ abs( a( 1, 2 ) )+abs( a( 2, 2 ) ), safmin )
200 ascale = one / anorm
201 a( 1, 1 ) = ascale*a( 1, 1 )
202 a( 1, 2 ) = ascale*a( 1, 2 )
203 a( 2, 1 ) = ascale*a( 2, 1 )
204 a( 2, 2 ) = ascale*a( 2, 2 )
205*
206* Scale B
207*
208 bnorm = max( abs( b( 1, 1 ) ), abs( b( 1, 2 ) )+abs( b( 2, 2 ) ),
209 $ safmin )
210 bscale = one / bnorm
211 b( 1, 1 ) = bscale*b( 1, 1 )
212 b( 1, 2 ) = bscale*b( 1, 2 )
213 b( 2, 2 ) = bscale*b( 2, 2 )
214*
215* Check if A can be deflated
216*
217 IF( abs( a( 2, 1 ) ).LE.ulp ) THEN
218 csl = one
219 snl = zero
220 csr = one
221 snr = zero
222 a( 2, 1 ) = zero
223 b( 2, 1 ) = zero
224 wi = zero
225*
226* Check if B is singular
227*
228 ELSE IF( abs( b( 1, 1 ) ).LE.ulp ) THEN
229 CALL slartg( a( 1, 1 ), a( 2, 1 ), csl, snl, r )
230 csr = one
231 snr = zero
232 CALL srot( 2, a( 1, 1 ), lda, a( 2, 1 ), lda, csl, snl )
233 CALL srot( 2, b( 1, 1 ), ldb, b( 2, 1 ), ldb, csl, snl )
234 a( 2, 1 ) = zero
235 b( 1, 1 ) = zero
236 b( 2, 1 ) = zero
237 wi = zero
238*
239 ELSE IF( abs( b( 2, 2 ) ).LE.ulp ) THEN
240 CALL slartg( a( 2, 2 ), a( 2, 1 ), csr, snr, t )
241 snr = -snr
242 CALL srot( 2, a( 1, 1 ), 1, a( 1, 2 ), 1, csr, snr )
243 CALL srot( 2, b( 1, 1 ), 1, b( 1, 2 ), 1, csr, snr )
244 csl = one
245 snl = zero
246 a( 2, 1 ) = zero
247 b( 2, 1 ) = zero
248 b( 2, 2 ) = zero
249 wi = zero
250*
251 ELSE
252*
253* B is nonsingular, first compute the eigenvalues of (A,B)
254*
255 CALL slag2( a, lda, b, ldb, safmin, scale1, scale2, wr1,
256 $ wr2,
257 $ wi )
258*
259 IF( wi.EQ.zero ) THEN
260*
261* two real eigenvalues, compute s*A-w*B
262*
263 h1 = scale1*a( 1, 1 ) - wr1*b( 1, 1 )
264 h2 = scale1*a( 1, 2 ) - wr1*b( 1, 2 )
265 h3 = scale1*a( 2, 2 ) - wr1*b( 2, 2 )
266*
267 rr = slapy2( h1, h2 )
268 qq = slapy2( scale1*a( 2, 1 ), h3 )
269*
270 IF( rr.GT.qq ) THEN
271*
272* find right rotation matrix to zero 1,1 element of
273* (sA - wB)
274*
275 CALL slartg( h2, h1, csr, snr, t )
276*
277 ELSE
278*
279* find right rotation matrix to zero 2,1 element of
280* (sA - wB)
281*
282 CALL slartg( h3, scale1*a( 2, 1 ), csr, snr, t )
283*
284 END IF
285*
286 snr = -snr
287 CALL srot( 2, a( 1, 1 ), 1, a( 1, 2 ), 1, csr, snr )
288 CALL srot( 2, b( 1, 1 ), 1, b( 1, 2 ), 1, csr, snr )
289*
290* compute inf norms of A and B
291*
292 h1 = max( abs( a( 1, 1 ) )+abs( a( 1, 2 ) ),
293 $ abs( a( 2, 1 ) )+abs( a( 2, 2 ) ) )
294 h2 = max( abs( b( 1, 1 ) )+abs( b( 1, 2 ) ),
295 $ abs( b( 2, 1 ) )+abs( b( 2, 2 ) ) )
296*
297 IF( ( scale1*h1 ).GE.abs( wr1 )*h2 ) THEN
298*
299* find left rotation matrix Q to zero out B(2,1)
300*
301 CALL slartg( b( 1, 1 ), b( 2, 1 ), csl, snl, r )
302*
303 ELSE
304*
305* find left rotation matrix Q to zero out A(2,1)
306*
307 CALL slartg( a( 1, 1 ), a( 2, 1 ), csl, snl, r )
308*
309 END IF
310*
311 CALL srot( 2, a( 1, 1 ), lda, a( 2, 1 ), lda, csl, snl )
312 CALL srot( 2, b( 1, 1 ), ldb, b( 2, 1 ), ldb, csl, snl )
313*
314 a( 2, 1 ) = zero
315 b( 2, 1 ) = zero
316*
317 ELSE
318*
319* a pair of complex conjugate eigenvalues
320* first compute the SVD of the matrix B
321*
322 CALL slasv2( b( 1, 1 ), b( 1, 2 ), b( 2, 2 ), r, t, snr,
323 $ csr, snl, csl )
324*
325* Form (A,B) := Q(A,B)Z**T where Q is left rotation matrix and
326* Z is right rotation matrix computed from SLASV2
327*
328 CALL srot( 2, a( 1, 1 ), lda, a( 2, 1 ), lda, csl, snl )
329 CALL srot( 2, b( 1, 1 ), ldb, b( 2, 1 ), ldb, csl, snl )
330 CALL srot( 2, a( 1, 1 ), 1, a( 1, 2 ), 1, csr, snr )
331 CALL srot( 2, b( 1, 1 ), 1, b( 1, 2 ), 1, csr, snr )
332*
333 b( 2, 1 ) = zero
334 b( 1, 2 ) = zero
335*
336 END IF
337*
338 END IF
339*
340* Unscaling
341*
342 a( 1, 1 ) = anorm*a( 1, 1 )
343 a( 2, 1 ) = anorm*a( 2, 1 )
344 a( 1, 2 ) = anorm*a( 1, 2 )
345 a( 2, 2 ) = anorm*a( 2, 2 )
346 b( 1, 1 ) = bnorm*b( 1, 1 )
347 b( 2, 1 ) = bnorm*b( 2, 1 )
348 b( 1, 2 ) = bnorm*b( 1, 2 )
349 b( 2, 2 ) = bnorm*b( 2, 2 )
350*
351 IF( wi.EQ.zero ) THEN
352 alphar( 1 ) = a( 1, 1 )
353 alphar( 2 ) = a( 2, 2 )
354 alphai( 1 ) = zero
355 alphai( 2 ) = zero
356 beta( 1 ) = b( 1, 1 )
357 beta( 2 ) = b( 2, 2 )
358 ELSE
359 alphar( 1 ) = anorm*wr1 / scale1 / bnorm
360 alphai( 1 ) = anorm*wi / scale1 / bnorm
361 alphar( 2 ) = alphar( 1 )
362 alphai( 2 ) = -alphai( 1 )
363 beta( 1 ) = one
364 beta( 2 ) = one
365 END IF
366*
367 RETURN
368*
369* End of SLAGV2
370*
371 END
subroutine slag2(a, lda, b, ldb, safmin, scale1, scale2, wr1, wr2, wi)
SLAG2 computes the eigenvalues of a 2-by-2 generalized eigenvalue problem, with scaling as necessary ...
Definition slag2.f:154
subroutine slagv2(a, lda, b, ldb, alphar, alphai, beta, csl, snl, csr, snr)
SLAGV2 computes the Generalized Schur factorization of a real 2-by-2 matrix pencil (A,...
Definition slagv2.f:156
subroutine slartg(f, g, c, s, r)
SLARTG generates a plane rotation with real cosine and real sine.
Definition slartg.f90:111
subroutine slasv2(f, g, h, ssmin, ssmax, snr, csr, snl, csl)
SLASV2 computes the singular value decomposition of a 2-by-2 triangular matrix.
Definition slasv2.f:134
subroutine srot(n, sx, incx, sy, incy, c, s)
SROT
Definition srot.f:92