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.

dgesvdx.f 27 kB

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