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.

sgesvdx.f 28 kB

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