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.

cgesvdx.f 28 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834
  1. *> \brief <b> CGESVDX computes the singular value decomposition (SVD) 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 CGESVDX + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgesvdx.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgesvdx.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgesvdx.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE CGESVDX( JOBU, JOBVT, RANGE, M, N, A, LDA, VL, VU,
  22. * $ IL, IU, NS, S, U, LDU, VT, LDVT, WORK,
  23. * $ LWORK, RWORK, IWORK, INFO )
  24. *
  25. *
  26. * .. Scalar Arguments ..
  27. * CHARACTER JOBU, JOBVT, RANGE
  28. * INTEGER IL, INFO, IU, LDA, LDU, LDVT, LWORK, M, N, NS
  29. * REAL VL, VU
  30. * ..
  31. * .. Array Arguments ..
  32. * INTEGER IWORK( * )
  33. * REAL S( * ), RWORK( * )
  34. * COMPLEX A( LDA, * ), U( LDU, * ), VT( LDVT, * ),
  35. * $ WORK( * )
  36. * ..
  37. *
  38. *
  39. *> \par Purpose:
  40. * =============
  41. *>
  42. *> \verbatim
  43. *>
  44. *> CGESVDX computes the singular value decomposition (SVD) of a complex
  45. *> M-by-N matrix A, optionally computing the left and/or right singular
  46. *> vectors. The SVD is written
  47. *>
  48. *> A = U * SIGMA * transpose(V)
  49. *>
  50. *> where SIGMA is an M-by-N matrix which is zero except for its
  51. *> min(m,n) diagonal elements, U is an M-by-M unitary matrix, and
  52. *> V is an N-by-N unitary matrix. The diagonal elements of SIGMA
  53. *> are the singular values of A; they are real and non-negative, and
  54. *> are returned in descending order. The first min(m,n) columns of
  55. *> U and V are the left and right singular vectors of A.
  56. *>
  57. *> CGESVDX uses an eigenvalue problem for obtaining the SVD, which
  58. *> allows for the computation of a subset of singular values and
  59. *> vectors. See SBDSVDX for details.
  60. *>
  61. *> Note that the routine returns V**T, not V.
  62. *> \endverbatim
  63. *
  64. * Arguments:
  65. * ==========
  66. *
  67. *> \param[in] JOBU
  68. *> \verbatim
  69. *> JOBU is CHARACTER*1
  70. *> Specifies options for computing all or part of the matrix U:
  71. *> = 'V': the first min(m,n) columns of U (the left singular
  72. *> vectors) or as specified by RANGE are returned in
  73. *> the array U;
  74. *> = 'N': no columns of U (no left singular vectors) are
  75. *> computed.
  76. *> \endverbatim
  77. *>
  78. *> \param[in] JOBVT
  79. *> \verbatim
  80. *> JOBVT is CHARACTER*1
  81. *> Specifies options for computing all or part of the matrix
  82. *> V**T:
  83. *> = 'V': the first min(m,n) rows of V**T (the right singular
  84. *> vectors) or as specified by RANGE are returned in
  85. *> the array VT;
  86. *> = 'N': no rows of V**T (no right singular vectors) are
  87. *> computed.
  88. *> \endverbatim
  89. *>
  90. *> \param[in] RANGE
  91. *> \verbatim
  92. *> RANGE is CHARACTER*1
  93. *> = 'A': all singular values will be found.
  94. *> = 'V': all singular values in the half-open interval (VL,VU]
  95. *> will be found.
  96. *> = 'I': the IL-th through IU-th singular values will be found.
  97. *> \endverbatim
  98. *>
  99. *> \param[in] M
  100. *> \verbatim
  101. *> M is INTEGER
  102. *> The number of rows of the input matrix A. M >= 0.
  103. *> \endverbatim
  104. *>
  105. *> \param[in] N
  106. *> \verbatim
  107. *> N is INTEGER
  108. *> The number of columns of the input matrix A. N >= 0.
  109. *> \endverbatim
  110. *>
  111. *> \param[in,out] A
  112. *> \verbatim
  113. *> A is COMPLEX array, dimension (LDA,N)
  114. *> On entry, the M-by-N matrix A.
  115. *> On exit, the contents of A are destroyed.
  116. *> \endverbatim
  117. *>
  118. *> \param[in] LDA
  119. *> \verbatim
  120. *> LDA is INTEGER
  121. *> The leading dimension of the array A. LDA >= max(1,M).
  122. *> \endverbatim
  123. *>
  124. *> \param[in] VL
  125. *> \verbatim
  126. *> VL is REAL
  127. *> VL >=0.
  128. *> \endverbatim
  129. *>
  130. *> \param[in] VU
  131. *> \verbatim
  132. *> VU is REAL
  133. *> If RANGE='V', the lower and upper bounds of the interval to
  134. *> be searched for singular values. VU > VL.
  135. *> Not referenced if RANGE = 'A' or 'I'.
  136. *> \endverbatim
  137. *>
  138. *> \param[in] IL
  139. *> \verbatim
  140. *> IL is INTEGER
  141. *> \endverbatim
  142. *>
  143. *> \param[in] IU
  144. *> \verbatim
  145. *> IU is INTEGER
  146. *> If RANGE='I', the indices (in ascending order) of the
  147. *> smallest and largest singular values to be returned.
  148. *> 1 <= IL <= IU <= min(M,N), if min(M,N) > 0.
  149. *> Not referenced if RANGE = 'A' or 'V'.
  150. *> \endverbatim
  151. *>
  152. *> \param[out] NS
  153. *> \verbatim
  154. *> NS is INTEGER
  155. *> The total number of singular values found,
  156. *> 0 <= NS <= min(M,N).
  157. *> If RANGE = 'A', NS = min(M,N); if RANGE = 'I', NS = IU-IL+1.
  158. *> \endverbatim
  159. *>
  160. *> \param[out] S
  161. *> \verbatim
  162. *> S is REAL array, dimension (min(M,N))
  163. *> The singular values of A, sorted so that S(i) >= S(i+1).
  164. *> \endverbatim
  165. *>
  166. *> \param[out] U
  167. *> \verbatim
  168. *> U is COMPLEX array, dimension (LDU,UCOL)
  169. *> If JOBU = 'V', U contains columns of U (the left singular
  170. *> vectors, stored columnwise) as specified by RANGE; if
  171. *> JOBU = 'N', U is not referenced.
  172. *> Note: The user must ensure that UCOL >= NS; if RANGE = 'V',
  173. *> the exact value of NS is not known ILQFin advance and an upper
  174. *> bound must be used.
  175. *> \endverbatim
  176. *>
  177. *> \param[in] LDU
  178. *> \verbatim
  179. *> LDU is INTEGER
  180. *> The leading dimension of the array U. LDU >= 1; if
  181. *> JOBU = 'V', LDU >= M.
  182. *> \endverbatim
  183. *>
  184. *> \param[out] VT
  185. *> \verbatim
  186. *> VT is COMPLEX array, dimension (LDVT,N)
  187. *> If JOBVT = 'V', VT contains the rows of V**T (the right singular
  188. *> vectors, stored rowwise) as specified by RANGE; if JOBVT = 'N',
  189. *> VT is not referenced.
  190. *> Note: The user must ensure that LDVT >= NS; if RANGE = 'V',
  191. *> the exact value of NS is not known in advance and an upper
  192. *> bound must be used.
  193. *> \endverbatim
  194. *>
  195. *> \param[in] LDVT
  196. *> \verbatim
  197. *> LDVT is INTEGER
  198. *> The leading dimension of the array VT. LDVT >= 1; if
  199. *> JOBVT = 'V', LDVT >= NS (see above).
  200. *> \endverbatim
  201. *>
  202. *> \param[out] WORK
  203. *> \verbatim
  204. *> WORK is COMPLEX array, dimension (MAX(1,LWORK))
  205. *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK;
  206. *> \endverbatim
  207. *>
  208. *> \param[in] LWORK
  209. *> \verbatim
  210. *> LWORK is INTEGER
  211. *> The dimension of the array WORK.
  212. *> LWORK >= MAX(1,MIN(M,N)*(MIN(M,N)+4)) for the paths (see
  213. *> comments inside the code):
  214. *> - PATH 1 (M much larger than N)
  215. *> - PATH 1t (N much larger than M)
  216. *> LWORK >= MAX(1,MIN(M,N)*2+MAX(M,N)) for the other paths.
  217. *> For good performance, LWORK should generally be larger.
  218. *>
  219. *> If LWORK = -1, then a workspace query is assumed; the routine
  220. *> only calculates the optimal size of the WORK array, returns
  221. *> this value as the first entry of the WORK array, and no error
  222. *> message related to LWORK is issued by XERBLA.
  223. *> \endverbatim
  224. *>
  225. *> \param[out] RWORK
  226. *> \verbatim
  227. *> RWORK is REAL array, dimension (MAX(1,LRWORK))
  228. *> LRWORK >= MIN(M,N)*(MIN(M,N)*2+15*MIN(M,N)).
  229. *> \endverbatim
  230. *>
  231. *> \param[out] IWORK
  232. *> \verbatim
  233. *> IWORK is INTEGER array, dimension (12*MIN(M,N))
  234. *> If INFO = 0, the first NS elements of IWORK are zero. If INFO > 0,
  235. *> then IWORK contains the indices of the eigenvectors that failed
  236. *> to converge in SBDSVDX/SSTEVX.
  237. *> \endverbatim
  238. *>
  239. *> \param[out] INFO
  240. *> \verbatim
  241. *> INFO is INTEGER
  242. *> = 0: successful exit
  243. *> < 0: if INFO = -i, the i-th argument had an illegal value
  244. *> > 0: if INFO = i, then i eigenvectors failed to converge
  245. *> in SBDSVDX/SSTEVX.
  246. *> if INFO = N*2 + 1, an internal error occurred in
  247. *> SBDSVDX
  248. *> \endverbatim
  249. *
  250. * Authors:
  251. * ========
  252. *
  253. *> \author Univ. of Tennessee
  254. *> \author Univ. of California Berkeley
  255. *> \author Univ. of Colorado Denver
  256. *> \author NAG Ltd.
  257. *
  258. *> \date November 2015
  259. *
  260. *> \ingroup complexGEsing
  261. *
  262. * =====================================================================
  263. SUBROUTINE CGESVDX( JOBU, JOBVT, RANGE, M, N, A, LDA, VL, VU,
  264. $ IL, IU, NS, S, U, LDU, VT, LDVT, WORK,
  265. $ LWORK, RWORK, IWORK, INFO )
  266. *
  267. * -- LAPACK driver routine (version 3.6.0) --
  268. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  269. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  270. * November 2015
  271. *
  272. * .. Scalar Arguments ..
  273. CHARACTER JOBU, JOBVT, RANGE
  274. INTEGER IL, INFO, IU, LDA, LDU, LDVT, LWORK, M, N, NS
  275. REAL VL, VU
  276. * ..
  277. * .. Array Arguments ..
  278. INTEGER IWORK( * )
  279. REAL S( * ), RWORK( * )
  280. COMPLEX A( LDA, * ), U( LDU, * ), VT( LDVT, * ),
  281. $ WORK( * )
  282. * ..
  283. *
  284. * =====================================================================
  285. *
  286. * .. Parameters ..
  287. COMPLEX CZERO, CONE
  288. PARAMETER ( CZERO = ( 0.0E0, 0.0E0 ),
  289. $ CONE = ( 1.0E0, 0.0E0 ) )
  290. REAL ZERO, ONE
  291. PARAMETER ( ZERO = 0.0E0, ONE = 1.0E0 )
  292. * ..
  293. * .. Local Scalars ..
  294. CHARACTER JOBZ, RNGTGK
  295. LOGICAL ALLS, INDS, LQUERY, VALS, WANTU, WANTVT
  296. INTEGER I, ID, IE, IERR, ILQF, ILTGK, IQRF, ISCL,
  297. $ ITAU, ITAUP, ITAUQ, ITEMP, ITGKZ, IUTGK,
  298. $ J, K, MAXWRK, MINMN, MINWRK, MNTHR
  299. REAL ABSTOL, ANRM, BIGNUM, EPS, SMLNUM
  300. * ..
  301. * .. Local Arrays ..
  302. REAL DUM( 1 )
  303. * ..
  304. * .. External Subroutines ..
  305. EXTERNAL CGEBRD, CGELQF, CGEQRF, CLASCL, CLASET,
  306. $ SLASCL, XERBLA
  307. * ..
  308. * .. External Functions ..
  309. LOGICAL LSAME
  310. INTEGER ILAENV
  311. REAL SLAMCH, CLANGE
  312. EXTERNAL LSAME, ILAENV, SLAMCH, CLANGE
  313. * ..
  314. * .. Intrinsic Functions ..
  315. INTRINSIC MAX, MIN, SQRT
  316. * ..
  317. * .. Executable Statements ..
  318. *
  319. * Test the input arguments.
  320. *
  321. NS = 0
  322. INFO = 0
  323. ABSTOL = 2*SLAMCH('S')
  324. LQUERY = ( LWORK.EQ.-1 )
  325. MINMN = MIN( M, N )
  326. WANTU = LSAME( JOBU, 'V' )
  327. WANTVT = LSAME( JOBVT, 'V' )
  328. IF( WANTU .OR. WANTVT ) THEN
  329. JOBZ = 'V'
  330. ELSE
  331. JOBZ = 'N'
  332. END IF
  333. ALLS = LSAME( RANGE, 'A' )
  334. VALS = LSAME( RANGE, 'V' )
  335. INDS = LSAME( RANGE, 'I' )
  336. *
  337. INFO = 0
  338. IF( .NOT.LSAME( JOBU, 'V' ) .AND.
  339. $ .NOT.LSAME( JOBU, 'N' ) ) THEN
  340. INFO = -1
  341. ELSE IF( .NOT.LSAME( JOBVT, 'V' ) .AND.
  342. $ .NOT.LSAME( JOBVT, 'N' ) ) THEN
  343. INFO = -2
  344. ELSE IF( .NOT.( ALLS .OR. VALS .OR. INDS ) ) THEN
  345. INFO = -3
  346. ELSE IF( M.LT.0 ) THEN
  347. INFO = -4
  348. ELSE IF( N.LT.0 ) THEN
  349. INFO = -5
  350. ELSE IF( M.GT.LDA ) THEN
  351. INFO = -7
  352. ELSE IF( MINMN.GT.0 ) THEN
  353. IF( VALS ) THEN
  354. IF( VL.LT.ZERO ) THEN
  355. INFO = -8
  356. ELSE IF( VU.LE.VL ) THEN
  357. INFO = -9
  358. END IF
  359. ELSE IF( INDS ) THEN
  360. IF( IL.LT.1 .OR. IL.GT.MAX( 1, MINMN ) ) THEN
  361. INFO = -10
  362. ELSE IF( IU.LT.MIN( MINMN, IL ) .OR. IU.GT.MINMN ) THEN
  363. INFO = -11
  364. END IF
  365. END IF
  366. IF( INFO.EQ.0 ) THEN
  367. IF( WANTU .AND. LDU.LT.M ) THEN
  368. INFO = -15
  369. ELSE IF( WANTVT .AND. LDVT.LT.MINMN ) THEN
  370. INFO = -16
  371. END IF
  372. END IF
  373. END IF
  374. *
  375. * Compute workspace
  376. * (Note: Comments in the code beginning "Workspace:" describe the
  377. * minimal amount of workspace needed at that point in the code,
  378. * as well as the preferred amount for good performance.
  379. * NB refers to the optimal block size for the immediately
  380. * following subroutine, as returned by ILAENV.)
  381. *
  382. IF( INFO.EQ.0 ) THEN
  383. MINWRK = 1
  384. MAXWRK = 1
  385. IF( MINMN.GT.0 ) THEN
  386. IF( M.GE.N ) THEN
  387. MNTHR = ILAENV( 6, 'CGESVD', JOBU // JOBVT, M, N, 0, 0 )
  388. IF( M.GE.MNTHR ) THEN
  389. *
  390. * Path 1 (M much larger than N)
  391. *
  392. MAXWRK = N + N*
  393. $ ILAENV( 1, 'SGEQRF', ' ', M, N, -1, -1 )
  394. MAXWRK = MAX( MAXWRK, N*N + N + 2*N*
  395. $ ILAENV( 1, 'SGEBRD', ' ', N, N, -1, -1 ) )
  396. MINWRK = N*(N+4)
  397. ELSE
  398. *
  399. * Path 2 (M at least N, but not much larger)
  400. *
  401. MAXWRK = 2*N + ( M+N )*
  402. $ ILAENV( 1, 'CGEBRD', ' ', M, N, -1, -1 )
  403. MINWRK = 2*N + M
  404. END IF
  405. ELSE
  406. MNTHR = ILAENV( 6, 'CGESVD', JOBU // JOBVT, M, N, 0, 0 )
  407. IF( N.GE.MNTHR ) THEN
  408. *
  409. * Path 1t (N much larger than M)
  410. *
  411. MAXWRK = M + M*
  412. $ ILAENV( 1, 'CGELQF', ' ', M, N, -1, -1 )
  413. MAXWRK = MAX( MAXWRK, M*M + M + 2*M*
  414. $ ILAENV( 1, 'CGEBRD', ' ', M, M, -1, -1 ) )
  415. MINWRK = M*(M+4)
  416. ELSE
  417. *
  418. * Path 2t (N greater than M, but not much larger)
  419. *
  420. MAXWRK = M*(M*2+19) + ( M+N )*
  421. $ ILAENV( 1, 'CGEBRD', ' ', M, N, -1, -1 )
  422. MINWRK = 2*M + N
  423. END IF
  424. END IF
  425. END IF
  426. MAXWRK = MAX( MAXWRK, MINWRK )
  427. WORK( 1 ) = CMPLX( REAL( MAXWRK ), ZERO )
  428. *
  429. IF( LWORK.LT.MINWRK .AND. .NOT.LQUERY ) THEN
  430. INFO = -19
  431. END IF
  432. END IF
  433. *
  434. IF( INFO.NE.0 ) THEN
  435. CALL XERBLA( 'CGESVDX', -INFO )
  436. RETURN
  437. ELSE IF( LQUERY ) THEN
  438. RETURN
  439. END IF
  440. *
  441. * Quick return if possible
  442. *
  443. IF( M.EQ.0 .OR. N.EQ.0 ) THEN
  444. RETURN
  445. END IF
  446. *
  447. * Set singular values indices accord to RANGE='A'.
  448. *
  449. ALLS = LSAME( RANGE, 'A' )
  450. INDS = LSAME( RANGE, 'I' )
  451. IF( ALLS ) THEN
  452. RNGTGK = 'I'
  453. ILTGK = 1
  454. IUTGK = MIN( M, N )
  455. ELSE IF( INDS ) THEN
  456. RNGTGK = 'I'
  457. ILTGK = IL
  458. IUTGK = IU
  459. ELSE
  460. RNGTGK = 'V'
  461. ILTGK = 0
  462. IUTGK = 0
  463. END IF
  464. *
  465. * Get machine constants
  466. *
  467. EPS = SLAMCH( 'P' )
  468. SMLNUM = SQRT( SLAMCH( 'S' ) ) / EPS
  469. BIGNUM = ONE / SMLNUM
  470. *
  471. * Scale A if max element outside range [SMLNUM,BIGNUM]
  472. *
  473. ANRM = CLANGE( 'M', M, N, A, LDA, DUM )
  474. ISCL = 0
  475. IF( ANRM.GT.ZERO .AND. ANRM.LT.SMLNUM ) THEN
  476. ISCL = 1
  477. CALL CLASCL( 'G', 0, 0, ANRM, SMLNUM, M, N, A, LDA, INFO )
  478. ELSE IF( ANRM.GT.BIGNUM ) THEN
  479. ISCL = 1
  480. CALL CLASCL( 'G', 0, 0, ANRM, BIGNUM, M, N, A, LDA, INFO )
  481. END IF
  482. *
  483. IF( M.GE.N ) THEN
  484. *
  485. * A has at least as many rows as columns. If A has sufficiently
  486. * more rows than columns, first reduce A using the QR
  487. * decomposition.
  488. *
  489. IF( M.GE.MNTHR ) THEN
  490. *
  491. * Path 1 (M much larger than N):
  492. * A = Q * R = Q * ( QB * B * PB**T )
  493. * = Q * ( QB * ( UB * S * VB**T ) * PB**T )
  494. * U = Q * QB * UB; V**T = VB**T * PB**T
  495. *
  496. * Compute A=Q*R
  497. * (Workspace: need 2*N, prefer N+N*NB)
  498. *
  499. ITAU = 1
  500. ITEMP = ITAU + N
  501. CALL CGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( ITEMP ),
  502. $ LWORK-ITEMP+1, INFO )
  503. *
  504. * Copy R into WORK and bidiagonalize it:
  505. * (Workspace: need N*N+3*N, prefer N*N+N+2*N*NB)
  506. *
  507. IQRF = ITEMP
  508. ITAUQ = ITEMP + N*N
  509. ITAUP = ITAUQ + N
  510. ITEMP = ITAUP + N
  511. ID = 1
  512. IE = ID + N
  513. ITGKZ = IE + N
  514. CALL CLACPY( 'U', N, N, A, LDA, WORK( IQRF ), N )
  515. CALL CLASET( 'L', N-1, N-1, CZERO, CZERO,
  516. $ WORK( IQRF+1 ), N )
  517. CALL CGEBRD( N, N, WORK( IQRF ), N, RWORK( ID ),
  518. $ RWORK( IE ), WORK( ITAUQ ), WORK( ITAUP ),
  519. $ WORK( ITEMP ), LWORK-ITEMP+1, INFO )
  520. ITEMP = ITGKZ + N*(N*2+1)
  521. *
  522. * Solve eigenvalue problem TGK*Z=Z*S.
  523. * (Workspace: need 2*N*N+14*N)
  524. *
  525. CALL SBDSVDX( 'U', JOBZ, RNGTGK, N, RWORK( ID ),
  526. $ RWORK( IE ), VL, VU, ILTGK, IUTGK, NS, S,
  527. $ RWORK( ITGKZ ), N*2, RWORK( ITEMP ),
  528. $ IWORK, INFO)
  529. *
  530. * If needed, compute left singular vectors.
  531. *
  532. IF( WANTU ) THEN
  533. K = ITGKZ
  534. DO I = 1, NS
  535. DO J = 1, N
  536. U( J, I ) = CMPLX( RWORK( K ), ZERO )
  537. K = K + 1
  538. END DO
  539. K = K + N
  540. END DO
  541. CALL CLASET( 'A', M-N, N, CZERO, CZERO, U( N+1,1 ), LDU )
  542. *
  543. * Call CUNMBR to compute QB*UB.
  544. * (Workspace in WORK( ITEMP ): need N, prefer N*NB)
  545. *
  546. CALL CUNMBR( 'Q', 'L', 'N', N, NS, N, WORK( IQRF ), N,
  547. $ WORK( ITAUQ ), U, LDU, WORK( ITEMP ),
  548. $ LWORK-ITEMP+1, INFO )
  549. *
  550. * Call CUNMQR to compute Q*(QB*UB).
  551. * (Workspace in WORK( ITEMP ): need N, prefer N*NB)
  552. *
  553. CALL CUNMQR( 'L', 'N', M, NS, N, A, LDA,
  554. $ WORK( ITAU ), U, LDU, WORK( ITEMP ),
  555. $ LWORK-ITEMP+1, INFO )
  556. END IF
  557. *
  558. * If needed, compute right singular vectors.
  559. *
  560. IF( WANTVT) THEN
  561. K = ITGKZ + N
  562. DO I = 1, NS
  563. DO J = 1, N
  564. VT( I, J ) = CMPLX( RWORK( K ), ZERO )
  565. K = K + 1
  566. END DO
  567. K = K + N
  568. END DO
  569. *
  570. * Call CUNMBR to compute VB**T * PB**T
  571. * (Workspace in WORK( ITEMP ): need N, prefer N*NB)
  572. *
  573. CALL CUNMBR( 'P', 'R', 'C', NS, N, N, WORK( IQRF ), N,
  574. $ WORK( ITAUP ), VT, LDVT, WORK( ITEMP ),
  575. $ LWORK-ITEMP+1, INFO )
  576. END IF
  577. ELSE
  578. *
  579. * Path 2 (M at least N, but not much larger)
  580. * Reduce A to bidiagonal form without QR decomposition
  581. * A = QB * B * PB**T = QB * ( UB * S * VB**T ) * PB**T
  582. * U = QB * UB; V**T = VB**T * PB**T
  583. *
  584. * Bidiagonalize A
  585. * (Workspace: need 2*N+M, prefer 2*N+(M+N)*NB)
  586. *
  587. ITAUQ = 1
  588. ITAUP = ITAUQ + N
  589. ITEMP = ITAUP + N
  590. ID = 1
  591. IE = ID + N
  592. ITGKZ = IE + N
  593. CALL CGEBRD( M, N, A, LDA, RWORK( ID ), RWORK( IE ),
  594. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( ITEMP ),
  595. $ LWORK-ITEMP+1, INFO )
  596. ITEMP = ITGKZ + N*(N*2+1)
  597. *
  598. * Solve eigenvalue problem TGK*Z=Z*S.
  599. * (Workspace: need 2*N*N+14*N)
  600. *
  601. CALL SBDSVDX( 'U', JOBZ, RNGTGK, N, RWORK( ID ),
  602. $ RWORK( IE ), VL, VU, ILTGK, IUTGK, NS, S,
  603. $ RWORK( ITGKZ ), N*2, RWORK( ITEMP ),
  604. $ IWORK, INFO)
  605. *
  606. * If needed, compute left singular vectors.
  607. *
  608. IF( WANTU ) THEN
  609. K = ITGKZ
  610. DO I = 1, NS
  611. DO J = 1, N
  612. U( J, I ) = CMPLX( RWORK( K ), ZERO )
  613. K = K + 1
  614. END DO
  615. K = K + N
  616. END DO
  617. CALL CLASET( 'A', M-N, N, CZERO, CZERO, U( N+1,1 ), LDU )
  618. *
  619. * Call CUNMBR to compute QB*UB.
  620. * (Workspace in WORK( ITEMP ): need N, prefer N*NB)
  621. *
  622. CALL CUNMBR( 'Q', 'L', 'N', M, NS, N, A, LDA,
  623. $ WORK( ITAUQ ), U, LDU, WORK( ITEMP ),
  624. $ LWORK-ITEMP+1, IERR )
  625. END IF
  626. *
  627. * If needed, compute right singular vectors.
  628. *
  629. IF( WANTVT) THEN
  630. K = ITGKZ + N
  631. DO I = 1, NS
  632. DO J = 1, N
  633. VT( I, J ) = CMPLX( RWORK( K ), ZERO )
  634. K = K + 1
  635. END DO
  636. K = K + N
  637. END DO
  638. *
  639. * Call CUNMBR to compute VB**T * PB**T
  640. * (Workspace in WORK( ITEMP ): need N, prefer N*NB)
  641. *
  642. CALL CUNMBR( 'P', 'R', 'C', NS, N, N, A, LDA,
  643. $ WORK( ITAUP ), VT, LDVT, WORK( ITEMP ),
  644. $ LWORK-ITEMP+1, IERR )
  645. END IF
  646. END IF
  647. ELSE
  648. *
  649. * A has more columns than rows. If A has sufficiently more
  650. * columns than rows, first reduce A using the LQ decomposition.
  651. *
  652. IF( N.GE.MNTHR ) THEN
  653. *
  654. * Path 1t (N much larger than M):
  655. * A = L * Q = ( QB * B * PB**T ) * Q
  656. * = ( QB * ( UB * S * VB**T ) * PB**T ) * Q
  657. * U = QB * UB ; V**T = VB**T * PB**T * Q
  658. *
  659. * Compute A=L*Q
  660. * (Workspace: need 2*M, prefer M+M*NB)
  661. *
  662. ITAU = 1
  663. ITEMP = ITAU + M
  664. CALL CGELQF( M, N, A, LDA, WORK( ITAU ), WORK( ITEMP ),
  665. $ LWORK-ITEMP+1, INFO )
  666. * Copy L into WORK and bidiagonalize it:
  667. * (Workspace in WORK( ITEMP ): need M*M+3*M, prefer M*M+M+2*M*NB)
  668. *
  669. ILQF = ITEMP
  670. ITAUQ = ILQF + M*M
  671. ITAUP = ITAUQ + M
  672. ITEMP = ITAUP + M
  673. ID = 1
  674. IE = ID + M
  675. ITGKZ = IE + M
  676. CALL CLACPY( 'L', M, M, A, LDA, WORK( ILQF ), M )
  677. CALL CLASET( 'U', M-1, M-1, CZERO, CZERO,
  678. $ WORK( ILQF+M ), M )
  679. CALL CGEBRD( M, M, WORK( ILQF ), M, RWORK( ID ),
  680. $ RWORK( IE ), WORK( ITAUQ ), WORK( ITAUP ),
  681. $ WORK( ITEMP ), LWORK-ITEMP+1, INFO )
  682. ITEMP = ITGKZ + M*(M*2+1)
  683. *
  684. * Solve eigenvalue problem TGK*Z=Z*S.
  685. * (Workspace: need 2*M*M+14*M)
  686. *
  687. CALL SBDSVDX( 'U', JOBZ, RNGTGK, M, RWORK( ID ),
  688. $ RWORK( IE ), VL, VU, ILTGK, IUTGK, NS, S,
  689. $ RWORK( ITGKZ ), M*2, RWORK( ITEMP ),
  690. $ IWORK, INFO)
  691. *
  692. * If needed, compute left singular vectors.
  693. *
  694. IF( WANTU ) THEN
  695. K = ITGKZ
  696. DO I = 1, NS
  697. DO J = 1, M
  698. U( J, I ) = CMPLX( RWORK( K ), ZERO )
  699. K = K + 1
  700. END DO
  701. K = K + M
  702. END DO
  703. *
  704. * Call CUNMBR to compute QB*UB.
  705. * (Workspace in WORK( ITEMP ): need M, prefer M*NB)
  706. *
  707. CALL CUNMBR( 'Q', 'L', 'N', M, NS, M, WORK( ILQF ), M,
  708. $ WORK( ITAUQ ), U, LDU, WORK( ITEMP ),
  709. $ LWORK-ITEMP+1, INFO )
  710. END IF
  711. *
  712. * If needed, compute right singular vectors.
  713. *
  714. IF( WANTVT) THEN
  715. K = ITGKZ + M
  716. DO I = 1, NS
  717. DO J = 1, M
  718. VT( I, J ) = CMPLX( RWORK( K ), ZERO )
  719. K = K + 1
  720. END DO
  721. K = K + M
  722. END DO
  723. CALL CLASET( 'A', M, N-M, CZERO, CZERO,
  724. $ VT( 1,M+1 ), LDVT )
  725. *
  726. * Call CUNMBR to compute (VB**T)*(PB**T)
  727. * (Workspace in WORK( ITEMP ): need M, prefer M*NB)
  728. *
  729. CALL CUNMBR( 'P', 'R', 'C', NS, M, M, WORK( ILQF ), M,
  730. $ WORK( ITAUP ), VT, LDVT, WORK( ITEMP ),
  731. $ LWORK-ITEMP+1, INFO )
  732. *
  733. * Call CUNMLQ to compute ((VB**T)*(PB**T))*Q.
  734. * (Workspace in WORK( ITEMP ): need M, prefer M*NB)
  735. *
  736. CALL CUNMLQ( 'R', 'N', NS, N, M, A, LDA,
  737. $ WORK( ITAU ), VT, LDVT, WORK( ITEMP ),
  738. $ LWORK-ITEMP+1, INFO )
  739. END IF
  740. ELSE
  741. *
  742. * Path 2t (N greater than M, but not much larger)
  743. * Reduce to bidiagonal form without LQ decomposition
  744. * A = QB * B * PB**T = QB * ( UB * S * VB**T ) * PB**T
  745. * U = QB * UB; V**T = VB**T * PB**T
  746. *
  747. * Bidiagonalize A
  748. * (Workspace: need 2*M+N, prefer 2*M+(M+N)*NB)
  749. *
  750. ITAUQ = 1
  751. ITAUP = ITAUQ + M
  752. ITEMP = ITAUP + M
  753. ID = 1
  754. IE = ID + M
  755. ITGKZ = IE + M
  756. CALL CGEBRD( M, N, A, LDA, RWORK( ID ), RWORK( IE ),
  757. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( ITEMP ),
  758. $ LWORK-ITEMP+1, INFO )
  759. ITEMP = ITGKZ + M*(M*2+1)
  760. *
  761. * Solve eigenvalue problem TGK*Z=Z*S.
  762. * (Workspace: need 2*M*M+14*M)
  763. *
  764. CALL SBDSVDX( 'L', JOBZ, RNGTGK, M, RWORK( ID ),
  765. $ RWORK( IE ), VL, VU, ILTGK, IUTGK, NS, S,
  766. $ RWORK( ITGKZ ), M*2, RWORK( ITEMP ),
  767. $ IWORK, INFO)
  768. *
  769. * If needed, compute left singular vectors.
  770. *
  771. IF( WANTU ) THEN
  772. K = ITGKZ
  773. DO I = 1, NS
  774. DO J = 1, M
  775. U( J, I ) = CMPLX( RWORK( K ), ZERO )
  776. K = K + 1
  777. END DO
  778. K = K + M
  779. END DO
  780. *
  781. * Call CUNMBR to compute QB*UB.
  782. * (Workspace in WORK( ITEMP ): need M, prefer M*NB)
  783. *
  784. CALL CUNMBR( 'Q', 'L', 'N', M, NS, N, A, LDA,
  785. $ WORK( ITAUQ ), U, LDU, WORK( ITEMP ),
  786. $ LWORK-ITEMP+1, INFO )
  787. END IF
  788. *
  789. * If needed, compute right singular vectors.
  790. *
  791. IF( WANTVT) THEN
  792. K = ITGKZ + M
  793. DO I = 1, NS
  794. DO J = 1, M
  795. VT( I, J ) = CMPLX( RWORK( K ), ZERO )
  796. K = K + 1
  797. END DO
  798. K = K + M
  799. END DO
  800. CALL CLASET( 'A', M, N-M, CZERO, CZERO,
  801. $ VT( 1,M+1 ), LDVT )
  802. *
  803. * Call CUNMBR to compute VB**T * PB**T
  804. * (Workspace in WORK( ITEMP ): need M, prefer M*NB)
  805. *
  806. CALL CUNMBR( 'P', 'R', 'C', NS, N, M, A, LDA,
  807. $ WORK( ITAUP ), VT, LDVT, WORK( ITEMP ),
  808. $ LWORK-ITEMP+1, INFO )
  809. END IF
  810. END IF
  811. END IF
  812. *
  813. * Undo scaling if necessary
  814. *
  815. IF( ISCL.EQ.1 ) THEN
  816. IF( ANRM.GT.BIGNUM )
  817. $ CALL SLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN, 1,
  818. $ S, MINMN, INFO )
  819. IF( ANRM.LT.SMLNUM )
  820. $ CALL SLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN, 1,
  821. $ S, MINMN, INFO )
  822. END IF
  823. *
  824. * Return optimal workspace in WORK(1)
  825. *
  826. WORK( 1 ) = CMPLX( REAL( MAXWRK ), ZERO )
  827. *
  828. RETURN
  829. *
  830. * End of CGESVDX
  831. *
  832. END