You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

ctrevc.f 15 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483
  1. *> \brief \b CTREVC
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download CTREVC + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ctrevc.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ctrevc.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ctrevc.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE CTREVC( SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR,
  22. * LDVR, MM, M, WORK, RWORK, INFO )
  23. *
  24. * .. Scalar Arguments ..
  25. * CHARACTER HOWMNY, SIDE
  26. * INTEGER INFO, LDT, LDVL, LDVR, M, MM, N
  27. * ..
  28. * .. Array Arguments ..
  29. * LOGICAL SELECT( * )
  30. * REAL RWORK( * )
  31. * COMPLEX T( LDT, * ), VL( LDVL, * ), VR( LDVR, * ),
  32. * $ WORK( * )
  33. * ..
  34. *
  35. *
  36. *> \par Purpose:
  37. * =============
  38. *>
  39. *> \verbatim
  40. *>
  41. *> CTREVC computes some or all of the right and/or left eigenvectors of
  42. *> a complex upper triangular matrix T.
  43. *> Matrices of this type are produced by the Schur factorization of
  44. *> a complex general matrix: A = Q*T*Q**H, as computed by CHSEQR.
  45. *>
  46. *> The right eigenvector x and the left eigenvector y of T corresponding
  47. *> to an eigenvalue w are defined by:
  48. *>
  49. *> T*x = w*x, (y**H)*T = w*(y**H)
  50. *>
  51. *> where y**H denotes the conjugate transpose of the vector y.
  52. *> The eigenvalues are not input to this routine, but are read directly
  53. *> from the diagonal of T.
  54. *>
  55. *> This routine returns the matrices X and/or Y of right and left
  56. *> eigenvectors of T, or the products Q*X and/or Q*Y, where Q is an
  57. *> input matrix. If Q is the unitary factor that reduces a matrix A to
  58. *> Schur form T, then Q*X and Q*Y are the matrices of right and left
  59. *> eigenvectors of A.
  60. *> \endverbatim
  61. *
  62. * Arguments:
  63. * ==========
  64. *
  65. *> \param[in] SIDE
  66. *> \verbatim
  67. *> SIDE is CHARACTER*1
  68. *> = 'R': compute right eigenvectors only;
  69. *> = 'L': compute left eigenvectors only;
  70. *> = 'B': compute both right and left eigenvectors.
  71. *> \endverbatim
  72. *>
  73. *> \param[in] HOWMNY
  74. *> \verbatim
  75. *> HOWMNY is CHARACTER*1
  76. *> = 'A': compute all right and/or left eigenvectors;
  77. *> = 'B': compute all right and/or left eigenvectors,
  78. *> backtransformed using the matrices supplied in
  79. *> VR and/or VL;
  80. *> = 'S': compute selected right and/or left eigenvectors,
  81. *> as indicated by the logical array SELECT.
  82. *> \endverbatim
  83. *>
  84. *> \param[in] SELECT
  85. *> \verbatim
  86. *> SELECT is LOGICAL array, dimension (N)
  87. *> If HOWMNY = 'S', SELECT specifies the eigenvectors to be
  88. *> computed.
  89. *> The eigenvector corresponding to the j-th eigenvalue is
  90. *> computed if SELECT(j) = .TRUE..
  91. *> Not referenced if HOWMNY = 'A' or 'B'.
  92. *> \endverbatim
  93. *>
  94. *> \param[in] N
  95. *> \verbatim
  96. *> N is INTEGER
  97. *> The order of the matrix T. N >= 0.
  98. *> \endverbatim
  99. *>
  100. *> \param[in,out] T
  101. *> \verbatim
  102. *> T is COMPLEX array, dimension (LDT,N)
  103. *> The upper triangular matrix T. T is modified, but restored
  104. *> on exit.
  105. *> \endverbatim
  106. *>
  107. *> \param[in] LDT
  108. *> \verbatim
  109. *> LDT is INTEGER
  110. *> The leading dimension of the array T. LDT >= max(1,N).
  111. *> \endverbatim
  112. *>
  113. *> \param[in,out] VL
  114. *> \verbatim
  115. *> VL is COMPLEX array, dimension (LDVL,MM)
  116. *> On entry, if SIDE = 'L' or 'B' and HOWMNY = 'B', VL must
  117. *> contain an N-by-N matrix Q (usually the unitary matrix Q of
  118. *> Schur vectors returned by CHSEQR).
  119. *> On exit, if SIDE = 'L' or 'B', VL contains:
  120. *> if HOWMNY = 'A', the matrix Y of left eigenvectors of T;
  121. *> if HOWMNY = 'B', the matrix Q*Y;
  122. *> if HOWMNY = 'S', the left eigenvectors of T specified by
  123. *> SELECT, stored consecutively in the columns
  124. *> of VL, in the same order as their
  125. *> eigenvalues.
  126. *> Not referenced if SIDE = 'R'.
  127. *> \endverbatim
  128. *>
  129. *> \param[in] LDVL
  130. *> \verbatim
  131. *> LDVL is INTEGER
  132. *> The leading dimension of the array VL. LDVL >= 1, and if
  133. *> SIDE = 'L' or 'B', LDVL >= N.
  134. *> \endverbatim
  135. *>
  136. *> \param[in,out] VR
  137. *> \verbatim
  138. *> VR is COMPLEX array, dimension (LDVR,MM)
  139. *> On entry, if SIDE = 'R' or 'B' and HOWMNY = 'B', VR must
  140. *> contain an N-by-N matrix Q (usually the unitary matrix Q of
  141. *> Schur vectors returned by CHSEQR).
  142. *> On exit, if SIDE = 'R' or 'B', VR contains:
  143. *> if HOWMNY = 'A', the matrix X of right eigenvectors of T;
  144. *> if HOWMNY = 'B', the matrix Q*X;
  145. *> if HOWMNY = 'S', the right eigenvectors of T specified by
  146. *> SELECT, stored consecutively in the columns
  147. *> of VR, in the same order as their
  148. *> eigenvalues.
  149. *> Not referenced if SIDE = 'L'.
  150. *> \endverbatim
  151. *>
  152. *> \param[in] LDVR
  153. *> \verbatim
  154. *> LDVR is INTEGER
  155. *> The leading dimension of the array VR. LDVR >= 1, and if
  156. *> SIDE = 'R' or 'B'; LDVR >= N.
  157. *> \endverbatim
  158. *>
  159. *> \param[in] MM
  160. *> \verbatim
  161. *> MM is INTEGER
  162. *> The number of columns in the arrays VL and/or VR. MM >= M.
  163. *> \endverbatim
  164. *>
  165. *> \param[out] M
  166. *> \verbatim
  167. *> M is INTEGER
  168. *> The number of columns in the arrays VL and/or VR actually
  169. *> used to store the eigenvectors. If HOWMNY = 'A' or 'B', M
  170. *> is set to N. Each selected eigenvector occupies one
  171. *> column.
  172. *> \endverbatim
  173. *>
  174. *> \param[out] WORK
  175. *> \verbatim
  176. *> WORK is COMPLEX array, dimension (2*N)
  177. *> \endverbatim
  178. *>
  179. *> \param[out] RWORK
  180. *> \verbatim
  181. *> RWORK is REAL array, dimension (N)
  182. *> \endverbatim
  183. *>
  184. *> \param[out] INFO
  185. *> \verbatim
  186. *> INFO is INTEGER
  187. *> = 0: successful exit
  188. *> < 0: if INFO = -i, the i-th argument had an illegal value
  189. *> \endverbatim
  190. *
  191. * Authors:
  192. * ========
  193. *
  194. *> \author Univ. of Tennessee
  195. *> \author Univ. of California Berkeley
  196. *> \author Univ. of Colorado Denver
  197. *> \author NAG Ltd.
  198. *
  199. *> \ingroup complexOTHERcomputational
  200. *
  201. *> \par Further Details:
  202. * =====================
  203. *>
  204. *> \verbatim
  205. *>
  206. *> The algorithm used in this program is basically backward (forward)
  207. *> substitution, with scaling to make the the code robust against
  208. *> possible overflow.
  209. *>
  210. *> Each eigenvector is normalized so that the element of largest
  211. *> magnitude has magnitude 1; here the magnitude of a complex number
  212. *> (x,y) is taken to be |x| + |y|.
  213. *> \endverbatim
  214. *>
  215. * =====================================================================
  216. SUBROUTINE CTREVC( SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL, VR,
  217. $ LDVR, MM, M, WORK, RWORK, INFO )
  218. *
  219. * -- LAPACK computational routine --
  220. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  221. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  222. *
  223. * .. Scalar Arguments ..
  224. CHARACTER HOWMNY, SIDE
  225. INTEGER INFO, LDT, LDVL, LDVR, M, MM, N
  226. * ..
  227. * .. Array Arguments ..
  228. LOGICAL SELECT( * )
  229. REAL RWORK( * )
  230. COMPLEX T( LDT, * ), VL( LDVL, * ), VR( LDVR, * ),
  231. $ WORK( * )
  232. * ..
  233. *
  234. * =====================================================================
  235. *
  236. * .. Parameters ..
  237. REAL ZERO, ONE
  238. PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 )
  239. COMPLEX CMZERO, CMONE
  240. PARAMETER ( CMZERO = ( 0.0E+0, 0.0E+0 ),
  241. $ CMONE = ( 1.0E+0, 0.0E+0 ) )
  242. * ..
  243. * .. Local Scalars ..
  244. LOGICAL ALLV, BOTHV, LEFTV, OVER, RIGHTV, SOMEV
  245. INTEGER I, II, IS, J, K, KI
  246. REAL OVFL, REMAX, SCALE, SMIN, SMLNUM, ULP, UNFL
  247. COMPLEX CDUM
  248. * ..
  249. * .. External Functions ..
  250. LOGICAL LSAME
  251. INTEGER ICAMAX
  252. REAL SCASUM, SLAMCH
  253. EXTERNAL LSAME, ICAMAX, SCASUM, SLAMCH
  254. * ..
  255. * .. External Subroutines ..
  256. EXTERNAL CCOPY, CGEMV, CLATRS, CSSCAL, SLABAD, XERBLA
  257. * ..
  258. * .. Intrinsic Functions ..
  259. INTRINSIC ABS, AIMAG, CMPLX, CONJG, MAX, REAL
  260. * ..
  261. * .. Statement Functions ..
  262. REAL CABS1
  263. * ..
  264. * .. Statement Function definitions ..
  265. CABS1( CDUM ) = ABS( REAL( CDUM ) ) + ABS( AIMAG( CDUM ) )
  266. * ..
  267. * .. Executable Statements ..
  268. *
  269. * Decode and test the input parameters
  270. *
  271. BOTHV = LSAME( SIDE, 'B' )
  272. RIGHTV = LSAME( SIDE, 'R' ) .OR. BOTHV
  273. LEFTV = LSAME( SIDE, 'L' ) .OR. BOTHV
  274. *
  275. ALLV = LSAME( HOWMNY, 'A' )
  276. OVER = LSAME( HOWMNY, 'B' )
  277. SOMEV = LSAME( HOWMNY, 'S' )
  278. *
  279. * Set M to the number of columns required to store the selected
  280. * eigenvectors.
  281. *
  282. IF( SOMEV ) THEN
  283. M = 0
  284. DO 10 J = 1, N
  285. IF( SELECT( J ) )
  286. $ M = M + 1
  287. 10 CONTINUE
  288. ELSE
  289. M = N
  290. END IF
  291. *
  292. INFO = 0
  293. IF( .NOT.RIGHTV .AND. .NOT.LEFTV ) THEN
  294. INFO = -1
  295. ELSE IF( .NOT.ALLV .AND. .NOT.OVER .AND. .NOT.SOMEV ) THEN
  296. INFO = -2
  297. ELSE IF( N.LT.0 ) THEN
  298. INFO = -4
  299. ELSE IF( LDT.LT.MAX( 1, N ) ) THEN
  300. INFO = -6
  301. ELSE IF( LDVL.LT.1 .OR. ( LEFTV .AND. LDVL.LT.N ) ) THEN
  302. INFO = -8
  303. ELSE IF( LDVR.LT.1 .OR. ( RIGHTV .AND. LDVR.LT.N ) ) THEN
  304. INFO = -10
  305. ELSE IF( MM.LT.M ) THEN
  306. INFO = -11
  307. END IF
  308. IF( INFO.NE.0 ) THEN
  309. CALL XERBLA( 'CTREVC', -INFO )
  310. RETURN
  311. END IF
  312. *
  313. * Quick return if possible.
  314. *
  315. IF( N.EQ.0 )
  316. $ RETURN
  317. *
  318. * Set the constants to control overflow.
  319. *
  320. UNFL = SLAMCH( 'Safe minimum' )
  321. OVFL = ONE / UNFL
  322. CALL SLABAD( UNFL, OVFL )
  323. ULP = SLAMCH( 'Precision' )
  324. SMLNUM = UNFL*( N / ULP )
  325. *
  326. * Store the diagonal elements of T in working array WORK.
  327. *
  328. DO 20 I = 1, N
  329. WORK( I+N ) = T( I, I )
  330. 20 CONTINUE
  331. *
  332. * Compute 1-norm of each column of strictly upper triangular
  333. * part of T to control overflow in triangular solver.
  334. *
  335. RWORK( 1 ) = ZERO
  336. DO 30 J = 2, N
  337. RWORK( J ) = SCASUM( J-1, T( 1, J ), 1 )
  338. 30 CONTINUE
  339. *
  340. IF( RIGHTV ) THEN
  341. *
  342. * Compute right eigenvectors.
  343. *
  344. IS = M
  345. DO 80 KI = N, 1, -1
  346. *
  347. IF( SOMEV ) THEN
  348. IF( .NOT.SELECT( KI ) )
  349. $ GO TO 80
  350. END IF
  351. SMIN = MAX( ULP*( CABS1( T( KI, KI ) ) ), SMLNUM )
  352. *
  353. WORK( 1 ) = CMONE
  354. *
  355. * Form right-hand side.
  356. *
  357. DO 40 K = 1, KI - 1
  358. WORK( K ) = -T( K, KI )
  359. 40 CONTINUE
  360. *
  361. * Solve the triangular system:
  362. * (T(1:KI-1,1:KI-1) - T(KI,KI))*X = SCALE*WORK.
  363. *
  364. DO 50 K = 1, KI - 1
  365. T( K, K ) = T( K, K ) - T( KI, KI )
  366. IF( CABS1( T( K, K ) ).LT.SMIN )
  367. $ T( K, K ) = SMIN
  368. 50 CONTINUE
  369. *
  370. IF( KI.GT.1 ) THEN
  371. CALL CLATRS( 'Upper', 'No transpose', 'Non-unit', 'Y',
  372. $ KI-1, T, LDT, WORK( 1 ), SCALE, RWORK,
  373. $ INFO )
  374. WORK( KI ) = SCALE
  375. END IF
  376. *
  377. * Copy the vector x or Q*x to VR and normalize.
  378. *
  379. IF( .NOT.OVER ) THEN
  380. CALL CCOPY( KI, WORK( 1 ), 1, VR( 1, IS ), 1 )
  381. *
  382. II = ICAMAX( KI, VR( 1, IS ), 1 )
  383. REMAX = ONE / CABS1( VR( II, IS ) )
  384. CALL CSSCAL( KI, REMAX, VR( 1, IS ), 1 )
  385. *
  386. DO 60 K = KI + 1, N
  387. VR( K, IS ) = CMZERO
  388. 60 CONTINUE
  389. ELSE
  390. IF( KI.GT.1 )
  391. $ CALL CGEMV( 'N', N, KI-1, CMONE, VR, LDVR, WORK( 1 ),
  392. $ 1, CMPLX( SCALE ), VR( 1, KI ), 1 )
  393. *
  394. II = ICAMAX( N, VR( 1, KI ), 1 )
  395. REMAX = ONE / CABS1( VR( II, KI ) )
  396. CALL CSSCAL( N, REMAX, VR( 1, KI ), 1 )
  397. END IF
  398. *
  399. * Set back the original diagonal elements of T.
  400. *
  401. DO 70 K = 1, KI - 1
  402. T( K, K ) = WORK( K+N )
  403. 70 CONTINUE
  404. *
  405. IS = IS - 1
  406. 80 CONTINUE
  407. END IF
  408. *
  409. IF( LEFTV ) THEN
  410. *
  411. * Compute left eigenvectors.
  412. *
  413. IS = 1
  414. DO 130 KI = 1, N
  415. *
  416. IF( SOMEV ) THEN
  417. IF( .NOT.SELECT( KI ) )
  418. $ GO TO 130
  419. END IF
  420. SMIN = MAX( ULP*( CABS1( T( KI, KI ) ) ), SMLNUM )
  421. *
  422. WORK( N ) = CMONE
  423. *
  424. * Form right-hand side.
  425. *
  426. DO 90 K = KI + 1, N
  427. WORK( K ) = -CONJG( T( KI, K ) )
  428. 90 CONTINUE
  429. *
  430. * Solve the triangular system:
  431. * (T(KI+1:N,KI+1:N) - T(KI,KI))**H*X = SCALE*WORK.
  432. *
  433. DO 100 K = KI + 1, N
  434. T( K, K ) = T( K, K ) - T( KI, KI )
  435. IF( CABS1( T( K, K ) ).LT.SMIN )
  436. $ T( K, K ) = SMIN
  437. 100 CONTINUE
  438. *
  439. IF( KI.LT.N ) THEN
  440. CALL CLATRS( 'Upper', 'Conjugate transpose', 'Non-unit',
  441. $ 'Y', N-KI, T( KI+1, KI+1 ), LDT,
  442. $ WORK( KI+1 ), SCALE, RWORK, INFO )
  443. WORK( KI ) = SCALE
  444. END IF
  445. *
  446. * Copy the vector x or Q*x to VL and normalize.
  447. *
  448. IF( .NOT.OVER ) THEN
  449. CALL CCOPY( N-KI+1, WORK( KI ), 1, VL( KI, IS ), 1 )
  450. *
  451. II = ICAMAX( N-KI+1, VL( KI, IS ), 1 ) + KI - 1
  452. REMAX = ONE / CABS1( VL( II, IS ) )
  453. CALL CSSCAL( N-KI+1, REMAX, VL( KI, IS ), 1 )
  454. *
  455. DO 110 K = 1, KI - 1
  456. VL( K, IS ) = CMZERO
  457. 110 CONTINUE
  458. ELSE
  459. IF( KI.LT.N )
  460. $ CALL CGEMV( 'N', N, N-KI, CMONE, VL( 1, KI+1 ), LDVL,
  461. $ WORK( KI+1 ), 1, CMPLX( SCALE ),
  462. $ VL( 1, KI ), 1 )
  463. *
  464. II = ICAMAX( N, VL( 1, KI ), 1 )
  465. REMAX = ONE / CABS1( VL( II, KI ) )
  466. CALL CSSCAL( N, REMAX, VL( 1, KI ), 1 )
  467. END IF
  468. *
  469. * Set back the original diagonal elements of T.
  470. *
  471. DO 120 K = KI + 1, N
  472. T( K, K ) = WORK( K+N )
  473. 120 CONTINUE
  474. *
  475. IS = IS + 1
  476. 130 CONTINUE
  477. END IF
  478. *
  479. RETURN
  480. *
  481. * End of CTREVC
  482. *
  483. END