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.

zgesvdx.f 28 kB

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