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.

dgesdd.f 60 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545
  1. *> \brief \b DGESDD
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download DGESDD + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgesdd.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgesdd.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgesdd.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE DGESDD( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
  22. * WORK, LWORK, IWORK, INFO )
  23. *
  24. * .. Scalar Arguments ..
  25. * CHARACTER JOBZ
  26. * INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
  27. * ..
  28. * .. Array Arguments ..
  29. * INTEGER IWORK( * )
  30. * DOUBLE PRECISION A( LDA, * ), S( * ), U( LDU, * ),
  31. * $ VT( LDVT, * ), WORK( * )
  32. * ..
  33. *
  34. *
  35. *> \par Purpose:
  36. * =============
  37. *>
  38. *> \verbatim
  39. *>
  40. *> DGESDD computes the singular value decomposition (SVD) of a real
  41. *> M-by-N matrix A, optionally computing the left and right singular
  42. *> vectors. If singular vectors are desired, it uses a
  43. *> divide-and-conquer algorithm.
  44. *>
  45. *> 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. *> Note that the routine returns VT = V**T, not V.
  57. *>
  58. *> \endverbatim
  59. *
  60. * Arguments:
  61. * ==========
  62. *
  63. *> \param[in] JOBZ
  64. *> \verbatim
  65. *> JOBZ is CHARACTER*1
  66. *> Specifies options for computing all or part of the matrix U:
  67. *> = 'A': all M columns of U and all N rows of V**T are
  68. *> returned in the arrays U and VT;
  69. *> = 'S': the first min(M,N) columns of U and the first
  70. *> min(M,N) rows of V**T are returned in the arrays U
  71. *> and VT;
  72. *> = 'O': If M >= N, the first N columns of U are overwritten
  73. *> on the array A and all rows of V**T are returned in
  74. *> the array VT;
  75. *> otherwise, all columns of U are returned in the
  76. *> array U and the first M rows of V**T are overwritten
  77. *> in the array A;
  78. *> = 'N': no columns of U or rows of V**T are computed.
  79. *> \endverbatim
  80. *>
  81. *> \param[in] M
  82. *> \verbatim
  83. *> M is INTEGER
  84. *> The number of rows of the input matrix A. M >= 0.
  85. *> \endverbatim
  86. *>
  87. *> \param[in] N
  88. *> \verbatim
  89. *> N is INTEGER
  90. *> The number of columns of the input matrix A. N >= 0.
  91. *> \endverbatim
  92. *>
  93. *> \param[in,out] A
  94. *> \verbatim
  95. *> A is DOUBLE PRECISION array, dimension (LDA,N)
  96. *> On entry, the M-by-N matrix A.
  97. *> On exit,
  98. *> if JOBZ = 'O', A is overwritten with the first N columns
  99. *> of U (the left singular vectors, stored
  100. *> columnwise) if M >= N;
  101. *> A is overwritten with the first M rows
  102. *> of V**T (the right singular vectors, stored
  103. *> rowwise) otherwise.
  104. *> if JOBZ .ne. 'O', the contents of A are destroyed.
  105. *> \endverbatim
  106. *>
  107. *> \param[in] LDA
  108. *> \verbatim
  109. *> LDA is INTEGER
  110. *> The leading dimension of the array A. LDA >= max(1,M).
  111. *> \endverbatim
  112. *>
  113. *> \param[out] S
  114. *> \verbatim
  115. *> S is DOUBLE PRECISION array, dimension (min(M,N))
  116. *> The singular values of A, sorted so that S(i) >= S(i+1).
  117. *> \endverbatim
  118. *>
  119. *> \param[out] U
  120. *> \verbatim
  121. *> U is DOUBLE PRECISION array, dimension (LDU,UCOL)
  122. *> UCOL = M if JOBZ = 'A' or JOBZ = 'O' and M < N;
  123. *> UCOL = min(M,N) if JOBZ = 'S'.
  124. *> If JOBZ = 'A' or JOBZ = 'O' and M < N, U contains the M-by-M
  125. *> orthogonal matrix U;
  126. *> if JOBZ = 'S', U contains the first min(M,N) columns of U
  127. *> (the left singular vectors, stored columnwise);
  128. *> if JOBZ = 'O' and M >= N, or JOBZ = 'N', U is not referenced.
  129. *> \endverbatim
  130. *>
  131. *> \param[in] LDU
  132. *> \verbatim
  133. *> LDU is INTEGER
  134. *> The leading dimension of the array U. LDU >= 1; if
  135. *> JOBZ = 'S' or 'A' or JOBZ = 'O' and M < N, LDU >= M.
  136. *> \endverbatim
  137. *>
  138. *> \param[out] VT
  139. *> \verbatim
  140. *> VT is DOUBLE PRECISION array, dimension (LDVT,N)
  141. *> If JOBZ = 'A' or JOBZ = 'O' and M >= N, VT contains the
  142. *> N-by-N orthogonal matrix V**T;
  143. *> if JOBZ = 'S', VT contains the first min(M,N) rows of
  144. *> V**T (the right singular vectors, stored rowwise);
  145. *> if JOBZ = 'O' and M < N, or JOBZ = 'N', VT is not referenced.
  146. *> \endverbatim
  147. *>
  148. *> \param[in] LDVT
  149. *> \verbatim
  150. *> LDVT is INTEGER
  151. *> The leading dimension of the array VT. LDVT >= 1;
  152. *> if JOBZ = 'A' or JOBZ = 'O' and M >= N, LDVT >= N;
  153. *> if JOBZ = 'S', LDVT >= min(M,N).
  154. *> \endverbatim
  155. *>
  156. *> \param[out] WORK
  157. *> \verbatim
  158. *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
  159. *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK;
  160. *> \endverbatim
  161. *>
  162. *> \param[in] LWORK
  163. *> \verbatim
  164. *> LWORK is INTEGER
  165. *> The dimension of the array WORK. LWORK >= 1.
  166. *> If LWORK = -1, a workspace query is assumed. The optimal
  167. *> size for the WORK array is calculated and stored in WORK(1),
  168. *> and no other work except argument checking is performed.
  169. *>
  170. *> Let mx = max(M,N) and mn = min(M,N).
  171. *> If JOBZ = 'N', LWORK >= 3*mn + max( mx, 7*mn ).
  172. *> If JOBZ = 'O', LWORK >= 3*mn + max( mx, 5*mn*mn + 4*mn ).
  173. *> If JOBZ = 'S', LWORK >= 4*mn*mn + 7*mn.
  174. *> If JOBZ = 'A', LWORK >= 4*mn*mn + 6*mn + mx.
  175. *> These are not tight minimums in all cases; see comments inside code.
  176. *> For good performance, LWORK should generally be larger;
  177. *> a query is recommended.
  178. *> \endverbatim
  179. *>
  180. *> \param[out] IWORK
  181. *> \verbatim
  182. *> IWORK is INTEGER array, dimension (8*min(M,N))
  183. *> \endverbatim
  184. *>
  185. *> \param[out] INFO
  186. *> \verbatim
  187. *> INFO is INTEGER
  188. *> < 0: if INFO = -i, the i-th argument had an illegal value.
  189. *> = -4: if A had a NAN entry.
  190. *> > 0: DBDSDC did not converge, updating process failed.
  191. *> = 0: successful exit.
  192. *> \endverbatim
  193. *
  194. * Authors:
  195. * ========
  196. *
  197. *> \author Univ. of Tennessee
  198. *> \author Univ. of California Berkeley
  199. *> \author Univ. of Colorado Denver
  200. *> \author NAG Ltd.
  201. *
  202. *> \ingroup doubleGEsing
  203. *
  204. *> \par Contributors:
  205. * ==================
  206. *>
  207. *> Ming Gu and Huan Ren, Computer Science Division, University of
  208. *> California at Berkeley, USA
  209. *>
  210. * =====================================================================
  211. SUBROUTINE DGESDD( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
  212. $ WORK, LWORK, IWORK, INFO )
  213. implicit none
  214. *
  215. * -- LAPACK driver routine --
  216. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  217. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  218. *
  219. * .. Scalar Arguments ..
  220. CHARACTER JOBZ
  221. INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
  222. * ..
  223. * .. Array Arguments ..
  224. INTEGER IWORK( * )
  225. DOUBLE PRECISION A( LDA, * ), S( * ), U( LDU, * ),
  226. $ VT( LDVT, * ), WORK( * )
  227. * ..
  228. *
  229. * =====================================================================
  230. *
  231. * .. Parameters ..
  232. DOUBLE PRECISION ZERO, ONE
  233. PARAMETER ( ZERO = 0.0D0, ONE = 1.0D0 )
  234. * ..
  235. * .. Local Scalars ..
  236. LOGICAL LQUERY, WNTQA, WNTQAS, WNTQN, WNTQO, WNTQS
  237. INTEGER BDSPAC, BLK, CHUNK, I, IE, IERR, IL,
  238. $ IR, ISCL, ITAU, ITAUP, ITAUQ, IU, IVT, LDWKVT,
  239. $ LDWRKL, LDWRKR, LDWRKU, MAXWRK, MINMN, MINWRK,
  240. $ MNTHR, NWORK, WRKBL
  241. INTEGER LWORK_DGEBRD_MN, LWORK_DGEBRD_MM,
  242. $ LWORK_DGEBRD_NN, LWORK_DGELQF_MN,
  243. $ LWORK_DGEQRF_MN,
  244. $ LWORK_DORGBR_P_MM, LWORK_DORGBR_Q_NN,
  245. $ LWORK_DORGLQ_MN, LWORK_DORGLQ_NN,
  246. $ LWORK_DORGQR_MM, LWORK_DORGQR_MN,
  247. $ LWORK_DORMBR_PRT_MM, LWORK_DORMBR_QLN_MM,
  248. $ LWORK_DORMBR_PRT_MN, LWORK_DORMBR_QLN_MN,
  249. $ LWORK_DORMBR_PRT_NN, LWORK_DORMBR_QLN_NN
  250. DOUBLE PRECISION ANRM, BIGNUM, EPS, SMLNUM
  251. * ..
  252. * .. Local Arrays ..
  253. INTEGER IDUM( 1 )
  254. DOUBLE PRECISION DUM( 1 )
  255. * ..
  256. * .. External Subroutines ..
  257. EXTERNAL DBDSDC, DGEBRD, DGELQF, DGEMM, DGEQRF, DLACPY,
  258. $ DLASCL, DLASET, DORGBR, DORGLQ, DORGQR, DORMBR,
  259. $ XERBLA
  260. * ..
  261. * .. External Functions ..
  262. LOGICAL LSAME, DISNAN
  263. DOUBLE PRECISION DLAMCH, DLANGE, DROUNDUP_LWORK
  264. EXTERNAL DLAMCH, DLANGE, LSAME, DISNAN,
  265. $ DROUNDUP_LWORK
  266. * ..
  267. * .. Intrinsic Functions ..
  268. INTRINSIC INT, MAX, MIN, SQRT
  269. * ..
  270. * .. Executable Statements ..
  271. *
  272. * Test the input arguments
  273. *
  274. INFO = 0
  275. MINMN = MIN( M, N )
  276. WNTQA = LSAME( JOBZ, 'A' )
  277. WNTQS = LSAME( JOBZ, 'S' )
  278. WNTQAS = WNTQA .OR. WNTQS
  279. WNTQO = LSAME( JOBZ, 'O' )
  280. WNTQN = LSAME( JOBZ, 'N' )
  281. LQUERY = ( LWORK.EQ.-1 )
  282. *
  283. IF( .NOT.( WNTQA .OR. WNTQS .OR. WNTQO .OR. WNTQN ) ) THEN
  284. INFO = -1
  285. ELSE IF( M.LT.0 ) THEN
  286. INFO = -2
  287. ELSE IF( N.LT.0 ) THEN
  288. INFO = -3
  289. ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
  290. INFO = -5
  291. ELSE IF( LDU.LT.1 .OR. ( WNTQAS .AND. LDU.LT.M ) .OR.
  292. $ ( WNTQO .AND. M.LT.N .AND. LDU.LT.M ) ) THEN
  293. INFO = -8
  294. ELSE IF( LDVT.LT.1 .OR. ( WNTQA .AND. LDVT.LT.N ) .OR.
  295. $ ( WNTQS .AND. LDVT.LT.MINMN ) .OR.
  296. $ ( WNTQO .AND. M.GE.N .AND. LDVT.LT.N ) ) THEN
  297. INFO = -10
  298. END IF
  299. *
  300. * Compute workspace
  301. * Note: Comments in the code beginning "Workspace:" describe the
  302. * minimal amount of workspace allocated at that point in the code,
  303. * as well as the preferred amount for good performance.
  304. * NB refers to the optimal block size for the immediately
  305. * following subroutine, as returned by ILAENV.
  306. *
  307. IF( INFO.EQ.0 ) THEN
  308. MINWRK = 1
  309. MAXWRK = 1
  310. BDSPAC = 0
  311. MNTHR = INT( MINMN*11.0D0 / 6.0D0 )
  312. IF( M.GE.N .AND. MINMN.GT.0 ) THEN
  313. *
  314. * Compute space needed for DBDSDC
  315. *
  316. IF( WNTQN ) THEN
  317. * dbdsdc needs only 4*N (or 6*N for uplo=L for LAPACK <= 3.6)
  318. * keep 7*N for backwards compatibility.
  319. BDSPAC = 7*N
  320. ELSE
  321. BDSPAC = 3*N*N + 4*N
  322. END IF
  323. *
  324. * Compute space preferred for each routine
  325. CALL DGEBRD( M, N, DUM(1), M, DUM(1), DUM(1), DUM(1),
  326. $ DUM(1), DUM(1), -1, IERR )
  327. LWORK_DGEBRD_MN = INT( DUM(1) )
  328. *
  329. CALL DGEBRD( N, N, DUM(1), N, DUM(1), DUM(1), DUM(1),
  330. $ DUM(1), DUM(1), -1, IERR )
  331. LWORK_DGEBRD_NN = INT( DUM(1) )
  332. *
  333. CALL DGEQRF( M, N, DUM(1), M, DUM(1), DUM(1), -1, IERR )
  334. LWORK_DGEQRF_MN = INT( DUM(1) )
  335. *
  336. CALL DORGBR( 'Q', N, N, N, DUM(1), N, DUM(1), DUM(1), -1,
  337. $ IERR )
  338. LWORK_DORGBR_Q_NN = INT( DUM(1) )
  339. *
  340. CALL DORGQR( M, M, N, DUM(1), M, DUM(1), DUM(1), -1, IERR )
  341. LWORK_DORGQR_MM = INT( DUM(1) )
  342. *
  343. CALL DORGQR( M, N, N, DUM(1), M, DUM(1), DUM(1), -1, IERR )
  344. LWORK_DORGQR_MN = INT( DUM(1) )
  345. *
  346. CALL DORMBR( 'P', 'R', 'T', N, N, N, DUM(1), N,
  347. $ DUM(1), DUM(1), N, DUM(1), -1, IERR )
  348. LWORK_DORMBR_PRT_NN = INT( DUM(1) )
  349. *
  350. CALL DORMBR( 'Q', 'L', 'N', N, N, N, DUM(1), N,
  351. $ DUM(1), DUM(1), N, DUM(1), -1, IERR )
  352. LWORK_DORMBR_QLN_NN = INT( DUM(1) )
  353. *
  354. CALL DORMBR( 'Q', 'L', 'N', M, N, N, DUM(1), M,
  355. $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
  356. LWORK_DORMBR_QLN_MN = INT( DUM(1) )
  357. *
  358. CALL DORMBR( 'Q', 'L', 'N', M, M, N, DUM(1), M,
  359. $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
  360. LWORK_DORMBR_QLN_MM = INT( DUM(1) )
  361. *
  362. IF( M.GE.MNTHR ) THEN
  363. IF( WNTQN ) THEN
  364. *
  365. * Path 1 (M >> N, JOBZ='N')
  366. *
  367. WRKBL = N + LWORK_DGEQRF_MN
  368. WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD_NN )
  369. MAXWRK = MAX( WRKBL, BDSPAC + N )
  370. MINWRK = BDSPAC + N
  371. ELSE IF( WNTQO ) THEN
  372. *
  373. * Path 2 (M >> N, JOBZ='O')
  374. *
  375. WRKBL = N + LWORK_DGEQRF_MN
  376. WRKBL = MAX( WRKBL, N + LWORK_DORGQR_MN )
  377. WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD_NN )
  378. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_QLN_NN )
  379. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_PRT_NN )
  380. WRKBL = MAX( WRKBL, 3*N + BDSPAC )
  381. MAXWRK = WRKBL + 2*N*N
  382. MINWRK = BDSPAC + 2*N*N + 3*N
  383. ELSE IF( WNTQS ) THEN
  384. *
  385. * Path 3 (M >> N, JOBZ='S')
  386. *
  387. WRKBL = N + LWORK_DGEQRF_MN
  388. WRKBL = MAX( WRKBL, N + LWORK_DORGQR_MN )
  389. WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD_NN )
  390. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_QLN_NN )
  391. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_PRT_NN )
  392. WRKBL = MAX( WRKBL, 3*N + BDSPAC )
  393. MAXWRK = WRKBL + N*N
  394. MINWRK = BDSPAC + N*N + 3*N
  395. ELSE IF( WNTQA ) THEN
  396. *
  397. * Path 4 (M >> N, JOBZ='A')
  398. *
  399. WRKBL = N + LWORK_DGEQRF_MN
  400. WRKBL = MAX( WRKBL, N + LWORK_DORGQR_MM )
  401. WRKBL = MAX( WRKBL, 3*N + LWORK_DGEBRD_NN )
  402. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_QLN_NN )
  403. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_PRT_NN )
  404. WRKBL = MAX( WRKBL, 3*N + BDSPAC )
  405. MAXWRK = WRKBL + N*N
  406. MINWRK = N*N + MAX( 3*N + BDSPAC, N + M )
  407. END IF
  408. ELSE
  409. *
  410. * Path 5 (M >= N, but not much larger)
  411. *
  412. WRKBL = 3*N + LWORK_DGEBRD_MN
  413. IF( WNTQN ) THEN
  414. * Path 5n (M >= N, jobz='N')
  415. MAXWRK = MAX( WRKBL, 3*N + BDSPAC )
  416. MINWRK = 3*N + MAX( M, BDSPAC )
  417. ELSE IF( WNTQO ) THEN
  418. * Path 5o (M >= N, jobz='O')
  419. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_PRT_NN )
  420. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_QLN_MN )
  421. WRKBL = MAX( WRKBL, 3*N + BDSPAC )
  422. MAXWRK = WRKBL + M*N
  423. MINWRK = 3*N + MAX( M, N*N + BDSPAC )
  424. ELSE IF( WNTQS ) THEN
  425. * Path 5s (M >= N, jobz='S')
  426. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_QLN_MN )
  427. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_PRT_NN )
  428. MAXWRK = MAX( WRKBL, 3*N + BDSPAC )
  429. MINWRK = 3*N + MAX( M, BDSPAC )
  430. ELSE IF( WNTQA ) THEN
  431. * Path 5a (M >= N, jobz='A')
  432. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_QLN_MM )
  433. WRKBL = MAX( WRKBL, 3*N + LWORK_DORMBR_PRT_NN )
  434. MAXWRK = MAX( WRKBL, 3*N + BDSPAC )
  435. MINWRK = 3*N + MAX( M, BDSPAC )
  436. END IF
  437. END IF
  438. ELSE IF( MINMN.GT.0 ) THEN
  439. *
  440. * Compute space needed for DBDSDC
  441. *
  442. IF( WNTQN ) THEN
  443. * dbdsdc needs only 4*N (or 6*N for uplo=L for LAPACK <= 3.6)
  444. * keep 7*N for backwards compatibility.
  445. BDSPAC = 7*M
  446. ELSE
  447. BDSPAC = 3*M*M + 4*M
  448. END IF
  449. *
  450. * Compute space preferred for each routine
  451. CALL DGEBRD( M, N, DUM(1), M, DUM(1), DUM(1), DUM(1),
  452. $ DUM(1), DUM(1), -1, IERR )
  453. LWORK_DGEBRD_MN = INT( DUM(1) )
  454. *
  455. CALL DGEBRD( M, M, A, M, S, DUM(1), DUM(1),
  456. $ DUM(1), DUM(1), -1, IERR )
  457. LWORK_DGEBRD_MM = INT( DUM(1) )
  458. *
  459. CALL DGELQF( M, N, A, M, DUM(1), DUM(1), -1, IERR )
  460. LWORK_DGELQF_MN = INT( DUM(1) )
  461. *
  462. CALL DORGLQ( N, N, M, DUM(1), N, DUM(1), DUM(1), -1, IERR )
  463. LWORK_DORGLQ_NN = INT( DUM(1) )
  464. *
  465. CALL DORGLQ( M, N, M, A, M, DUM(1), DUM(1), -1, IERR )
  466. LWORK_DORGLQ_MN = INT( DUM(1) )
  467. *
  468. CALL DORGBR( 'P', M, M, M, A, N, DUM(1), DUM(1), -1, IERR )
  469. LWORK_DORGBR_P_MM = INT( DUM(1) )
  470. *
  471. CALL DORMBR( 'P', 'R', 'T', M, M, M, DUM(1), M,
  472. $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
  473. LWORK_DORMBR_PRT_MM = INT( DUM(1) )
  474. *
  475. CALL DORMBR( 'P', 'R', 'T', M, N, M, DUM(1), M,
  476. $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
  477. LWORK_DORMBR_PRT_MN = INT( DUM(1) )
  478. *
  479. CALL DORMBR( 'P', 'R', 'T', N, N, M, DUM(1), N,
  480. $ DUM(1), DUM(1), N, DUM(1), -1, IERR )
  481. LWORK_DORMBR_PRT_NN = INT( DUM(1) )
  482. *
  483. CALL DORMBR( 'Q', 'L', 'N', M, M, M, DUM(1), M,
  484. $ DUM(1), DUM(1), M, DUM(1), -1, IERR )
  485. LWORK_DORMBR_QLN_MM = INT( DUM(1) )
  486. *
  487. IF( N.GE.MNTHR ) THEN
  488. IF( WNTQN ) THEN
  489. *
  490. * Path 1t (N >> M, JOBZ='N')
  491. *
  492. WRKBL = M + LWORK_DGELQF_MN
  493. WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD_MM )
  494. MAXWRK = MAX( WRKBL, BDSPAC + M )
  495. MINWRK = BDSPAC + M
  496. ELSE IF( WNTQO ) THEN
  497. *
  498. * Path 2t (N >> M, JOBZ='O')
  499. *
  500. WRKBL = M + LWORK_DGELQF_MN
  501. WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_MN )
  502. WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD_MM )
  503. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_QLN_MM )
  504. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_PRT_MM )
  505. WRKBL = MAX( WRKBL, 3*M + BDSPAC )
  506. MAXWRK = WRKBL + 2*M*M
  507. MINWRK = BDSPAC + 2*M*M + 3*M
  508. ELSE IF( WNTQS ) THEN
  509. *
  510. * Path 3t (N >> M, JOBZ='S')
  511. *
  512. WRKBL = M + LWORK_DGELQF_MN
  513. WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_MN )
  514. WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD_MM )
  515. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_QLN_MM )
  516. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_PRT_MM )
  517. WRKBL = MAX( WRKBL, 3*M + BDSPAC )
  518. MAXWRK = WRKBL + M*M
  519. MINWRK = BDSPAC + M*M + 3*M
  520. ELSE IF( WNTQA ) THEN
  521. *
  522. * Path 4t (N >> M, JOBZ='A')
  523. *
  524. WRKBL = M + LWORK_DGELQF_MN
  525. WRKBL = MAX( WRKBL, M + LWORK_DORGLQ_NN )
  526. WRKBL = MAX( WRKBL, 3*M + LWORK_DGEBRD_MM )
  527. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_QLN_MM )
  528. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_PRT_MM )
  529. WRKBL = MAX( WRKBL, 3*M + BDSPAC )
  530. MAXWRK = WRKBL + M*M
  531. MINWRK = M*M + MAX( 3*M + BDSPAC, M + N )
  532. END IF
  533. ELSE
  534. *
  535. * Path 5t (N > M, but not much larger)
  536. *
  537. WRKBL = 3*M + LWORK_DGEBRD_MN
  538. IF( WNTQN ) THEN
  539. * Path 5tn (N > M, jobz='N')
  540. MAXWRK = MAX( WRKBL, 3*M + BDSPAC )
  541. MINWRK = 3*M + MAX( N, BDSPAC )
  542. ELSE IF( WNTQO ) THEN
  543. * Path 5to (N > M, jobz='O')
  544. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_QLN_MM )
  545. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_PRT_MN )
  546. WRKBL = MAX( WRKBL, 3*M + BDSPAC )
  547. MAXWRK = WRKBL + M*N
  548. MINWRK = 3*M + MAX( N, M*M + BDSPAC )
  549. ELSE IF( WNTQS ) THEN
  550. * Path 5ts (N > M, jobz='S')
  551. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_QLN_MM )
  552. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_PRT_MN )
  553. MAXWRK = MAX( WRKBL, 3*M + BDSPAC )
  554. MINWRK = 3*M + MAX( N, BDSPAC )
  555. ELSE IF( WNTQA ) THEN
  556. * Path 5ta (N > M, jobz='A')
  557. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_QLN_MM )
  558. WRKBL = MAX( WRKBL, 3*M + LWORK_DORMBR_PRT_NN )
  559. MAXWRK = MAX( WRKBL, 3*M + BDSPAC )
  560. MINWRK = 3*M + MAX( N, BDSPAC )
  561. END IF
  562. END IF
  563. END IF
  564. MAXWRK = MAX( MAXWRK, MINWRK )
  565. WORK( 1 ) = DROUNDUP_LWORK( MAXWRK )
  566. *
  567. IF( LWORK.LT.MINWRK .AND. .NOT.LQUERY ) THEN
  568. INFO = -12
  569. END IF
  570. END IF
  571. *
  572. IF( INFO.NE.0 ) THEN
  573. CALL XERBLA( 'DGESDD', -INFO )
  574. RETURN
  575. ELSE IF( LQUERY ) THEN
  576. RETURN
  577. END IF
  578. *
  579. * Quick return if possible
  580. *
  581. IF( M.EQ.0 .OR. N.EQ.0 ) THEN
  582. RETURN
  583. END IF
  584. *
  585. * Get machine constants
  586. *
  587. EPS = DLAMCH( 'P' )
  588. SMLNUM = SQRT( DLAMCH( 'S' ) ) / EPS
  589. BIGNUM = ONE / SMLNUM
  590. *
  591. * Scale A if max element outside range [SMLNUM,BIGNUM]
  592. *
  593. ANRM = DLANGE( 'M', M, N, A, LDA, DUM )
  594. IF( DISNAN( ANRM ) ) THEN
  595. INFO = -4
  596. RETURN
  597. END IF
  598. ISCL = 0
  599. IF( ANRM.GT.ZERO .AND. ANRM.LT.SMLNUM ) THEN
  600. ISCL = 1
  601. CALL DLASCL( 'G', 0, 0, ANRM, SMLNUM, M, N, A, LDA, IERR )
  602. ELSE IF( ANRM.GT.BIGNUM ) THEN
  603. ISCL = 1
  604. CALL DLASCL( 'G', 0, 0, ANRM, BIGNUM, M, N, A, LDA, IERR )
  605. END IF
  606. *
  607. IF( M.GE.N ) THEN
  608. *
  609. * A has at least as many rows as columns. If A has sufficiently
  610. * more rows than columns, first reduce using the QR
  611. * decomposition (if sufficient workspace available)
  612. *
  613. IF( M.GE.MNTHR ) THEN
  614. *
  615. IF( WNTQN ) THEN
  616. *
  617. * Path 1 (M >> N, JOBZ='N')
  618. * No singular vectors to be computed
  619. *
  620. ITAU = 1
  621. NWORK = ITAU + N
  622. *
  623. * Compute A=Q*R
  624. * Workspace: need N [tau] + N [work]
  625. * Workspace: prefer N [tau] + N*NB [work]
  626. *
  627. CALL DGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  628. $ LWORK - NWORK + 1, IERR )
  629. *
  630. * Zero out below R
  631. *
  632. CALL DLASET( 'L', N-1, N-1, ZERO, ZERO, A( 2, 1 ), LDA )
  633. IE = 1
  634. ITAUQ = IE + N
  635. ITAUP = ITAUQ + N
  636. NWORK = ITAUP + N
  637. *
  638. * Bidiagonalize R in A
  639. * Workspace: need 3*N [e, tauq, taup] + N [work]
  640. * Workspace: prefer 3*N [e, tauq, taup] + 2*N*NB [work]
  641. *
  642. CALL DGEBRD( N, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
  643. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  644. $ IERR )
  645. NWORK = IE + N
  646. *
  647. * Perform bidiagonal SVD, computing singular values only
  648. * Workspace: need N [e] + BDSPAC
  649. *
  650. CALL DBDSDC( 'U', 'N', N, S, WORK( IE ), DUM, 1, DUM, 1,
  651. $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
  652. *
  653. ELSE IF( WNTQO ) THEN
  654. *
  655. * Path 2 (M >> N, JOBZ = 'O')
  656. * N left singular vectors to be overwritten on A and
  657. * N right singular vectors to be computed in VT
  658. *
  659. IR = 1
  660. *
  661. * WORK(IR) is LDWRKR by N
  662. *
  663. IF( LWORK .GE. LDA*N + N*N + 3*N + BDSPAC ) THEN
  664. LDWRKR = LDA
  665. ELSE
  666. LDWRKR = ( LWORK - N*N - 3*N - BDSPAC ) / N
  667. END IF
  668. ITAU = IR + LDWRKR*N
  669. NWORK = ITAU + N
  670. *
  671. * Compute A=Q*R
  672. * Workspace: need N*N [R] + N [tau] + N [work]
  673. * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
  674. *
  675. CALL DGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  676. $ LWORK - NWORK + 1, IERR )
  677. *
  678. * Copy R to WORK(IR), zeroing out below it
  679. *
  680. CALL DLACPY( 'U', N, N, A, LDA, WORK( IR ), LDWRKR )
  681. CALL DLASET( 'L', N - 1, N - 1, ZERO, ZERO, WORK(IR+1),
  682. $ LDWRKR )
  683. *
  684. * Generate Q in A
  685. * Workspace: need N*N [R] + N [tau] + N [work]
  686. * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
  687. *
  688. CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
  689. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  690. IE = ITAU
  691. ITAUQ = IE + N
  692. ITAUP = ITAUQ + N
  693. NWORK = ITAUP + N
  694. *
  695. * Bidiagonalize R in WORK(IR)
  696. * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N [work]
  697. * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + 2*N*NB [work]
  698. *
  699. CALL DGEBRD( N, N, WORK( IR ), LDWRKR, S, WORK( IE ),
  700. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  701. $ LWORK - NWORK + 1, IERR )
  702. *
  703. * WORK(IU) is N by N
  704. *
  705. IU = NWORK
  706. NWORK = IU + N*N
  707. *
  708. * Perform bidiagonal SVD, computing left singular vectors
  709. * of bidiagonal matrix in WORK(IU) and computing right
  710. * singular vectors of bidiagonal matrix in VT
  711. * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N*N [U] + BDSPAC
  712. *
  713. CALL DBDSDC( 'U', 'I', N, S, WORK( IE ), WORK( IU ), N,
  714. $ VT, LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  715. $ INFO )
  716. *
  717. * Overwrite WORK(IU) by left singular vectors of R
  718. * and VT by right singular vectors of R
  719. * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N*N [U] + N [work]
  720. * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + N*N [U] + N*NB [work]
  721. *
  722. CALL DORMBR( 'Q', 'L', 'N', N, N, N, WORK( IR ), LDWRKR,
  723. $ WORK( ITAUQ ), WORK( IU ), N, WORK( NWORK ),
  724. $ LWORK - NWORK + 1, IERR )
  725. CALL DORMBR( 'P', 'R', 'T', N, N, N, WORK( IR ), LDWRKR,
  726. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  727. $ LWORK - NWORK + 1, IERR )
  728. *
  729. * Multiply Q in A by left singular vectors of R in
  730. * WORK(IU), storing result in WORK(IR) and copying to A
  731. * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N*N [U]
  732. * Workspace: prefer M*N [R] + 3*N [e, tauq, taup] + N*N [U]
  733. *
  734. DO 10 I = 1, M, LDWRKR
  735. CHUNK = MIN( M - I + 1, LDWRKR )
  736. CALL DGEMM( 'N', 'N', CHUNK, N, N, ONE, A( I, 1 ),
  737. $ LDA, WORK( IU ), N, ZERO, WORK( IR ),
  738. $ LDWRKR )
  739. CALL DLACPY( 'F', CHUNK, N, WORK( IR ), LDWRKR,
  740. $ A( I, 1 ), LDA )
  741. 10 CONTINUE
  742. *
  743. ELSE IF( WNTQS ) THEN
  744. *
  745. * Path 3 (M >> N, JOBZ='S')
  746. * N left singular vectors to be computed in U and
  747. * N right singular vectors to be computed in VT
  748. *
  749. IR = 1
  750. *
  751. * WORK(IR) is N by N
  752. *
  753. LDWRKR = N
  754. ITAU = IR + LDWRKR*N
  755. NWORK = ITAU + N
  756. *
  757. * Compute A=Q*R
  758. * Workspace: need N*N [R] + N [tau] + N [work]
  759. * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
  760. *
  761. CALL DGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  762. $ LWORK - NWORK + 1, IERR )
  763. *
  764. * Copy R to WORK(IR), zeroing out below it
  765. *
  766. CALL DLACPY( 'U', N, N, A, LDA, WORK( IR ), LDWRKR )
  767. CALL DLASET( 'L', N - 1, N - 1, ZERO, ZERO, WORK(IR+1),
  768. $ LDWRKR )
  769. *
  770. * Generate Q in A
  771. * Workspace: need N*N [R] + N [tau] + N [work]
  772. * Workspace: prefer N*N [R] + N [tau] + N*NB [work]
  773. *
  774. CALL DORGQR( M, N, N, A, LDA, WORK( ITAU ),
  775. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  776. IE = ITAU
  777. ITAUQ = IE + N
  778. ITAUP = ITAUQ + N
  779. NWORK = ITAUP + N
  780. *
  781. * Bidiagonalize R in WORK(IR)
  782. * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N [work]
  783. * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + 2*N*NB [work]
  784. *
  785. CALL DGEBRD( N, N, WORK( IR ), LDWRKR, S, WORK( IE ),
  786. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  787. $ LWORK - NWORK + 1, IERR )
  788. *
  789. * Perform bidiagonal SVD, computing left singular vectors
  790. * of bidiagoal matrix in U and computing right singular
  791. * vectors of bidiagonal matrix in VT
  792. * Workspace: need N*N [R] + 3*N [e, tauq, taup] + BDSPAC
  793. *
  794. CALL DBDSDC( 'U', 'I', N, S, WORK( IE ), U, LDU, VT,
  795. $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  796. $ INFO )
  797. *
  798. * Overwrite U by left singular vectors of R and VT
  799. * by right singular vectors of R
  800. * Workspace: need N*N [R] + 3*N [e, tauq, taup] + N [work]
  801. * Workspace: prefer N*N [R] + 3*N [e, tauq, taup] + N*NB [work]
  802. *
  803. CALL DORMBR( 'Q', 'L', 'N', N, N, N, WORK( IR ), LDWRKR,
  804. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  805. $ LWORK - NWORK + 1, IERR )
  806. *
  807. CALL DORMBR( 'P', 'R', 'T', N, N, N, WORK( IR ), LDWRKR,
  808. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  809. $ LWORK - NWORK + 1, IERR )
  810. *
  811. * Multiply Q in A by left singular vectors of R in
  812. * WORK(IR), storing result in U
  813. * Workspace: need N*N [R]
  814. *
  815. CALL DLACPY( 'F', N, N, U, LDU, WORK( IR ), LDWRKR )
  816. CALL DGEMM( 'N', 'N', M, N, N, ONE, A, LDA, WORK( IR ),
  817. $ LDWRKR, ZERO, U, LDU )
  818. *
  819. ELSE IF( WNTQA ) THEN
  820. *
  821. * Path 4 (M >> N, JOBZ='A')
  822. * M left singular vectors to be computed in U and
  823. * N right singular vectors to be computed in VT
  824. *
  825. IU = 1
  826. *
  827. * WORK(IU) is N by N
  828. *
  829. LDWRKU = N
  830. ITAU = IU + LDWRKU*N
  831. NWORK = ITAU + N
  832. *
  833. * Compute A=Q*R, copying result to U
  834. * Workspace: need N*N [U] + N [tau] + N [work]
  835. * Workspace: prefer N*N [U] + N [tau] + N*NB [work]
  836. *
  837. CALL DGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  838. $ LWORK - NWORK + 1, IERR )
  839. CALL DLACPY( 'L', M, N, A, LDA, U, LDU )
  840. *
  841. * Generate Q in U
  842. * Workspace: need N*N [U] + N [tau] + M [work]
  843. * Workspace: prefer N*N [U] + N [tau] + M*NB [work]
  844. CALL DORGQR( M, M, N, U, LDU, WORK( ITAU ),
  845. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  846. *
  847. * Produce R in A, zeroing out other entries
  848. *
  849. CALL DLASET( 'L', N-1, N-1, ZERO, ZERO, A( 2, 1 ), LDA )
  850. IE = ITAU
  851. ITAUQ = IE + N
  852. ITAUP = ITAUQ + N
  853. NWORK = ITAUP + N
  854. *
  855. * Bidiagonalize R in A
  856. * Workspace: need N*N [U] + 3*N [e, tauq, taup] + N [work]
  857. * Workspace: prefer N*N [U] + 3*N [e, tauq, taup] + 2*N*NB [work]
  858. *
  859. CALL DGEBRD( N, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
  860. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  861. $ IERR )
  862. *
  863. * Perform bidiagonal SVD, computing left singular vectors
  864. * of bidiagonal matrix in WORK(IU) and computing right
  865. * singular vectors of bidiagonal matrix in VT
  866. * Workspace: need N*N [U] + 3*N [e, tauq, taup] + BDSPAC
  867. *
  868. CALL DBDSDC( 'U', 'I', N, S, WORK( IE ), WORK( IU ), N,
  869. $ VT, LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  870. $ INFO )
  871. *
  872. * Overwrite WORK(IU) by left singular vectors of R and VT
  873. * by right singular vectors of R
  874. * Workspace: need N*N [U] + 3*N [e, tauq, taup] + N [work]
  875. * Workspace: prefer N*N [U] + 3*N [e, tauq, taup] + N*NB [work]
  876. *
  877. CALL DORMBR( 'Q', 'L', 'N', N, N, N, A, LDA,
  878. $ WORK( ITAUQ ), WORK( IU ), LDWRKU,
  879. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  880. CALL DORMBR( 'P', 'R', 'T', N, N, N, A, LDA,
  881. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  882. $ LWORK - NWORK + 1, IERR )
  883. *
  884. * Multiply Q in U by left singular vectors of R in
  885. * WORK(IU), storing result in A
  886. * Workspace: need N*N [U]
  887. *
  888. CALL DGEMM( 'N', 'N', M, N, N, ONE, U, LDU, WORK( IU ),
  889. $ LDWRKU, ZERO, A, LDA )
  890. *
  891. * Copy left singular vectors of A from A to U
  892. *
  893. CALL DLACPY( 'F', M, N, A, LDA, U, LDU )
  894. *
  895. END IF
  896. *
  897. ELSE
  898. *
  899. * M .LT. MNTHR
  900. *
  901. * Path 5 (M >= N, but not much larger)
  902. * Reduce to bidiagonal form without QR decomposition
  903. *
  904. IE = 1
  905. ITAUQ = IE + N
  906. ITAUP = ITAUQ + N
  907. NWORK = ITAUP + N
  908. *
  909. * Bidiagonalize A
  910. * Workspace: need 3*N [e, tauq, taup] + M [work]
  911. * Workspace: prefer 3*N [e, tauq, taup] + (M+N)*NB [work]
  912. *
  913. CALL DGEBRD( M, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
  914. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  915. $ IERR )
  916. IF( WNTQN ) THEN
  917. *
  918. * Path 5n (M >= N, JOBZ='N')
  919. * Perform bidiagonal SVD, only computing singular values
  920. * Workspace: need 3*N [e, tauq, taup] + BDSPAC
  921. *
  922. CALL DBDSDC( 'U', 'N', N, S, WORK( IE ), DUM, 1, DUM, 1,
  923. $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
  924. ELSE IF( WNTQO ) THEN
  925. * Path 5o (M >= N, JOBZ='O')
  926. IU = NWORK
  927. IF( LWORK .GE. M*N + 3*N + BDSPAC ) THEN
  928. *
  929. * WORK( IU ) is M by N
  930. *
  931. LDWRKU = M
  932. NWORK = IU + LDWRKU*N
  933. CALL DLASET( 'F', M, N, ZERO, ZERO, WORK( IU ),
  934. $ LDWRKU )
  935. * IR is unused; silence compile warnings
  936. IR = -1
  937. ELSE
  938. *
  939. * WORK( IU ) is N by N
  940. *
  941. LDWRKU = N
  942. NWORK = IU + LDWRKU*N
  943. *
  944. * WORK(IR) is LDWRKR by N
  945. *
  946. IR = NWORK
  947. LDWRKR = ( LWORK - N*N - 3*N ) / N
  948. END IF
  949. NWORK = IU + LDWRKU*N
  950. *
  951. * Perform bidiagonal SVD, computing left singular vectors
  952. * of bidiagonal matrix in WORK(IU) and computing right
  953. * singular vectors of bidiagonal matrix in VT
  954. * Workspace: need 3*N [e, tauq, taup] + N*N [U] + BDSPAC
  955. *
  956. CALL DBDSDC( 'U', 'I', N, S, WORK( IE ), WORK( IU ),
  957. $ LDWRKU, VT, LDVT, DUM, IDUM, WORK( NWORK ),
  958. $ IWORK, INFO )
  959. *
  960. * Overwrite VT by right singular vectors of A
  961. * Workspace: need 3*N [e, tauq, taup] + N*N [U] + N [work]
  962. * Workspace: prefer 3*N [e, tauq, taup] + N*N [U] + N*NB [work]
  963. *
  964. CALL DORMBR( 'P', 'R', 'T', N, N, N, A, LDA,
  965. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  966. $ LWORK - NWORK + 1, IERR )
  967. *
  968. IF( LWORK .GE. M*N + 3*N + BDSPAC ) THEN
  969. *
  970. * Path 5o-fast
  971. * Overwrite WORK(IU) by left singular vectors of A
  972. * Workspace: need 3*N [e, tauq, taup] + M*N [U] + N [work]
  973. * Workspace: prefer 3*N [e, tauq, taup] + M*N [U] + N*NB [work]
  974. *
  975. CALL DORMBR( 'Q', 'L', 'N', M, N, N, A, LDA,
  976. $ WORK( ITAUQ ), WORK( IU ), LDWRKU,
  977. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  978. *
  979. * Copy left singular vectors of A from WORK(IU) to A
  980. *
  981. CALL DLACPY( 'F', M, N, WORK( IU ), LDWRKU, A, LDA )
  982. ELSE
  983. *
  984. * Path 5o-slow
  985. * Generate Q in A
  986. * Workspace: need 3*N [e, tauq, taup] + N*N [U] + N [work]
  987. * Workspace: prefer 3*N [e, tauq, taup] + N*N [U] + N*NB [work]
  988. *
  989. CALL DORGBR( 'Q', M, N, N, A, LDA, WORK( ITAUQ ),
  990. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  991. *
  992. * Multiply Q in A by left singular vectors of
  993. * bidiagonal matrix in WORK(IU), storing result in
  994. * WORK(IR) and copying to A
  995. * Workspace: need 3*N [e, tauq, taup] + N*N [U] + NB*N [R]
  996. * Workspace: prefer 3*N [e, tauq, taup] + N*N [U] + M*N [R]
  997. *
  998. DO 20 I = 1, M, LDWRKR
  999. CHUNK = MIN( M - I + 1, LDWRKR )
  1000. CALL DGEMM( 'N', 'N', CHUNK, N, N, ONE, A( I, 1 ),
  1001. $ LDA, WORK( IU ), LDWRKU, ZERO,
  1002. $ WORK( IR ), LDWRKR )
  1003. CALL DLACPY( 'F', CHUNK, N, WORK( IR ), LDWRKR,
  1004. $ A( I, 1 ), LDA )
  1005. 20 CONTINUE
  1006. END IF
  1007. *
  1008. ELSE IF( WNTQS ) THEN
  1009. *
  1010. * Path 5s (M >= N, JOBZ='S')
  1011. * Perform bidiagonal SVD, computing left singular vectors
  1012. * of bidiagonal matrix in U and computing right singular
  1013. * vectors of bidiagonal matrix in VT
  1014. * Workspace: need 3*N [e, tauq, taup] + BDSPAC
  1015. *
  1016. CALL DLASET( 'F', M, N, ZERO, ZERO, U, LDU )
  1017. CALL DBDSDC( 'U', 'I', N, S, WORK( IE ), U, LDU, VT,
  1018. $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  1019. $ INFO )
  1020. *
  1021. * Overwrite U by left singular vectors of A and VT
  1022. * by right singular vectors of A
  1023. * Workspace: need 3*N [e, tauq, taup] + N [work]
  1024. * Workspace: prefer 3*N [e, tauq, taup] + N*NB [work]
  1025. *
  1026. CALL DORMBR( 'Q', 'L', 'N', M, N, N, A, LDA,
  1027. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1028. $ LWORK - NWORK + 1, IERR )
  1029. CALL DORMBR( 'P', 'R', 'T', N, N, N, A, LDA,
  1030. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1031. $ LWORK - NWORK + 1, IERR )
  1032. ELSE IF( WNTQA ) THEN
  1033. *
  1034. * Path 5a (M >= N, JOBZ='A')
  1035. * Perform bidiagonal SVD, computing left singular vectors
  1036. * of bidiagonal matrix in U and computing right singular
  1037. * vectors of bidiagonal matrix in VT
  1038. * Workspace: need 3*N [e, tauq, taup] + BDSPAC
  1039. *
  1040. CALL DLASET( 'F', M, M, ZERO, ZERO, U, LDU )
  1041. CALL DBDSDC( 'U', 'I', N, S, WORK( IE ), U, LDU, VT,
  1042. $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  1043. $ INFO )
  1044. *
  1045. * Set the right corner of U to identity matrix
  1046. *
  1047. IF( M.GT.N ) THEN
  1048. CALL DLASET( 'F', M - N, M - N, ZERO, ONE, U(N+1,N+1),
  1049. $ LDU )
  1050. END IF
  1051. *
  1052. * Overwrite U by left singular vectors of A and VT
  1053. * by right singular vectors of A
  1054. * Workspace: need 3*N [e, tauq, taup] + M [work]
  1055. * Workspace: prefer 3*N [e, tauq, taup] + M*NB [work]
  1056. *
  1057. CALL DORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  1058. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1059. $ LWORK - NWORK + 1, IERR )
  1060. CALL DORMBR( 'P', 'R', 'T', N, N, M, A, LDA,
  1061. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1062. $ LWORK - NWORK + 1, IERR )
  1063. END IF
  1064. *
  1065. END IF
  1066. *
  1067. ELSE
  1068. *
  1069. * A has more columns than rows. If A has sufficiently more
  1070. * columns than rows, first reduce using the LQ decomposition (if
  1071. * sufficient workspace available)
  1072. *
  1073. IF( N.GE.MNTHR ) THEN
  1074. *
  1075. IF( WNTQN ) THEN
  1076. *
  1077. * Path 1t (N >> M, JOBZ='N')
  1078. * No singular vectors to be computed
  1079. *
  1080. ITAU = 1
  1081. NWORK = ITAU + M
  1082. *
  1083. * Compute A=L*Q
  1084. * Workspace: need M [tau] + M [work]
  1085. * Workspace: prefer M [tau] + M*NB [work]
  1086. *
  1087. CALL DGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1088. $ LWORK - NWORK + 1, IERR )
  1089. *
  1090. * Zero out above L
  1091. *
  1092. CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ), LDA )
  1093. IE = 1
  1094. ITAUQ = IE + M
  1095. ITAUP = ITAUQ + M
  1096. NWORK = ITAUP + M
  1097. *
  1098. * Bidiagonalize L in A
  1099. * Workspace: need 3*M [e, tauq, taup] + M [work]
  1100. * Workspace: prefer 3*M [e, tauq, taup] + 2*M*NB [work]
  1101. *
  1102. CALL DGEBRD( M, M, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
  1103. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1104. $ IERR )
  1105. NWORK = IE + M
  1106. *
  1107. * Perform bidiagonal SVD, computing singular values only
  1108. * Workspace: need M [e] + BDSPAC
  1109. *
  1110. CALL DBDSDC( 'U', 'N', M, S, WORK( IE ), DUM, 1, DUM, 1,
  1111. $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
  1112. *
  1113. ELSE IF( WNTQO ) THEN
  1114. *
  1115. * Path 2t (N >> M, JOBZ='O')
  1116. * M right singular vectors to be overwritten on A and
  1117. * M left singular vectors to be computed in U
  1118. *
  1119. IVT = 1
  1120. *
  1121. * WORK(IVT) is M by M
  1122. * WORK(IL) is M by M; it is later resized to M by chunk for gemm
  1123. *
  1124. IL = IVT + M*M
  1125. IF( LWORK .GE. M*N + M*M + 3*M + BDSPAC ) THEN
  1126. LDWRKL = M
  1127. CHUNK = N
  1128. ELSE
  1129. LDWRKL = M
  1130. CHUNK = ( LWORK - M*M ) / M
  1131. END IF
  1132. ITAU = IL + LDWRKL*M
  1133. NWORK = ITAU + M
  1134. *
  1135. * Compute A=L*Q
  1136. * Workspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
  1137. * Workspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
  1138. *
  1139. CALL DGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1140. $ LWORK - NWORK + 1, IERR )
  1141. *
  1142. * Copy L to WORK(IL), zeroing about above it
  1143. *
  1144. CALL DLACPY( 'L', M, M, A, LDA, WORK( IL ), LDWRKL )
  1145. CALL DLASET( 'U', M - 1, M - 1, ZERO, ZERO,
  1146. $ WORK( IL + LDWRKL ), LDWRKL )
  1147. *
  1148. * Generate Q in A
  1149. * Workspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
  1150. * Workspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
  1151. *
  1152. CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
  1153. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  1154. IE = ITAU
  1155. ITAUQ = IE + M
  1156. ITAUP = ITAUQ + M
  1157. NWORK = ITAUP + M
  1158. *
  1159. * Bidiagonalize L in WORK(IL)
  1160. * Workspace: need M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + M [work]
  1161. * Workspace: prefer M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + 2*M*NB [work]
  1162. *
  1163. CALL DGEBRD( M, M, WORK( IL ), LDWRKL, S, WORK( IE ),
  1164. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  1165. $ LWORK - NWORK + 1, IERR )
  1166. *
  1167. * Perform bidiagonal SVD, computing left singular vectors
  1168. * of bidiagonal matrix in U, and computing right singular
  1169. * vectors of bidiagonal matrix in WORK(IVT)
  1170. * Workspace: need M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + BDSPAC
  1171. *
  1172. CALL DBDSDC( 'U', 'I', M, S, WORK( IE ), U, LDU,
  1173. $ WORK( IVT ), M, DUM, IDUM, WORK( NWORK ),
  1174. $ IWORK, INFO )
  1175. *
  1176. * Overwrite U by left singular vectors of L and WORK(IVT)
  1177. * by right singular vectors of L
  1178. * Workspace: need M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + M [work]
  1179. * Workspace: prefer M*M [VT] + M*M [L] + 3*M [e, tauq, taup] + M*NB [work]
  1180. *
  1181. CALL DORMBR( 'Q', 'L', 'N', M, M, M, WORK( IL ), LDWRKL,
  1182. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1183. $ LWORK - NWORK + 1, IERR )
  1184. CALL DORMBR( 'P', 'R', 'T', M, M, M, WORK( IL ), LDWRKL,
  1185. $ WORK( ITAUP ), WORK( IVT ), M,
  1186. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  1187. *
  1188. * Multiply right singular vectors of L in WORK(IVT) by Q
  1189. * in A, storing result in WORK(IL) and copying to A
  1190. * Workspace: need M*M [VT] + M*M [L]
  1191. * Workspace: prefer M*M [VT] + M*N [L]
  1192. * At this point, L is resized as M by chunk.
  1193. *
  1194. DO 30 I = 1, N, CHUNK
  1195. BLK = MIN( N - I + 1, CHUNK )
  1196. CALL DGEMM( 'N', 'N', M, BLK, M, ONE, WORK( IVT ), M,
  1197. $ A( 1, I ), LDA, ZERO, WORK( IL ), LDWRKL )
  1198. CALL DLACPY( 'F', M, BLK, WORK( IL ), LDWRKL,
  1199. $ A( 1, I ), LDA )
  1200. 30 CONTINUE
  1201. *
  1202. ELSE IF( WNTQS ) THEN
  1203. *
  1204. * Path 3t (N >> M, JOBZ='S')
  1205. * M right singular vectors to be computed in VT and
  1206. * M left singular vectors to be computed in U
  1207. *
  1208. IL = 1
  1209. *
  1210. * WORK(IL) is M by M
  1211. *
  1212. LDWRKL = M
  1213. ITAU = IL + LDWRKL*M
  1214. NWORK = ITAU + M
  1215. *
  1216. * Compute A=L*Q
  1217. * Workspace: need M*M [L] + M [tau] + M [work]
  1218. * Workspace: prefer M*M [L] + M [tau] + M*NB [work]
  1219. *
  1220. CALL DGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1221. $ LWORK - NWORK + 1, IERR )
  1222. *
  1223. * Copy L to WORK(IL), zeroing out above it
  1224. *
  1225. CALL DLACPY( 'L', M, M, A, LDA, WORK( IL ), LDWRKL )
  1226. CALL DLASET( 'U', M - 1, M - 1, ZERO, ZERO,
  1227. $ WORK( IL + LDWRKL ), LDWRKL )
  1228. *
  1229. * Generate Q in A
  1230. * Workspace: need M*M [L] + M [tau] + M [work]
  1231. * Workspace: prefer M*M [L] + M [tau] + M*NB [work]
  1232. *
  1233. CALL DORGLQ( M, N, M, A, LDA, WORK( ITAU ),
  1234. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  1235. IE = ITAU
  1236. ITAUQ = IE + M
  1237. ITAUP = ITAUQ + M
  1238. NWORK = ITAUP + M
  1239. *
  1240. * Bidiagonalize L in WORK(IU).
  1241. * Workspace: need M*M [L] + 3*M [e, tauq, taup] + M [work]
  1242. * Workspace: prefer M*M [L] + 3*M [e, tauq, taup] + 2*M*NB [work]
  1243. *
  1244. CALL DGEBRD( M, M, WORK( IL ), LDWRKL, S, WORK( IE ),
  1245. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  1246. $ LWORK - NWORK + 1, IERR )
  1247. *
  1248. * Perform bidiagonal SVD, computing left singular vectors
  1249. * of bidiagonal matrix in U and computing right singular
  1250. * vectors of bidiagonal matrix in VT
  1251. * Workspace: need M*M [L] + 3*M [e, tauq, taup] + BDSPAC
  1252. *
  1253. CALL DBDSDC( 'U', 'I', M, S, WORK( IE ), U, LDU, VT,
  1254. $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  1255. $ INFO )
  1256. *
  1257. * Overwrite U by left singular vectors of L and VT
  1258. * by right singular vectors of L
  1259. * Workspace: need M*M [L] + 3*M [e, tauq, taup] + M [work]
  1260. * Workspace: prefer M*M [L] + 3*M [e, tauq, taup] + M*NB [work]
  1261. *
  1262. CALL DORMBR( 'Q', 'L', 'N', M, M, M, WORK( IL ), LDWRKL,
  1263. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1264. $ LWORK - NWORK + 1, IERR )
  1265. CALL DORMBR( 'P', 'R', 'T', M, M, M, WORK( IL ), LDWRKL,
  1266. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1267. $ LWORK - NWORK + 1, IERR )
  1268. *
  1269. * Multiply right singular vectors of L in WORK(IL) by
  1270. * Q in A, storing result in VT
  1271. * Workspace: need M*M [L]
  1272. *
  1273. CALL DLACPY( 'F', M, M, VT, LDVT, WORK( IL ), LDWRKL )
  1274. CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IL ), LDWRKL,
  1275. $ A, LDA, ZERO, VT, LDVT )
  1276. *
  1277. ELSE IF( WNTQA ) THEN
  1278. *
  1279. * Path 4t (N >> M, JOBZ='A')
  1280. * N right singular vectors to be computed in VT and
  1281. * M left singular vectors to be computed in U
  1282. *
  1283. IVT = 1
  1284. *
  1285. * WORK(IVT) is M by M
  1286. *
  1287. LDWKVT = M
  1288. ITAU = IVT + LDWKVT*M
  1289. NWORK = ITAU + M
  1290. *
  1291. * Compute A=L*Q, copying result to VT
  1292. * Workspace: need M*M [VT] + M [tau] + M [work]
  1293. * Workspace: prefer M*M [VT] + M [tau] + M*NB [work]
  1294. *
  1295. CALL DGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1296. $ LWORK - NWORK + 1, IERR )
  1297. CALL DLACPY( 'U', M, N, A, LDA, VT, LDVT )
  1298. *
  1299. * Generate Q in VT
  1300. * Workspace: need M*M [VT] + M [tau] + N [work]
  1301. * Workspace: prefer M*M [VT] + M [tau] + N*NB [work]
  1302. *
  1303. CALL DORGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
  1304. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  1305. *
  1306. * Produce L in A, zeroing out other entries
  1307. *
  1308. CALL DLASET( 'U', M-1, M-1, ZERO, ZERO, A( 1, 2 ), LDA )
  1309. IE = ITAU
  1310. ITAUQ = IE + M
  1311. ITAUP = ITAUQ + M
  1312. NWORK = ITAUP + M
  1313. *
  1314. * Bidiagonalize L in A
  1315. * Workspace: need M*M [VT] + 3*M [e, tauq, taup] + M [work]
  1316. * Workspace: prefer M*M [VT] + 3*M [e, tauq, taup] + 2*M*NB [work]
  1317. *
  1318. CALL DGEBRD( M, M, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
  1319. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1320. $ IERR )
  1321. *
  1322. * Perform bidiagonal SVD, computing left singular vectors
  1323. * of bidiagonal matrix in U and computing right singular
  1324. * vectors of bidiagonal matrix in WORK(IVT)
  1325. * Workspace: need M*M [VT] + 3*M [e, tauq, taup] + BDSPAC
  1326. *
  1327. CALL DBDSDC( 'U', 'I', M, S, WORK( IE ), U, LDU,
  1328. $ WORK( IVT ), LDWKVT, DUM, IDUM,
  1329. $ WORK( NWORK ), IWORK, INFO )
  1330. *
  1331. * Overwrite U by left singular vectors of L and WORK(IVT)
  1332. * by right singular vectors of L
  1333. * Workspace: need M*M [VT] + 3*M [e, tauq, taup]+ M [work]
  1334. * Workspace: prefer M*M [VT] + 3*M [e, tauq, taup]+ M*NB [work]
  1335. *
  1336. CALL DORMBR( 'Q', 'L', 'N', M, M, M, A, LDA,
  1337. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1338. $ LWORK - NWORK + 1, IERR )
  1339. CALL DORMBR( 'P', 'R', 'T', M, M, M, A, LDA,
  1340. $ WORK( ITAUP ), WORK( IVT ), LDWKVT,
  1341. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  1342. *
  1343. * Multiply right singular vectors of L in WORK(IVT) by
  1344. * Q in VT, storing result in A
  1345. * Workspace: need M*M [VT]
  1346. *
  1347. CALL DGEMM( 'N', 'N', M, N, M, ONE, WORK( IVT ), LDWKVT,
  1348. $ VT, LDVT, ZERO, A, LDA )
  1349. *
  1350. * Copy right singular vectors of A from A to VT
  1351. *
  1352. CALL DLACPY( 'F', M, N, A, LDA, VT, LDVT )
  1353. *
  1354. END IF
  1355. *
  1356. ELSE
  1357. *
  1358. * N .LT. MNTHR
  1359. *
  1360. * Path 5t (N > M, but not much larger)
  1361. * Reduce to bidiagonal form without LQ decomposition
  1362. *
  1363. IE = 1
  1364. ITAUQ = IE + M
  1365. ITAUP = ITAUQ + M
  1366. NWORK = ITAUP + M
  1367. *
  1368. * Bidiagonalize A
  1369. * Workspace: need 3*M [e, tauq, taup] + N [work]
  1370. * Workspace: prefer 3*M [e, tauq, taup] + (M+N)*NB [work]
  1371. *
  1372. CALL DGEBRD( M, N, A, LDA, S, WORK( IE ), WORK( ITAUQ ),
  1373. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1374. $ IERR )
  1375. IF( WNTQN ) THEN
  1376. *
  1377. * Path 5tn (N > M, JOBZ='N')
  1378. * Perform bidiagonal SVD, only computing singular values
  1379. * Workspace: need 3*M [e, tauq, taup] + BDSPAC
  1380. *
  1381. CALL DBDSDC( 'L', 'N', M, S, WORK( IE ), DUM, 1, DUM, 1,
  1382. $ DUM, IDUM, WORK( NWORK ), IWORK, INFO )
  1383. ELSE IF( WNTQO ) THEN
  1384. * Path 5to (N > M, JOBZ='O')
  1385. LDWKVT = M
  1386. IVT = NWORK
  1387. IF( LWORK .GE. M*N + 3*M + BDSPAC ) THEN
  1388. *
  1389. * WORK( IVT ) is M by N
  1390. *
  1391. CALL DLASET( 'F', M, N, ZERO, ZERO, WORK( IVT ),
  1392. $ LDWKVT )
  1393. NWORK = IVT + LDWKVT*N
  1394. * IL is unused; silence compile warnings
  1395. IL = -1
  1396. ELSE
  1397. *
  1398. * WORK( IVT ) is M by M
  1399. *
  1400. NWORK = IVT + LDWKVT*M
  1401. IL = NWORK
  1402. *
  1403. * WORK(IL) is M by CHUNK
  1404. *
  1405. CHUNK = ( LWORK - M*M - 3*M ) / M
  1406. END IF
  1407. *
  1408. * Perform bidiagonal SVD, computing left singular vectors
  1409. * of bidiagonal matrix in U and computing right singular
  1410. * vectors of bidiagonal matrix in WORK(IVT)
  1411. * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + BDSPAC
  1412. *
  1413. CALL DBDSDC( 'L', 'I', M, S, WORK( IE ), U, LDU,
  1414. $ WORK( IVT ), LDWKVT, DUM, IDUM,
  1415. $ WORK( NWORK ), IWORK, INFO )
  1416. *
  1417. * Overwrite U by left singular vectors of A
  1418. * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + M [work]
  1419. * Workspace: prefer 3*M [e, tauq, taup] + M*M [VT] + M*NB [work]
  1420. *
  1421. CALL DORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  1422. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1423. $ LWORK - NWORK + 1, IERR )
  1424. *
  1425. IF( LWORK .GE. M*N + 3*M + BDSPAC ) THEN
  1426. *
  1427. * Path 5to-fast
  1428. * Overwrite WORK(IVT) by left singular vectors of A
  1429. * Workspace: need 3*M [e, tauq, taup] + M*N [VT] + M [work]
  1430. * Workspace: prefer 3*M [e, tauq, taup] + M*N [VT] + M*NB [work]
  1431. *
  1432. CALL DORMBR( 'P', 'R', 'T', M, N, M, A, LDA,
  1433. $ WORK( ITAUP ), WORK( IVT ), LDWKVT,
  1434. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  1435. *
  1436. * Copy right singular vectors of A from WORK(IVT) to A
  1437. *
  1438. CALL DLACPY( 'F', M, N, WORK( IVT ), LDWKVT, A, LDA )
  1439. ELSE
  1440. *
  1441. * Path 5to-slow
  1442. * Generate P**T in A
  1443. * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + M [work]
  1444. * Workspace: prefer 3*M [e, tauq, taup] + M*M [VT] + M*NB [work]
  1445. *
  1446. CALL DORGBR( 'P', M, N, M, A, LDA, WORK( ITAUP ),
  1447. $ WORK( NWORK ), LWORK - NWORK + 1, IERR )
  1448. *
  1449. * Multiply Q in A by right singular vectors of
  1450. * bidiagonal matrix in WORK(IVT), storing result in
  1451. * WORK(IL) and copying to A
  1452. * Workspace: need 3*M [e, tauq, taup] + M*M [VT] + M*NB [L]
  1453. * Workspace: prefer 3*M [e, tauq, taup] + M*M [VT] + M*N [L]
  1454. *
  1455. DO 40 I = 1, N, CHUNK
  1456. BLK = MIN( N - I + 1, CHUNK )
  1457. CALL DGEMM( 'N', 'N', M, BLK, M, ONE, WORK( IVT ),
  1458. $ LDWKVT, A( 1, I ), LDA, ZERO,
  1459. $ WORK( IL ), M )
  1460. CALL DLACPY( 'F', M, BLK, WORK( IL ), M, A( 1, I ),
  1461. $ LDA )
  1462. 40 CONTINUE
  1463. END IF
  1464. ELSE IF( WNTQS ) THEN
  1465. *
  1466. * Path 5ts (N > M, JOBZ='S')
  1467. * Perform bidiagonal SVD, computing left singular vectors
  1468. * of bidiagonal matrix in U and computing right singular
  1469. * vectors of bidiagonal matrix in VT
  1470. * Workspace: need 3*M [e, tauq, taup] + BDSPAC
  1471. *
  1472. CALL DLASET( 'F', M, N, ZERO, ZERO, VT, LDVT )
  1473. CALL DBDSDC( 'L', 'I', M, S, WORK( IE ), U, LDU, VT,
  1474. $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  1475. $ INFO )
  1476. *
  1477. * Overwrite U by left singular vectors of A and VT
  1478. * by right singular vectors of A
  1479. * Workspace: need 3*M [e, tauq, taup] + M [work]
  1480. * Workspace: prefer 3*M [e, tauq, taup] + M*NB [work]
  1481. *
  1482. CALL DORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  1483. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1484. $ LWORK - NWORK + 1, IERR )
  1485. CALL DORMBR( 'P', 'R', 'T', M, N, M, A, LDA,
  1486. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1487. $ LWORK - NWORK + 1, IERR )
  1488. ELSE IF( WNTQA ) THEN
  1489. *
  1490. * Path 5ta (N > M, JOBZ='A')
  1491. * Perform bidiagonal SVD, computing left singular vectors
  1492. * of bidiagonal matrix in U and computing right singular
  1493. * vectors of bidiagonal matrix in VT
  1494. * Workspace: need 3*M [e, tauq, taup] + BDSPAC
  1495. *
  1496. CALL DLASET( 'F', N, N, ZERO, ZERO, VT, LDVT )
  1497. CALL DBDSDC( 'L', 'I', M, S, WORK( IE ), U, LDU, VT,
  1498. $ LDVT, DUM, IDUM, WORK( NWORK ), IWORK,
  1499. $ INFO )
  1500. *
  1501. * Set the right corner of VT to identity matrix
  1502. *
  1503. IF( N.GT.M ) THEN
  1504. CALL DLASET( 'F', N-M, N-M, ZERO, ONE, VT(M+1,M+1),
  1505. $ LDVT )
  1506. END IF
  1507. *
  1508. * Overwrite U by left singular vectors of A and VT
  1509. * by right singular vectors of A
  1510. * Workspace: need 3*M [e, tauq, taup] + N [work]
  1511. * Workspace: prefer 3*M [e, tauq, taup] + N*NB [work]
  1512. *
  1513. CALL DORMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  1514. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1515. $ LWORK - NWORK + 1, IERR )
  1516. CALL DORMBR( 'P', 'R', 'T', N, N, M, A, LDA,
  1517. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1518. $ LWORK - NWORK + 1, IERR )
  1519. END IF
  1520. *
  1521. END IF
  1522. *
  1523. END IF
  1524. *
  1525. * Undo scaling if necessary
  1526. *
  1527. IF( ISCL.EQ.1 ) THEN
  1528. IF( ANRM.GT.BIGNUM )
  1529. $ CALL DLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN, 1, S, MINMN,
  1530. $ IERR )
  1531. IF( ANRM.LT.SMLNUM )
  1532. $ CALL DLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN, 1, S, MINMN,
  1533. $ IERR )
  1534. END IF
  1535. *
  1536. * Return optimal workspace in WORK(1)
  1537. *
  1538. WORK( 1 ) = DROUNDUP_LWORK( MAXWRK )
  1539. *
  1540. RETURN
  1541. *
  1542. * End of DGESDD
  1543. *
  1544. END