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.

cgesdd.f 90 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217
  1. *> \brief \b CGESDD
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download CGESDD + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgesdd.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgesdd.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgesdd.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE CGESDD( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
  22. * WORK, LWORK, RWORK, 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. * REAL RWORK( * ), S( * )
  31. * COMPLEX A( LDA, * ), U( LDU, * ), VT( LDVT, * ),
  32. * $ WORK( * )
  33. * ..
  34. *
  35. *
  36. *> \par Purpose:
  37. * =============
  38. *>
  39. *> \verbatim
  40. *>
  41. *> CGESDD computes the singular value decomposition (SVD) of a complex
  42. *> M-by-N matrix A, optionally computing the left and/or right singular
  43. *> vectors, by using divide-and-conquer method. The SVD is written
  44. *>
  45. *> A = U * SIGMA * conjugate-transpose(V)
  46. *>
  47. *> where SIGMA is an M-by-N matrix which is zero except for its
  48. *> min(m,n) diagonal elements, U is an M-by-M unitary matrix, and
  49. *> V is an N-by-N unitary matrix. The diagonal elements of SIGMA
  50. *> are the singular values of A; they are real and non-negative, and
  51. *> are returned in descending order. The first min(m,n) columns of
  52. *> U and V are the left and right singular vectors of A.
  53. *>
  54. *> Note that the routine returns VT = V**H, not V.
  55. *>
  56. *> \endverbatim
  57. *
  58. * Arguments:
  59. * ==========
  60. *
  61. *> \param[in] JOBZ
  62. *> \verbatim
  63. *> JOBZ is CHARACTER*1
  64. *> Specifies options for computing all or part of the matrix U:
  65. *> = 'A': all M columns of U and all N rows of V**H are
  66. *> returned in the arrays U and VT;
  67. *> = 'S': the first min(M,N) columns of U and the first
  68. *> min(M,N) rows of V**H are returned in the arrays U
  69. *> and VT;
  70. *> = 'O': If M >= N, the first N columns of U are overwritten
  71. *> in the array A and all rows of V**H are returned in
  72. *> the array VT;
  73. *> otherwise, all columns of U are returned in the
  74. *> array U and the first M rows of V**H are overwritten
  75. *> in the array A;
  76. *> = 'N': no columns of U or rows of V**H are computed.
  77. *> \endverbatim
  78. *>
  79. *> \param[in] M
  80. *> \verbatim
  81. *> M is INTEGER
  82. *> The number of rows of the input matrix A. M >= 0.
  83. *> \endverbatim
  84. *>
  85. *> \param[in] N
  86. *> \verbatim
  87. *> N is INTEGER
  88. *> The number of columns of the input matrix A. N >= 0.
  89. *> \endverbatim
  90. *>
  91. *> \param[in,out] A
  92. *> \verbatim
  93. *> A is COMPLEX array, dimension (LDA,N)
  94. *> On entry, the M-by-N matrix A.
  95. *> On exit,
  96. *> if JOBZ = 'O', A is overwritten with the first N columns
  97. *> of U (the left singular vectors, stored
  98. *> columnwise) if M >= N;
  99. *> A is overwritten with the first M rows
  100. *> of V**H (the right singular vectors, stored
  101. *> rowwise) otherwise.
  102. *> if JOBZ .ne. 'O', the contents of A are destroyed.
  103. *> \endverbatim
  104. *>
  105. *> \param[in] LDA
  106. *> \verbatim
  107. *> LDA is INTEGER
  108. *> The leading dimension of the array A. LDA >= max(1,M).
  109. *> \endverbatim
  110. *>
  111. *> \param[out] S
  112. *> \verbatim
  113. *> S is REAL array, dimension (min(M,N))
  114. *> The singular values of A, sorted so that S(i) >= S(i+1).
  115. *> \endverbatim
  116. *>
  117. *> \param[out] U
  118. *> \verbatim
  119. *> U is COMPLEX array, dimension (LDU,UCOL)
  120. *> UCOL = M if JOBZ = 'A' or JOBZ = 'O' and M < N;
  121. *> UCOL = min(M,N) if JOBZ = 'S'.
  122. *> If JOBZ = 'A' or JOBZ = 'O' and M < N, U contains the M-by-M
  123. *> unitary matrix U;
  124. *> if JOBZ = 'S', U contains the first min(M,N) columns of U
  125. *> (the left singular vectors, stored columnwise);
  126. *> if JOBZ = 'O' and M >= N, or JOBZ = 'N', U is not referenced.
  127. *> \endverbatim
  128. *>
  129. *> \param[in] LDU
  130. *> \verbatim
  131. *> LDU is INTEGER
  132. *> The leading dimension of the array U. LDU >= 1;
  133. *> if JOBZ = 'S' or 'A' or JOBZ = 'O' and M < N, LDU >= M.
  134. *> \endverbatim
  135. *>
  136. *> \param[out] VT
  137. *> \verbatim
  138. *> VT is COMPLEX array, dimension (LDVT,N)
  139. *> If JOBZ = 'A' or JOBZ = 'O' and M >= N, VT contains the
  140. *> N-by-N unitary matrix V**H;
  141. *> if JOBZ = 'S', VT contains the first min(M,N) rows of
  142. *> V**H (the right singular vectors, stored rowwise);
  143. *> if JOBZ = 'O' and M < N, or JOBZ = 'N', VT is not referenced.
  144. *> \endverbatim
  145. *>
  146. *> \param[in] LDVT
  147. *> \verbatim
  148. *> LDVT is INTEGER
  149. *> The leading dimension of the array VT. LDVT >= 1;
  150. *> if JOBZ = 'A' or JOBZ = 'O' and M >= N, LDVT >= N;
  151. *> if JOBZ = 'S', LDVT >= min(M,N).
  152. *> \endverbatim
  153. *>
  154. *> \param[out] WORK
  155. *> \verbatim
  156. *> WORK is COMPLEX array, dimension (MAX(1,LWORK))
  157. *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
  158. *> \endverbatim
  159. *>
  160. *> \param[in] LWORK
  161. *> \verbatim
  162. *> LWORK is INTEGER
  163. *> The dimension of the array WORK. LWORK >= 1.
  164. *> If LWORK = -1, a workspace query is assumed. The optimal
  165. *> size for the WORK array is calculated and stored in WORK(1),
  166. *> and no other work except argument checking is performed.
  167. *>
  168. *> Let mx = max(M,N) and mn = min(M,N).
  169. *> If JOBZ = 'N', LWORK >= 2*mn + mx.
  170. *> If JOBZ = 'O', LWORK >= 2*mn*mn + 2*mn + mx.
  171. *> If JOBZ = 'S', LWORK >= mn*mn + 3*mn.
  172. *> If JOBZ = 'A', LWORK >= mn*mn + 2*mn + mx.
  173. *> These are not tight minimums in all cases; see comments inside code.
  174. *> For good performance, LWORK should generally be larger;
  175. *> a query is recommended.
  176. *> \endverbatim
  177. *>
  178. *> \param[out] RWORK
  179. *> \verbatim
  180. *> RWORK is REAL array, dimension (MAX(1,LRWORK))
  181. *> Let mx = max(M,N) and mn = min(M,N).
  182. *> If JOBZ = 'N', LRWORK >= 5*mn (LAPACK <= 3.6 needs 7*mn);
  183. *> else if mx >> mn, LRWORK >= 5*mn*mn + 5*mn;
  184. *> else LRWORK >= max( 5*mn*mn + 5*mn,
  185. *> 2*mx*mn + 2*mn*mn + mn ).
  186. *> \endverbatim
  187. *>
  188. *> \param[out] IWORK
  189. *> \verbatim
  190. *> IWORK is INTEGER array, dimension (8*min(M,N))
  191. *> \endverbatim
  192. *>
  193. *> \param[out] INFO
  194. *> \verbatim
  195. *> INFO is INTEGER
  196. *> < 0: if INFO = -i, the i-th argument had an illegal value.
  197. *> = -4: if A had a NAN entry.
  198. *> > 0: The updating process of SBDSDC did not converge.
  199. *> = 0: successful exit.
  200. *> \endverbatim
  201. *
  202. * Authors:
  203. * ========
  204. *
  205. *> \author Univ. of Tennessee
  206. *> \author Univ. of California Berkeley
  207. *> \author Univ. of Colorado Denver
  208. *> \author NAG Ltd.
  209. *
  210. *> \ingroup complexGEsing
  211. *
  212. *> \par Contributors:
  213. * ==================
  214. *>
  215. *> Ming Gu and Huan Ren, Computer Science Division, University of
  216. *> California at Berkeley, USA
  217. *>
  218. * =====================================================================
  219. SUBROUTINE CGESDD( JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT,
  220. $ WORK, LWORK, RWORK, IWORK, INFO )
  221. implicit none
  222. *
  223. * -- LAPACK driver routine --
  224. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  225. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  226. *
  227. * .. Scalar Arguments ..
  228. CHARACTER JOBZ
  229. INTEGER INFO, LDA, LDU, LDVT, LWORK, M, N
  230. * ..
  231. * .. Array Arguments ..
  232. INTEGER IWORK( * )
  233. REAL RWORK( * ), S( * )
  234. COMPLEX A( LDA, * ), U( LDU, * ), VT( LDVT, * ),
  235. $ WORK( * )
  236. * ..
  237. *
  238. * =====================================================================
  239. *
  240. * .. Parameters ..
  241. COMPLEX CZERO, CONE
  242. PARAMETER ( CZERO = ( 0.0E+0, 0.0E+0 ),
  243. $ CONE = ( 1.0E+0, 0.0E+0 ) )
  244. REAL ZERO, ONE
  245. PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 )
  246. * ..
  247. * .. Local Scalars ..
  248. LOGICAL LQUERY, WNTQA, WNTQAS, WNTQN, WNTQO, WNTQS
  249. INTEGER BLK, CHUNK, I, IE, IERR, IL, IR, IRU, IRVT,
  250. $ ISCL, ITAU, ITAUP, ITAUQ, IU, IVT, LDWKVT,
  251. $ LDWRKL, LDWRKR, LDWRKU, MAXWRK, MINMN, MINWRK,
  252. $ MNTHR1, MNTHR2, NRWORK, NWORK, WRKBL
  253. INTEGER LWORK_CGEBRD_MN, LWORK_CGEBRD_MM,
  254. $ LWORK_CGEBRD_NN, LWORK_CGELQF_MN,
  255. $ LWORK_CGEQRF_MN,
  256. $ LWORK_CUNGBR_P_MN, LWORK_CUNGBR_P_NN,
  257. $ LWORK_CUNGBR_Q_MN, LWORK_CUNGBR_Q_MM,
  258. $ LWORK_CUNGLQ_MN, LWORK_CUNGLQ_NN,
  259. $ LWORK_CUNGQR_MM, LWORK_CUNGQR_MN,
  260. $ LWORK_CUNMBR_PRC_MM, LWORK_CUNMBR_QLN_MM,
  261. $ LWORK_CUNMBR_PRC_MN, LWORK_CUNMBR_QLN_MN,
  262. $ LWORK_CUNMBR_PRC_NN, LWORK_CUNMBR_QLN_NN
  263. REAL ANRM, BIGNUM, EPS, SMLNUM
  264. * ..
  265. * .. Local Arrays ..
  266. INTEGER IDUM( 1 )
  267. REAL DUM( 1 )
  268. COMPLEX CDUM( 1 )
  269. * ..
  270. * .. External Subroutines ..
  271. EXTERNAL CGEBRD, CGELQF, CGEMM, CGEQRF, CLACP2, CLACPY,
  272. $ CLACRM, CLARCM, CLASCL, CLASET, CUNGBR, CUNGLQ,
  273. $ CUNGQR, CUNMBR, SBDSDC, SLASCL, XERBLA
  274. * ..
  275. * .. External Functions ..
  276. LOGICAL LSAME, SISNAN
  277. REAL SLAMCH, CLANGE, SROUNDUP_LWORK
  278. EXTERNAL LSAME, SLAMCH, CLANGE, SISNAN,
  279. $ SROUNDUP_LWORK
  280. * ..
  281. * .. Intrinsic Functions ..
  282. INTRINSIC INT, MAX, MIN, SQRT
  283. * ..
  284. * .. Executable Statements ..
  285. *
  286. * Test the input arguments
  287. *
  288. INFO = 0
  289. MINMN = MIN( M, N )
  290. MNTHR1 = INT( MINMN*17.0E0 / 9.0E0 )
  291. MNTHR2 = INT( MINMN*5.0E0 / 3.0E0 )
  292. WNTQA = LSAME( JOBZ, 'A' )
  293. WNTQS = LSAME( JOBZ, 'S' )
  294. WNTQAS = WNTQA .OR. WNTQS
  295. WNTQO = LSAME( JOBZ, 'O' )
  296. WNTQN = LSAME( JOBZ, 'N' )
  297. LQUERY = ( LWORK.EQ.-1 )
  298. MINWRK = 1
  299. MAXWRK = 1
  300. *
  301. IF( .NOT.( WNTQA .OR. WNTQS .OR. WNTQO .OR. WNTQN ) ) THEN
  302. INFO = -1
  303. ELSE IF( M.LT.0 ) THEN
  304. INFO = -2
  305. ELSE IF( N.LT.0 ) THEN
  306. INFO = -3
  307. ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
  308. INFO = -5
  309. ELSE IF( LDU.LT.1 .OR. ( WNTQAS .AND. LDU.LT.M ) .OR.
  310. $ ( WNTQO .AND. M.LT.N .AND. LDU.LT.M ) ) THEN
  311. INFO = -8
  312. ELSE IF( LDVT.LT.1 .OR. ( WNTQA .AND. LDVT.LT.N ) .OR.
  313. $ ( WNTQS .AND. LDVT.LT.MINMN ) .OR.
  314. $ ( WNTQO .AND. M.GE.N .AND. LDVT.LT.N ) ) THEN
  315. INFO = -10
  316. END IF
  317. *
  318. * Compute workspace
  319. * Note: Comments in the code beginning "Workspace:" describe the
  320. * minimal amount of workspace allocated at that point in the code,
  321. * as well as the preferred amount for good performance.
  322. * CWorkspace refers to complex workspace, and RWorkspace to
  323. * real workspace. NB refers to the optimal block size for the
  324. * immediately following subroutine, as returned by ILAENV.)
  325. *
  326. IF( INFO.EQ.0 ) THEN
  327. MINWRK = 1
  328. MAXWRK = 1
  329. IF( M.GE.N .AND. MINMN.GT.0 ) THEN
  330. *
  331. * There is no complex work space needed for bidiagonal SVD
  332. * The real work space needed for bidiagonal SVD (sbdsdc) is
  333. * BDSPAC = 3*N*N + 4*N for singular values and vectors;
  334. * BDSPAC = 4*N for singular values only;
  335. * not including e, RU, and RVT matrices.
  336. *
  337. * Compute space preferred for each routine
  338. CALL CGEBRD( M, N, CDUM(1), M, DUM(1), DUM(1), CDUM(1),
  339. $ CDUM(1), CDUM(1), -1, IERR )
  340. LWORK_CGEBRD_MN = INT( CDUM(1) )
  341. *
  342. CALL CGEBRD( N, N, CDUM(1), N, DUM(1), DUM(1), CDUM(1),
  343. $ CDUM(1), CDUM(1), -1, IERR )
  344. LWORK_CGEBRD_NN = INT( CDUM(1) )
  345. *
  346. CALL CGEQRF( M, N, CDUM(1), M, CDUM(1), CDUM(1), -1, IERR )
  347. LWORK_CGEQRF_MN = INT( CDUM(1) )
  348. *
  349. CALL CUNGBR( 'P', N, N, N, CDUM(1), N, CDUM(1), CDUM(1),
  350. $ -1, IERR )
  351. LWORK_CUNGBR_P_NN = INT( CDUM(1) )
  352. *
  353. CALL CUNGBR( 'Q', M, M, N, CDUM(1), M, CDUM(1), CDUM(1),
  354. $ -1, IERR )
  355. LWORK_CUNGBR_Q_MM = INT( CDUM(1) )
  356. *
  357. CALL CUNGBR( 'Q', M, N, N, CDUM(1), M, CDUM(1), CDUM(1),
  358. $ -1, IERR )
  359. LWORK_CUNGBR_Q_MN = INT( CDUM(1) )
  360. *
  361. CALL CUNGQR( M, M, N, CDUM(1), M, CDUM(1), CDUM(1),
  362. $ -1, IERR )
  363. LWORK_CUNGQR_MM = INT( CDUM(1) )
  364. *
  365. CALL CUNGQR( M, N, N, CDUM(1), M, CDUM(1), CDUM(1),
  366. $ -1, IERR )
  367. LWORK_CUNGQR_MN = INT( CDUM(1) )
  368. *
  369. CALL CUNMBR( 'P', 'R', 'C', N, N, N, CDUM(1), N, CDUM(1),
  370. $ CDUM(1), N, CDUM(1), -1, IERR )
  371. LWORK_CUNMBR_PRC_NN = INT( CDUM(1) )
  372. *
  373. CALL CUNMBR( 'Q', 'L', 'N', M, M, N, CDUM(1), M, CDUM(1),
  374. $ CDUM(1), M, CDUM(1), -1, IERR )
  375. LWORK_CUNMBR_QLN_MM = INT( CDUM(1) )
  376. *
  377. CALL CUNMBR( 'Q', 'L', 'N', M, N, N, CDUM(1), M, CDUM(1),
  378. $ CDUM(1), M, CDUM(1), -1, IERR )
  379. LWORK_CUNMBR_QLN_MN = INT( CDUM(1) )
  380. *
  381. CALL CUNMBR( 'Q', 'L', 'N', N, N, N, CDUM(1), N, CDUM(1),
  382. $ CDUM(1), N, CDUM(1), -1, IERR )
  383. LWORK_CUNMBR_QLN_NN = INT( CDUM(1) )
  384. *
  385. IF( M.GE.MNTHR1 ) THEN
  386. IF( WNTQN ) THEN
  387. *
  388. * Path 1 (M >> N, JOBZ='N')
  389. *
  390. MAXWRK = N + LWORK_CGEQRF_MN
  391. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CGEBRD_NN )
  392. MINWRK = 3*N
  393. ELSE IF( WNTQO ) THEN
  394. *
  395. * Path 2 (M >> N, JOBZ='O')
  396. *
  397. WRKBL = N + LWORK_CGEQRF_MN
  398. WRKBL = MAX( WRKBL, N + LWORK_CUNGQR_MN )
  399. WRKBL = MAX( WRKBL, 2*N + LWORK_CGEBRD_NN )
  400. WRKBL = MAX( WRKBL, 2*N + LWORK_CUNMBR_QLN_NN )
  401. WRKBL = MAX( WRKBL, 2*N + LWORK_CUNMBR_PRC_NN )
  402. MAXWRK = M*N + N*N + WRKBL
  403. MINWRK = 2*N*N + 3*N
  404. ELSE IF( WNTQS ) THEN
  405. *
  406. * Path 3 (M >> N, JOBZ='S')
  407. *
  408. WRKBL = N + LWORK_CGEQRF_MN
  409. WRKBL = MAX( WRKBL, N + LWORK_CUNGQR_MN )
  410. WRKBL = MAX( WRKBL, 2*N + LWORK_CGEBRD_NN )
  411. WRKBL = MAX( WRKBL, 2*N + LWORK_CUNMBR_QLN_NN )
  412. WRKBL = MAX( WRKBL, 2*N + LWORK_CUNMBR_PRC_NN )
  413. MAXWRK = N*N + WRKBL
  414. MINWRK = N*N + 3*N
  415. ELSE IF( WNTQA ) THEN
  416. *
  417. * Path 4 (M >> N, JOBZ='A')
  418. *
  419. WRKBL = N + LWORK_CGEQRF_MN
  420. WRKBL = MAX( WRKBL, N + LWORK_CUNGQR_MM )
  421. WRKBL = MAX( WRKBL, 2*N + LWORK_CGEBRD_NN )
  422. WRKBL = MAX( WRKBL, 2*N + LWORK_CUNMBR_QLN_NN )
  423. WRKBL = MAX( WRKBL, 2*N + LWORK_CUNMBR_PRC_NN )
  424. MAXWRK = N*N + WRKBL
  425. MINWRK = N*N + MAX( 3*N, N + M )
  426. END IF
  427. ELSE IF( M.GE.MNTHR2 ) THEN
  428. *
  429. * Path 5 (M >> N, but not as much as MNTHR1)
  430. *
  431. MAXWRK = 2*N + LWORK_CGEBRD_MN
  432. MINWRK = 2*N + M
  433. IF( WNTQO ) THEN
  434. * Path 5o (M >> N, JOBZ='O')
  435. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNGBR_P_NN )
  436. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNGBR_Q_MN )
  437. MAXWRK = MAXWRK + M*N
  438. MINWRK = MINWRK + N*N
  439. ELSE IF( WNTQS ) THEN
  440. * Path 5s (M >> N, JOBZ='S')
  441. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNGBR_P_NN )
  442. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNGBR_Q_MN )
  443. ELSE IF( WNTQA ) THEN
  444. * Path 5a (M >> N, JOBZ='A')
  445. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNGBR_P_NN )
  446. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNGBR_Q_MM )
  447. END IF
  448. ELSE
  449. *
  450. * Path 6 (M >= N, but not much larger)
  451. *
  452. MAXWRK = 2*N + LWORK_CGEBRD_MN
  453. MINWRK = 2*N + M
  454. IF( WNTQO ) THEN
  455. * Path 6o (M >= N, JOBZ='O')
  456. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNMBR_PRC_NN )
  457. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNMBR_QLN_MN )
  458. MAXWRK = MAXWRK + M*N
  459. MINWRK = MINWRK + N*N
  460. ELSE IF( WNTQS ) THEN
  461. * Path 6s (M >= N, JOBZ='S')
  462. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNMBR_QLN_MN )
  463. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNMBR_PRC_NN )
  464. ELSE IF( WNTQA ) THEN
  465. * Path 6a (M >= N, JOBZ='A')
  466. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNMBR_QLN_MM )
  467. MAXWRK = MAX( MAXWRK, 2*N + LWORK_CUNMBR_PRC_NN )
  468. END IF
  469. END IF
  470. ELSE IF( MINMN.GT.0 ) THEN
  471. *
  472. * There is no complex work space needed for bidiagonal SVD
  473. * The real work space needed for bidiagonal SVD (sbdsdc) is
  474. * BDSPAC = 3*M*M + 4*M for singular values and vectors;
  475. * BDSPAC = 4*M for singular values only;
  476. * not including e, RU, and RVT matrices.
  477. *
  478. * Compute space preferred for each routine
  479. CALL CGEBRD( M, N, CDUM(1), M, DUM(1), DUM(1), CDUM(1),
  480. $ CDUM(1), CDUM(1), -1, IERR )
  481. LWORK_CGEBRD_MN = INT( CDUM(1) )
  482. *
  483. CALL CGEBRD( M, M, CDUM(1), M, DUM(1), DUM(1), CDUM(1),
  484. $ CDUM(1), CDUM(1), -1, IERR )
  485. LWORK_CGEBRD_MM = INT( CDUM(1) )
  486. *
  487. CALL CGELQF( M, N, CDUM(1), M, CDUM(1), CDUM(1), -1, IERR )
  488. LWORK_CGELQF_MN = INT( CDUM(1) )
  489. *
  490. CALL CUNGBR( 'P', M, N, M, CDUM(1), M, CDUM(1), CDUM(1),
  491. $ -1, IERR )
  492. LWORK_CUNGBR_P_MN = INT( CDUM(1) )
  493. *
  494. CALL CUNGBR( 'P', N, N, M, CDUM(1), N, CDUM(1), CDUM(1),
  495. $ -1, IERR )
  496. LWORK_CUNGBR_P_NN = INT( CDUM(1) )
  497. *
  498. CALL CUNGBR( 'Q', M, M, N, CDUM(1), M, CDUM(1), CDUM(1),
  499. $ -1, IERR )
  500. LWORK_CUNGBR_Q_MM = INT( CDUM(1) )
  501. *
  502. CALL CUNGLQ( M, N, M, CDUM(1), M, CDUM(1), CDUM(1),
  503. $ -1, IERR )
  504. LWORK_CUNGLQ_MN = INT( CDUM(1) )
  505. *
  506. CALL CUNGLQ( N, N, M, CDUM(1), N, CDUM(1), CDUM(1),
  507. $ -1, IERR )
  508. LWORK_CUNGLQ_NN = INT( CDUM(1) )
  509. *
  510. CALL CUNMBR( 'P', 'R', 'C', M, M, M, CDUM(1), M, CDUM(1),
  511. $ CDUM(1), M, CDUM(1), -1, IERR )
  512. LWORK_CUNMBR_PRC_MM = INT( CDUM(1) )
  513. *
  514. CALL CUNMBR( 'P', 'R', 'C', M, N, M, CDUM(1), M, CDUM(1),
  515. $ CDUM(1), M, CDUM(1), -1, IERR )
  516. LWORK_CUNMBR_PRC_MN = INT( CDUM(1) )
  517. *
  518. CALL CUNMBR( 'P', 'R', 'C', N, N, M, CDUM(1), N, CDUM(1),
  519. $ CDUM(1), N, CDUM(1), -1, IERR )
  520. LWORK_CUNMBR_PRC_NN = INT( CDUM(1) )
  521. *
  522. CALL CUNMBR( 'Q', 'L', 'N', M, M, M, CDUM(1), M, CDUM(1),
  523. $ CDUM(1), M, CDUM(1), -1, IERR )
  524. LWORK_CUNMBR_QLN_MM = INT( CDUM(1) )
  525. *
  526. IF( N.GE.MNTHR1 ) THEN
  527. IF( WNTQN ) THEN
  528. *
  529. * Path 1t (N >> M, JOBZ='N')
  530. *
  531. MAXWRK = M + LWORK_CGELQF_MN
  532. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CGEBRD_MM )
  533. MINWRK = 3*M
  534. ELSE IF( WNTQO ) THEN
  535. *
  536. * Path 2t (N >> M, JOBZ='O')
  537. *
  538. WRKBL = M + LWORK_CGELQF_MN
  539. WRKBL = MAX( WRKBL, M + LWORK_CUNGLQ_MN )
  540. WRKBL = MAX( WRKBL, 2*M + LWORK_CGEBRD_MM )
  541. WRKBL = MAX( WRKBL, 2*M + LWORK_CUNMBR_QLN_MM )
  542. WRKBL = MAX( WRKBL, 2*M + LWORK_CUNMBR_PRC_MM )
  543. MAXWRK = M*N + M*M + WRKBL
  544. MINWRK = 2*M*M + 3*M
  545. ELSE IF( WNTQS ) THEN
  546. *
  547. * Path 3t (N >> M, JOBZ='S')
  548. *
  549. WRKBL = M + LWORK_CGELQF_MN
  550. WRKBL = MAX( WRKBL, M + LWORK_CUNGLQ_MN )
  551. WRKBL = MAX( WRKBL, 2*M + LWORK_CGEBRD_MM )
  552. WRKBL = MAX( WRKBL, 2*M + LWORK_CUNMBR_QLN_MM )
  553. WRKBL = MAX( WRKBL, 2*M + LWORK_CUNMBR_PRC_MM )
  554. MAXWRK = M*M + WRKBL
  555. MINWRK = M*M + 3*M
  556. ELSE IF( WNTQA ) THEN
  557. *
  558. * Path 4t (N >> M, JOBZ='A')
  559. *
  560. WRKBL = M + LWORK_CGELQF_MN
  561. WRKBL = MAX( WRKBL, M + LWORK_CUNGLQ_NN )
  562. WRKBL = MAX( WRKBL, 2*M + LWORK_CGEBRD_MM )
  563. WRKBL = MAX( WRKBL, 2*M + LWORK_CUNMBR_QLN_MM )
  564. WRKBL = MAX( WRKBL, 2*M + LWORK_CUNMBR_PRC_MM )
  565. MAXWRK = M*M + WRKBL
  566. MINWRK = M*M + MAX( 3*M, M + N )
  567. END IF
  568. ELSE IF( N.GE.MNTHR2 ) THEN
  569. *
  570. * Path 5t (N >> M, but not as much as MNTHR1)
  571. *
  572. MAXWRK = 2*M + LWORK_CGEBRD_MN
  573. MINWRK = 2*M + N
  574. IF( WNTQO ) THEN
  575. * Path 5to (N >> M, JOBZ='O')
  576. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNGBR_Q_MM )
  577. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNGBR_P_MN )
  578. MAXWRK = MAXWRK + M*N
  579. MINWRK = MINWRK + M*M
  580. ELSE IF( WNTQS ) THEN
  581. * Path 5ts (N >> M, JOBZ='S')
  582. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNGBR_Q_MM )
  583. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNGBR_P_MN )
  584. ELSE IF( WNTQA ) THEN
  585. * Path 5ta (N >> M, JOBZ='A')
  586. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNGBR_Q_MM )
  587. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNGBR_P_NN )
  588. END IF
  589. ELSE
  590. *
  591. * Path 6t (N > M, but not much larger)
  592. *
  593. MAXWRK = 2*M + LWORK_CGEBRD_MN
  594. MINWRK = 2*M + N
  595. IF( WNTQO ) THEN
  596. * Path 6to (N > M, JOBZ='O')
  597. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNMBR_QLN_MM )
  598. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNMBR_PRC_MN )
  599. MAXWRK = MAXWRK + M*N
  600. MINWRK = MINWRK + M*M
  601. ELSE IF( WNTQS ) THEN
  602. * Path 6ts (N > M, JOBZ='S')
  603. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNMBR_QLN_MM )
  604. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNMBR_PRC_MN )
  605. ELSE IF( WNTQA ) THEN
  606. * Path 6ta (N > M, JOBZ='A')
  607. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNMBR_QLN_MM )
  608. MAXWRK = MAX( MAXWRK, 2*M + LWORK_CUNMBR_PRC_NN )
  609. END IF
  610. END IF
  611. END IF
  612. MAXWRK = MAX( MAXWRK, MINWRK )
  613. END IF
  614. IF( INFO.EQ.0 ) THEN
  615. WORK( 1 ) = SROUNDUP_LWORK( MAXWRK )
  616. IF( LWORK.LT.MINWRK .AND. .NOT. LQUERY ) THEN
  617. INFO = -12
  618. END IF
  619. END IF
  620. *
  621. IF( INFO.NE.0 ) THEN
  622. CALL XERBLA( 'CGESDD', -INFO )
  623. RETURN
  624. ELSE IF( LQUERY ) THEN
  625. RETURN
  626. END IF
  627. *
  628. * Quick return if possible
  629. *
  630. IF( M.EQ.0 .OR. N.EQ.0 ) THEN
  631. RETURN
  632. END IF
  633. *
  634. * Get machine constants
  635. *
  636. EPS = SLAMCH( 'P' )
  637. SMLNUM = SQRT( SLAMCH( 'S' ) ) / EPS
  638. BIGNUM = ONE / SMLNUM
  639. *
  640. * Scale A if max element outside range [SMLNUM,BIGNUM]
  641. *
  642. ANRM = CLANGE( 'M', M, N, A, LDA, DUM )
  643. IF( SISNAN ( ANRM ) ) THEN
  644. INFO = -4
  645. RETURN
  646. END IF
  647. ISCL = 0
  648. IF( ANRM.GT.ZERO .AND. ANRM.LT.SMLNUM ) THEN
  649. ISCL = 1
  650. CALL CLASCL( 'G', 0, 0, ANRM, SMLNUM, M, N, A, LDA, IERR )
  651. ELSE IF( ANRM.GT.BIGNUM ) THEN
  652. ISCL = 1
  653. CALL CLASCL( 'G', 0, 0, ANRM, BIGNUM, M, N, A, LDA, IERR )
  654. END IF
  655. *
  656. IF( M.GE.N ) THEN
  657. *
  658. * A has at least as many rows as columns. If A has sufficiently
  659. * more rows than columns, first reduce using the QR
  660. * decomposition (if sufficient workspace available)
  661. *
  662. IF( M.GE.MNTHR1 ) THEN
  663. *
  664. IF( WNTQN ) THEN
  665. *
  666. * Path 1 (M >> N, JOBZ='N')
  667. * No singular vectors to be computed
  668. *
  669. ITAU = 1
  670. NWORK = ITAU + N
  671. *
  672. * Compute A=Q*R
  673. * CWorkspace: need N [tau] + N [work]
  674. * CWorkspace: prefer N [tau] + N*NB [work]
  675. * RWorkspace: need 0
  676. *
  677. CALL CGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  678. $ LWORK-NWORK+1, IERR )
  679. *
  680. * Zero out below R
  681. *
  682. CALL CLASET( 'L', N-1, N-1, CZERO, CZERO, A( 2, 1 ),
  683. $ LDA )
  684. IE = 1
  685. ITAUQ = 1
  686. ITAUP = ITAUQ + N
  687. NWORK = ITAUP + N
  688. *
  689. * Bidiagonalize R in A
  690. * CWorkspace: need 2*N [tauq, taup] + N [work]
  691. * CWorkspace: prefer 2*N [tauq, taup] + 2*N*NB [work]
  692. * RWorkspace: need N [e]
  693. *
  694. CALL CGEBRD( N, N, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  695. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  696. $ IERR )
  697. NRWORK = IE + N
  698. *
  699. * Perform bidiagonal SVD, compute singular values only
  700. * CWorkspace: need 0
  701. * RWorkspace: need N [e] + BDSPAC
  702. *
  703. CALL SBDSDC( 'U', 'N', N, S, RWORK( IE ), DUM,1,DUM,1,
  704. $ DUM, IDUM, RWORK( NRWORK ), IWORK, INFO )
  705. *
  706. ELSE IF( WNTQO ) THEN
  707. *
  708. * Path 2 (M >> N, JOBZ='O')
  709. * N left singular vectors to be overwritten on A and
  710. * N right singular vectors to be computed in VT
  711. *
  712. IU = 1
  713. *
  714. * WORK(IU) is N by N
  715. *
  716. LDWRKU = N
  717. IR = IU + LDWRKU*N
  718. IF( LWORK .GE. M*N + N*N + 3*N ) THEN
  719. *
  720. * WORK(IR) is M by N
  721. *
  722. LDWRKR = M
  723. ELSE
  724. LDWRKR = ( LWORK - N*N - 3*N ) / N
  725. END IF
  726. ITAU = IR + LDWRKR*N
  727. NWORK = ITAU + N
  728. *
  729. * Compute A=Q*R
  730. * CWorkspace: need N*N [U] + N*N [R] + N [tau] + N [work]
  731. * CWorkspace: prefer N*N [U] + N*N [R] + N [tau] + N*NB [work]
  732. * RWorkspace: need 0
  733. *
  734. CALL CGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  735. $ LWORK-NWORK+1, IERR )
  736. *
  737. * Copy R to WORK( IR ), zeroing out below it
  738. *
  739. CALL CLACPY( 'U', N, N, A, LDA, WORK( IR ), LDWRKR )
  740. CALL CLASET( 'L', N-1, N-1, CZERO, CZERO, WORK( IR+1 ),
  741. $ LDWRKR )
  742. *
  743. * Generate Q in A
  744. * CWorkspace: need N*N [U] + N*N [R] + N [tau] + N [work]
  745. * CWorkspace: prefer N*N [U] + N*N [R] + N [tau] + N*NB [work]
  746. * RWorkspace: need 0
  747. *
  748. CALL CUNGQR( M, N, N, A, LDA, WORK( ITAU ),
  749. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  750. IE = 1
  751. ITAUQ = ITAU
  752. ITAUP = ITAUQ + N
  753. NWORK = ITAUP + N
  754. *
  755. * Bidiagonalize R in WORK(IR)
  756. * CWorkspace: need N*N [U] + N*N [R] + 2*N [tauq, taup] + N [work]
  757. * CWorkspace: prefer N*N [U] + N*N [R] + 2*N [tauq, taup] + 2*N*NB [work]
  758. * RWorkspace: need N [e]
  759. *
  760. CALL CGEBRD( N, N, WORK( IR ), LDWRKR, S, RWORK( IE ),
  761. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  762. $ LWORK-NWORK+1, IERR )
  763. *
  764. * Perform bidiagonal SVD, computing left singular vectors
  765. * of R in WORK(IRU) and computing right singular vectors
  766. * of R in WORK(IRVT)
  767. * CWorkspace: need 0
  768. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  769. *
  770. IRU = IE + N
  771. IRVT = IRU + N*N
  772. NRWORK = IRVT + N*N
  773. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  774. $ N, RWORK( IRVT ), N, DUM, IDUM,
  775. $ RWORK( NRWORK ), IWORK, INFO )
  776. *
  777. * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
  778. * Overwrite WORK(IU) by the left singular vectors of R
  779. * CWorkspace: need N*N [U] + N*N [R] + 2*N [tauq, taup] + N [work]
  780. * CWorkspace: prefer N*N [U] + N*N [R] + 2*N [tauq, taup] + N*NB [work]
  781. * RWorkspace: need 0
  782. *
  783. CALL CLACP2( 'F', N, N, RWORK( IRU ), N, WORK( IU ),
  784. $ LDWRKU )
  785. CALL CUNMBR( 'Q', 'L', 'N', N, N, N, WORK( IR ), LDWRKR,
  786. $ WORK( ITAUQ ), WORK( IU ), LDWRKU,
  787. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  788. *
  789. * Copy real matrix RWORK(IRVT) to complex matrix VT
  790. * Overwrite VT by the right singular vectors of R
  791. * CWorkspace: need N*N [U] + N*N [R] + 2*N [tauq, taup] + N [work]
  792. * CWorkspace: prefer N*N [U] + N*N [R] + 2*N [tauq, taup] + N*NB [work]
  793. * RWorkspace: need 0
  794. *
  795. CALL CLACP2( 'F', N, N, RWORK( IRVT ), N, VT, LDVT )
  796. CALL CUNMBR( 'P', 'R', 'C', N, N, N, WORK( IR ), LDWRKR,
  797. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  798. $ LWORK-NWORK+1, IERR )
  799. *
  800. * Multiply Q in A by left singular vectors of R in
  801. * WORK(IU), storing result in WORK(IR) and copying to A
  802. * CWorkspace: need N*N [U] + N*N [R]
  803. * CWorkspace: prefer N*N [U] + M*N [R]
  804. * RWorkspace: need 0
  805. *
  806. DO 10 I = 1, M, LDWRKR
  807. CHUNK = MIN( M-I+1, LDWRKR )
  808. CALL CGEMM( 'N', 'N', CHUNK, N, N, CONE, A( I, 1 ),
  809. $ LDA, WORK( IU ), LDWRKU, CZERO,
  810. $ WORK( IR ), LDWRKR )
  811. CALL CLACPY( 'F', CHUNK, N, WORK( IR ), LDWRKR,
  812. $ A( I, 1 ), LDA )
  813. 10 CONTINUE
  814. *
  815. ELSE IF( WNTQS ) THEN
  816. *
  817. * Path 3 (M >> N, JOBZ='S')
  818. * N left singular vectors to be computed in U and
  819. * N right singular vectors to be computed in VT
  820. *
  821. IR = 1
  822. *
  823. * WORK(IR) is N by N
  824. *
  825. LDWRKR = N
  826. ITAU = IR + LDWRKR*N
  827. NWORK = ITAU + N
  828. *
  829. * Compute A=Q*R
  830. * CWorkspace: need N*N [R] + N [tau] + N [work]
  831. * CWorkspace: prefer N*N [R] + N [tau] + N*NB [work]
  832. * RWorkspace: need 0
  833. *
  834. CALL CGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  835. $ LWORK-NWORK+1, IERR )
  836. *
  837. * Copy R to WORK(IR), zeroing out below it
  838. *
  839. CALL CLACPY( 'U', N, N, A, LDA, WORK( IR ), LDWRKR )
  840. CALL CLASET( 'L', N-1, N-1, CZERO, CZERO, WORK( IR+1 ),
  841. $ LDWRKR )
  842. *
  843. * Generate Q in A
  844. * CWorkspace: need N*N [R] + N [tau] + N [work]
  845. * CWorkspace: prefer N*N [R] + N [tau] + N*NB [work]
  846. * RWorkspace: need 0
  847. *
  848. CALL CUNGQR( M, N, N, A, LDA, WORK( ITAU ),
  849. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  850. IE = 1
  851. ITAUQ = ITAU
  852. ITAUP = ITAUQ + N
  853. NWORK = ITAUP + N
  854. *
  855. * Bidiagonalize R in WORK(IR)
  856. * CWorkspace: need N*N [R] + 2*N [tauq, taup] + N [work]
  857. * CWorkspace: prefer N*N [R] + 2*N [tauq, taup] + 2*N*NB [work]
  858. * RWorkspace: need N [e]
  859. *
  860. CALL CGEBRD( N, N, WORK( IR ), LDWRKR, S, RWORK( IE ),
  861. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  862. $ LWORK-NWORK+1, IERR )
  863. *
  864. * Perform bidiagonal SVD, computing left singular vectors
  865. * of bidiagonal matrix in RWORK(IRU) and computing right
  866. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  867. * CWorkspace: need 0
  868. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  869. *
  870. IRU = IE + N
  871. IRVT = IRU + N*N
  872. NRWORK = IRVT + N*N
  873. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  874. $ N, RWORK( IRVT ), N, DUM, IDUM,
  875. $ RWORK( NRWORK ), IWORK, INFO )
  876. *
  877. * Copy real matrix RWORK(IRU) to complex matrix U
  878. * Overwrite U by left singular vectors of R
  879. * CWorkspace: need N*N [R] + 2*N [tauq, taup] + N [work]
  880. * CWorkspace: prefer N*N [R] + 2*N [tauq, taup] + N*NB [work]
  881. * RWorkspace: need 0
  882. *
  883. CALL CLACP2( 'F', N, N, RWORK( IRU ), N, U, LDU )
  884. CALL CUNMBR( 'Q', 'L', 'N', N, N, N, WORK( IR ), LDWRKR,
  885. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  886. $ LWORK-NWORK+1, IERR )
  887. *
  888. * Copy real matrix RWORK(IRVT) to complex matrix VT
  889. * Overwrite VT by right singular vectors of R
  890. * CWorkspace: need N*N [R] + 2*N [tauq, taup] + N [work]
  891. * CWorkspace: prefer N*N [R] + 2*N [tauq, taup] + N*NB [work]
  892. * RWorkspace: need 0
  893. *
  894. CALL CLACP2( 'F', N, N, RWORK( IRVT ), N, VT, LDVT )
  895. CALL CUNMBR( 'P', 'R', 'C', N, N, N, WORK( IR ), LDWRKR,
  896. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  897. $ LWORK-NWORK+1, IERR )
  898. *
  899. * Multiply Q in A by left singular vectors of R in
  900. * WORK(IR), storing result in U
  901. * CWorkspace: need N*N [R]
  902. * RWorkspace: need 0
  903. *
  904. CALL CLACPY( 'F', N, N, U, LDU, WORK( IR ), LDWRKR )
  905. CALL CGEMM( 'N', 'N', M, N, N, CONE, A, LDA, WORK( IR ),
  906. $ LDWRKR, CZERO, U, LDU )
  907. *
  908. ELSE IF( WNTQA ) THEN
  909. *
  910. * Path 4 (M >> N, JOBZ='A')
  911. * M left singular vectors to be computed in U and
  912. * N right singular vectors to be computed in VT
  913. *
  914. IU = 1
  915. *
  916. * WORK(IU) is N by N
  917. *
  918. LDWRKU = N
  919. ITAU = IU + LDWRKU*N
  920. NWORK = ITAU + N
  921. *
  922. * Compute A=Q*R, copying result to U
  923. * CWorkspace: need N*N [U] + N [tau] + N [work]
  924. * CWorkspace: prefer N*N [U] + N [tau] + N*NB [work]
  925. * RWorkspace: need 0
  926. *
  927. CALL CGEQRF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  928. $ LWORK-NWORK+1, IERR )
  929. CALL CLACPY( 'L', M, N, A, LDA, U, LDU )
  930. *
  931. * Generate Q in U
  932. * CWorkspace: need N*N [U] + N [tau] + M [work]
  933. * CWorkspace: prefer N*N [U] + N [tau] + M*NB [work]
  934. * RWorkspace: need 0
  935. *
  936. CALL CUNGQR( M, M, N, U, LDU, WORK( ITAU ),
  937. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  938. *
  939. * Produce R in A, zeroing out below it
  940. *
  941. CALL CLASET( 'L', N-1, N-1, CZERO, CZERO, A( 2, 1 ),
  942. $ LDA )
  943. IE = 1
  944. ITAUQ = ITAU
  945. ITAUP = ITAUQ + N
  946. NWORK = ITAUP + N
  947. *
  948. * Bidiagonalize R in A
  949. * CWorkspace: need N*N [U] + 2*N [tauq, taup] + N [work]
  950. * CWorkspace: prefer N*N [U] + 2*N [tauq, taup] + 2*N*NB [work]
  951. * RWorkspace: need N [e]
  952. *
  953. CALL CGEBRD( N, N, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  954. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  955. $ IERR )
  956. IRU = IE + N
  957. IRVT = IRU + N*N
  958. NRWORK = IRVT + N*N
  959. *
  960. * Perform bidiagonal SVD, computing left singular vectors
  961. * of bidiagonal matrix in RWORK(IRU) and computing right
  962. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  963. * CWorkspace: need 0
  964. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  965. *
  966. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  967. $ N, RWORK( IRVT ), N, DUM, IDUM,
  968. $ RWORK( NRWORK ), IWORK, INFO )
  969. *
  970. * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
  971. * Overwrite WORK(IU) by left singular vectors of R
  972. * CWorkspace: need N*N [U] + 2*N [tauq, taup] + N [work]
  973. * CWorkspace: prefer N*N [U] + 2*N [tauq, taup] + N*NB [work]
  974. * RWorkspace: need 0
  975. *
  976. CALL CLACP2( 'F', N, N, RWORK( IRU ), N, WORK( IU ),
  977. $ LDWRKU )
  978. CALL CUNMBR( 'Q', 'L', 'N', N, N, N, A, LDA,
  979. $ WORK( ITAUQ ), WORK( IU ), LDWRKU,
  980. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  981. *
  982. * Copy real matrix RWORK(IRVT) to complex matrix VT
  983. * Overwrite VT by right singular vectors of R
  984. * CWorkspace: need N*N [U] + 2*N [tauq, taup] + N [work]
  985. * CWorkspace: prefer N*N [U] + 2*N [tauq, taup] + N*NB [work]
  986. * RWorkspace: need 0
  987. *
  988. CALL CLACP2( 'F', N, N, RWORK( IRVT ), N, VT, LDVT )
  989. CALL CUNMBR( 'P', 'R', 'C', N, N, N, A, LDA,
  990. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  991. $ LWORK-NWORK+1, IERR )
  992. *
  993. * Multiply Q in U by left singular vectors of R in
  994. * WORK(IU), storing result in A
  995. * CWorkspace: need N*N [U]
  996. * RWorkspace: need 0
  997. *
  998. CALL CGEMM( 'N', 'N', M, N, N, CONE, U, LDU, WORK( IU ),
  999. $ LDWRKU, CZERO, A, LDA )
  1000. *
  1001. * Copy left singular vectors of A from A to U
  1002. *
  1003. CALL CLACPY( 'F', M, N, A, LDA, U, LDU )
  1004. *
  1005. END IF
  1006. *
  1007. ELSE IF( M.GE.MNTHR2 ) THEN
  1008. *
  1009. * MNTHR2 <= M < MNTHR1
  1010. *
  1011. * Path 5 (M >> N, but not as much as MNTHR1)
  1012. * Reduce to bidiagonal form without QR decomposition, use
  1013. * CUNGBR and matrix multiplication to compute singular vectors
  1014. *
  1015. IE = 1
  1016. NRWORK = IE + N
  1017. ITAUQ = 1
  1018. ITAUP = ITAUQ + N
  1019. NWORK = ITAUP + N
  1020. *
  1021. * Bidiagonalize A
  1022. * CWorkspace: need 2*N [tauq, taup] + M [work]
  1023. * CWorkspace: prefer 2*N [tauq, taup] + (M+N)*NB [work]
  1024. * RWorkspace: need N [e]
  1025. *
  1026. CALL CGEBRD( M, N, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  1027. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1028. $ IERR )
  1029. IF( WNTQN ) THEN
  1030. *
  1031. * Path 5n (M >> N, JOBZ='N')
  1032. * Compute singular values only
  1033. * CWorkspace: need 0
  1034. * RWorkspace: need N [e] + BDSPAC
  1035. *
  1036. CALL SBDSDC( 'U', 'N', N, S, RWORK( IE ), DUM, 1,DUM,1,
  1037. $ DUM, IDUM, RWORK( NRWORK ), IWORK, INFO )
  1038. ELSE IF( WNTQO ) THEN
  1039. IU = NWORK
  1040. IRU = NRWORK
  1041. IRVT = IRU + N*N
  1042. NRWORK = IRVT + N*N
  1043. *
  1044. * Path 5o (M >> N, JOBZ='O')
  1045. * Copy A to VT, generate P**H
  1046. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1047. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1048. * RWorkspace: need 0
  1049. *
  1050. CALL CLACPY( 'U', N, N, A, LDA, VT, LDVT )
  1051. CALL CUNGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
  1052. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1053. *
  1054. * Generate Q in A
  1055. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1056. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1057. * RWorkspace: need 0
  1058. *
  1059. CALL CUNGBR( 'Q', M, N, N, A, LDA, WORK( ITAUQ ),
  1060. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1061. *
  1062. IF( LWORK .GE. M*N + 3*N ) THEN
  1063. *
  1064. * WORK( IU ) is M by N
  1065. *
  1066. LDWRKU = M
  1067. ELSE
  1068. *
  1069. * WORK(IU) is LDWRKU by N
  1070. *
  1071. LDWRKU = ( LWORK - 3*N ) / N
  1072. END IF
  1073. NWORK = IU + LDWRKU*N
  1074. *
  1075. * Perform bidiagonal SVD, computing left singular vectors
  1076. * of bidiagonal matrix in RWORK(IRU) and computing right
  1077. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1078. * CWorkspace: need 0
  1079. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  1080. *
  1081. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  1082. $ N, RWORK( IRVT ), N, DUM, IDUM,
  1083. $ RWORK( NRWORK ), IWORK, INFO )
  1084. *
  1085. * Multiply real matrix RWORK(IRVT) by P**H in VT,
  1086. * storing the result in WORK(IU), copying to VT
  1087. * CWorkspace: need 2*N [tauq, taup] + N*N [U]
  1088. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + 2*N*N [rwork]
  1089. *
  1090. CALL CLARCM( N, N, RWORK( IRVT ), N, VT, LDVT,
  1091. $ WORK( IU ), LDWRKU, RWORK( NRWORK ) )
  1092. CALL CLACPY( 'F', N, N, WORK( IU ), LDWRKU, VT, LDVT )
  1093. *
  1094. * Multiply Q in A by real matrix RWORK(IRU), storing the
  1095. * result in WORK(IU), copying to A
  1096. * CWorkspace: need 2*N [tauq, taup] + N*N [U]
  1097. * CWorkspace: prefer 2*N [tauq, taup] + M*N [U]
  1098. * RWorkspace: need N [e] + N*N [RU] + 2*N*N [rwork]
  1099. * RWorkspace: prefer N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
  1100. *
  1101. NRWORK = IRVT
  1102. DO 20 I = 1, M, LDWRKU
  1103. CHUNK = MIN( M-I+1, LDWRKU )
  1104. CALL CLACRM( CHUNK, N, A( I, 1 ), LDA, RWORK( IRU ),
  1105. $ N, WORK( IU ), LDWRKU, RWORK( NRWORK ) )
  1106. CALL CLACPY( 'F', CHUNK, N, WORK( IU ), LDWRKU,
  1107. $ A( I, 1 ), LDA )
  1108. 20 CONTINUE
  1109. *
  1110. ELSE IF( WNTQS ) THEN
  1111. *
  1112. * Path 5s (M >> N, JOBZ='S')
  1113. * Copy A to VT, generate P**H
  1114. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1115. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1116. * RWorkspace: need 0
  1117. *
  1118. CALL CLACPY( 'U', N, N, A, LDA, VT, LDVT )
  1119. CALL CUNGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
  1120. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1121. *
  1122. * Copy A to U, generate Q
  1123. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1124. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1125. * RWorkspace: need 0
  1126. *
  1127. CALL CLACPY( 'L', M, N, A, LDA, U, LDU )
  1128. CALL CUNGBR( 'Q', M, N, N, U, LDU, WORK( ITAUQ ),
  1129. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1130. *
  1131. * Perform bidiagonal SVD, computing left singular vectors
  1132. * of bidiagonal matrix in RWORK(IRU) and computing right
  1133. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1134. * CWorkspace: need 0
  1135. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  1136. *
  1137. IRU = NRWORK
  1138. IRVT = IRU + N*N
  1139. NRWORK = IRVT + N*N
  1140. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  1141. $ N, RWORK( IRVT ), N, DUM, IDUM,
  1142. $ RWORK( NRWORK ), IWORK, INFO )
  1143. *
  1144. * Multiply real matrix RWORK(IRVT) by P**H in VT,
  1145. * storing the result in A, copying to VT
  1146. * CWorkspace: need 0
  1147. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + 2*N*N [rwork]
  1148. *
  1149. CALL CLARCM( N, N, RWORK( IRVT ), N, VT, LDVT, A, LDA,
  1150. $ RWORK( NRWORK ) )
  1151. CALL CLACPY( 'F', N, N, A, LDA, VT, LDVT )
  1152. *
  1153. * Multiply Q in U by real matrix RWORK(IRU), storing the
  1154. * result in A, copying to U
  1155. * CWorkspace: need 0
  1156. * RWorkspace: need N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
  1157. *
  1158. NRWORK = IRVT
  1159. CALL CLACRM( M, N, U, LDU, RWORK( IRU ), N, A, LDA,
  1160. $ RWORK( NRWORK ) )
  1161. CALL CLACPY( 'F', M, N, A, LDA, U, LDU )
  1162. ELSE
  1163. *
  1164. * Path 5a (M >> N, JOBZ='A')
  1165. * Copy A to VT, generate P**H
  1166. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1167. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1168. * RWorkspace: need 0
  1169. *
  1170. CALL CLACPY( 'U', N, N, A, LDA, VT, LDVT )
  1171. CALL CUNGBR( 'P', N, N, N, VT, LDVT, WORK( ITAUP ),
  1172. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1173. *
  1174. * Copy A to U, generate Q
  1175. * CWorkspace: need 2*N [tauq, taup] + M [work]
  1176. * CWorkspace: prefer 2*N [tauq, taup] + M*NB [work]
  1177. * RWorkspace: need 0
  1178. *
  1179. CALL CLACPY( 'L', M, N, A, LDA, U, LDU )
  1180. CALL CUNGBR( 'Q', M, M, N, U, LDU, WORK( ITAUQ ),
  1181. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1182. *
  1183. * Perform bidiagonal SVD, computing left singular vectors
  1184. * of bidiagonal matrix in RWORK(IRU) and computing right
  1185. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1186. * CWorkspace: need 0
  1187. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  1188. *
  1189. IRU = NRWORK
  1190. IRVT = IRU + N*N
  1191. NRWORK = IRVT + N*N
  1192. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  1193. $ N, RWORK( IRVT ), N, DUM, IDUM,
  1194. $ RWORK( NRWORK ), IWORK, INFO )
  1195. *
  1196. * Multiply real matrix RWORK(IRVT) by P**H in VT,
  1197. * storing the result in A, copying to VT
  1198. * CWorkspace: need 0
  1199. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + 2*N*N [rwork]
  1200. *
  1201. CALL CLARCM( N, N, RWORK( IRVT ), N, VT, LDVT, A, LDA,
  1202. $ RWORK( NRWORK ) )
  1203. CALL CLACPY( 'F', N, N, A, LDA, VT, LDVT )
  1204. *
  1205. * Multiply Q in U by real matrix RWORK(IRU), storing the
  1206. * result in A, copying to U
  1207. * CWorkspace: need 0
  1208. * RWorkspace: need N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
  1209. *
  1210. NRWORK = IRVT
  1211. CALL CLACRM( M, N, U, LDU, RWORK( IRU ), N, A, LDA,
  1212. $ RWORK( NRWORK ) )
  1213. CALL CLACPY( 'F', M, N, A, LDA, U, LDU )
  1214. END IF
  1215. *
  1216. ELSE
  1217. *
  1218. * M .LT. MNTHR2
  1219. *
  1220. * Path 6 (M >= N, but not much larger)
  1221. * Reduce to bidiagonal form without QR decomposition
  1222. * Use CUNMBR to compute singular vectors
  1223. *
  1224. IE = 1
  1225. NRWORK = IE + N
  1226. ITAUQ = 1
  1227. ITAUP = ITAUQ + N
  1228. NWORK = ITAUP + N
  1229. *
  1230. * Bidiagonalize A
  1231. * CWorkspace: need 2*N [tauq, taup] + M [work]
  1232. * CWorkspace: prefer 2*N [tauq, taup] + (M+N)*NB [work]
  1233. * RWorkspace: need N [e]
  1234. *
  1235. CALL CGEBRD( M, N, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  1236. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1237. $ IERR )
  1238. IF( WNTQN ) THEN
  1239. *
  1240. * Path 6n (M >= N, JOBZ='N')
  1241. * Compute singular values only
  1242. * CWorkspace: need 0
  1243. * RWorkspace: need N [e] + BDSPAC
  1244. *
  1245. CALL SBDSDC( 'U', 'N', N, S, RWORK( IE ), DUM,1,DUM,1,
  1246. $ DUM, IDUM, RWORK( NRWORK ), IWORK, INFO )
  1247. ELSE IF( WNTQO ) THEN
  1248. IU = NWORK
  1249. IRU = NRWORK
  1250. IRVT = IRU + N*N
  1251. NRWORK = IRVT + N*N
  1252. IF( LWORK .GE. M*N + 3*N ) THEN
  1253. *
  1254. * WORK( IU ) is M by N
  1255. *
  1256. LDWRKU = M
  1257. ELSE
  1258. *
  1259. * WORK( IU ) is LDWRKU by N
  1260. *
  1261. LDWRKU = ( LWORK - 3*N ) / N
  1262. END IF
  1263. NWORK = IU + LDWRKU*N
  1264. *
  1265. * Path 6o (M >= N, JOBZ='O')
  1266. * Perform bidiagonal SVD, computing left singular vectors
  1267. * of bidiagonal matrix in RWORK(IRU) and computing right
  1268. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1269. * CWorkspace: need 0
  1270. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  1271. *
  1272. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  1273. $ N, RWORK( IRVT ), N, DUM, IDUM,
  1274. $ RWORK( NRWORK ), IWORK, INFO )
  1275. *
  1276. * Copy real matrix RWORK(IRVT) to complex matrix VT
  1277. * Overwrite VT by right singular vectors of A
  1278. * CWorkspace: need 2*N [tauq, taup] + N*N [U] + N [work]
  1279. * CWorkspace: prefer 2*N [tauq, taup] + N*N [U] + N*NB [work]
  1280. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
  1281. *
  1282. CALL CLACP2( 'F', N, N, RWORK( IRVT ), N, VT, LDVT )
  1283. CALL CUNMBR( 'P', 'R', 'C', N, N, N, A, LDA,
  1284. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1285. $ LWORK-NWORK+1, IERR )
  1286. *
  1287. IF( LWORK .GE. M*N + 3*N ) THEN
  1288. *
  1289. * Path 6o-fast
  1290. * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
  1291. * Overwrite WORK(IU) by left singular vectors of A, copying
  1292. * to A
  1293. * CWorkspace: need 2*N [tauq, taup] + M*N [U] + N [work]
  1294. * CWorkspace: prefer 2*N [tauq, taup] + M*N [U] + N*NB [work]
  1295. * RWorkspace: need N [e] + N*N [RU]
  1296. *
  1297. CALL CLASET( 'F', M, N, CZERO, CZERO, WORK( IU ),
  1298. $ LDWRKU )
  1299. CALL CLACP2( 'F', N, N, RWORK( IRU ), N, WORK( IU ),
  1300. $ LDWRKU )
  1301. CALL CUNMBR( 'Q', 'L', 'N', M, N, N, A, LDA,
  1302. $ WORK( ITAUQ ), WORK( IU ), LDWRKU,
  1303. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1304. CALL CLACPY( 'F', M, N, WORK( IU ), LDWRKU, A, LDA )
  1305. ELSE
  1306. *
  1307. * Path 6o-slow
  1308. * Generate Q in A
  1309. * CWorkspace: need 2*N [tauq, taup] + N*N [U] + N [work]
  1310. * CWorkspace: prefer 2*N [tauq, taup] + N*N [U] + N*NB [work]
  1311. * RWorkspace: need 0
  1312. *
  1313. CALL CUNGBR( 'Q', M, N, N, A, LDA, WORK( ITAUQ ),
  1314. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1315. *
  1316. * Multiply Q in A by real matrix RWORK(IRU), storing the
  1317. * result in WORK(IU), copying to A
  1318. * CWorkspace: need 2*N [tauq, taup] + N*N [U]
  1319. * CWorkspace: prefer 2*N [tauq, taup] + M*N [U]
  1320. * RWorkspace: need N [e] + N*N [RU] + 2*N*N [rwork]
  1321. * RWorkspace: prefer N [e] + N*N [RU] + 2*M*N [rwork] < N + 5*N*N since M < 2*N here
  1322. *
  1323. NRWORK = IRVT
  1324. DO 30 I = 1, M, LDWRKU
  1325. CHUNK = MIN( M-I+1, LDWRKU )
  1326. CALL CLACRM( CHUNK, N, A( I, 1 ), LDA,
  1327. $ RWORK( IRU ), N, WORK( IU ), LDWRKU,
  1328. $ RWORK( NRWORK ) )
  1329. CALL CLACPY( 'F', CHUNK, N, WORK( IU ), LDWRKU,
  1330. $ A( I, 1 ), LDA )
  1331. 30 CONTINUE
  1332. END IF
  1333. *
  1334. ELSE IF( WNTQS ) THEN
  1335. *
  1336. * Path 6s (M >= N, JOBZ='S')
  1337. * Perform bidiagonal SVD, computing left singular vectors
  1338. * of bidiagonal matrix in RWORK(IRU) and computing right
  1339. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1340. * CWorkspace: need 0
  1341. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  1342. *
  1343. IRU = NRWORK
  1344. IRVT = IRU + N*N
  1345. NRWORK = IRVT + N*N
  1346. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  1347. $ N, RWORK( IRVT ), N, DUM, IDUM,
  1348. $ RWORK( NRWORK ), IWORK, INFO )
  1349. *
  1350. * Copy real matrix RWORK(IRU) to complex matrix U
  1351. * Overwrite U by left singular vectors of A
  1352. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1353. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1354. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
  1355. *
  1356. CALL CLASET( 'F', M, N, CZERO, CZERO, U, LDU )
  1357. CALL CLACP2( 'F', N, N, RWORK( IRU ), N, U, LDU )
  1358. CALL CUNMBR( 'Q', 'L', 'N', M, N, N, A, LDA,
  1359. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1360. $ LWORK-NWORK+1, IERR )
  1361. *
  1362. * Copy real matrix RWORK(IRVT) to complex matrix VT
  1363. * Overwrite VT by right singular vectors of A
  1364. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1365. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1366. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
  1367. *
  1368. CALL CLACP2( 'F', N, N, RWORK( IRVT ), N, VT, LDVT )
  1369. CALL CUNMBR( 'P', 'R', 'C', N, N, N, A, LDA,
  1370. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1371. $ LWORK-NWORK+1, IERR )
  1372. ELSE
  1373. *
  1374. * Path 6a (M >= N, JOBZ='A')
  1375. * Perform bidiagonal SVD, computing left singular vectors
  1376. * of bidiagonal matrix in RWORK(IRU) and computing right
  1377. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1378. * CWorkspace: need 0
  1379. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT] + BDSPAC
  1380. *
  1381. IRU = NRWORK
  1382. IRVT = IRU + N*N
  1383. NRWORK = IRVT + N*N
  1384. CALL SBDSDC( 'U', 'I', N, S, RWORK( IE ), RWORK( IRU ),
  1385. $ N, RWORK( IRVT ), N, DUM, IDUM,
  1386. $ RWORK( NRWORK ), IWORK, INFO )
  1387. *
  1388. * Set the right corner of U to identity matrix
  1389. *
  1390. CALL CLASET( 'F', M, M, CZERO, CZERO, U, LDU )
  1391. IF( M.GT.N ) THEN
  1392. CALL CLASET( 'F', M-N, M-N, CZERO, CONE,
  1393. $ U( N+1, N+1 ), LDU )
  1394. END IF
  1395. *
  1396. * Copy real matrix RWORK(IRU) to complex matrix U
  1397. * Overwrite U by left singular vectors of A
  1398. * CWorkspace: need 2*N [tauq, taup] + M [work]
  1399. * CWorkspace: prefer 2*N [tauq, taup] + M*NB [work]
  1400. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
  1401. *
  1402. CALL CLACP2( 'F', N, N, RWORK( IRU ), N, U, LDU )
  1403. CALL CUNMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  1404. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1405. $ LWORK-NWORK+1, IERR )
  1406. *
  1407. * Copy real matrix RWORK(IRVT) to complex matrix VT
  1408. * Overwrite VT by right singular vectors of A
  1409. * CWorkspace: need 2*N [tauq, taup] + N [work]
  1410. * CWorkspace: prefer 2*N [tauq, taup] + N*NB [work]
  1411. * RWorkspace: need N [e] + N*N [RU] + N*N [RVT]
  1412. *
  1413. CALL CLACP2( 'F', N, N, RWORK( IRVT ), N, VT, LDVT )
  1414. CALL CUNMBR( 'P', 'R', 'C', N, N, N, A, LDA,
  1415. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1416. $ LWORK-NWORK+1, IERR )
  1417. END IF
  1418. *
  1419. END IF
  1420. *
  1421. ELSE
  1422. *
  1423. * A has more columns than rows. If A has sufficiently more
  1424. * columns than rows, first reduce using the LQ decomposition (if
  1425. * sufficient workspace available)
  1426. *
  1427. IF( N.GE.MNTHR1 ) THEN
  1428. *
  1429. IF( WNTQN ) THEN
  1430. *
  1431. * Path 1t (N >> M, JOBZ='N')
  1432. * No singular vectors to be computed
  1433. *
  1434. ITAU = 1
  1435. NWORK = ITAU + M
  1436. *
  1437. * Compute A=L*Q
  1438. * CWorkspace: need M [tau] + M [work]
  1439. * CWorkspace: prefer M [tau] + M*NB [work]
  1440. * RWorkspace: need 0
  1441. *
  1442. CALL CGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1443. $ LWORK-NWORK+1, IERR )
  1444. *
  1445. * Zero out above L
  1446. *
  1447. CALL CLASET( 'U', M-1, M-1, CZERO, CZERO, A( 1, 2 ),
  1448. $ LDA )
  1449. IE = 1
  1450. ITAUQ = 1
  1451. ITAUP = ITAUQ + M
  1452. NWORK = ITAUP + M
  1453. *
  1454. * Bidiagonalize L in A
  1455. * CWorkspace: need 2*M [tauq, taup] + M [work]
  1456. * CWorkspace: prefer 2*M [tauq, taup] + 2*M*NB [work]
  1457. * RWorkspace: need M [e]
  1458. *
  1459. CALL CGEBRD( M, M, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  1460. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1461. $ IERR )
  1462. NRWORK = IE + M
  1463. *
  1464. * Perform bidiagonal SVD, compute singular values only
  1465. * CWorkspace: need 0
  1466. * RWorkspace: need M [e] + BDSPAC
  1467. *
  1468. CALL SBDSDC( 'U', 'N', M, S, RWORK( IE ), DUM,1,DUM,1,
  1469. $ DUM, IDUM, RWORK( NRWORK ), IWORK, INFO )
  1470. *
  1471. ELSE IF( WNTQO ) THEN
  1472. *
  1473. * Path 2t (N >> M, JOBZ='O')
  1474. * M right singular vectors to be overwritten on A and
  1475. * M left singular vectors to be computed in U
  1476. *
  1477. IVT = 1
  1478. LDWKVT = M
  1479. *
  1480. * WORK(IVT) is M by M
  1481. *
  1482. IL = IVT + LDWKVT*M
  1483. IF( LWORK .GE. M*N + M*M + 3*M ) THEN
  1484. *
  1485. * WORK(IL) M by N
  1486. *
  1487. LDWRKL = M
  1488. CHUNK = N
  1489. ELSE
  1490. *
  1491. * WORK(IL) is M by CHUNK
  1492. *
  1493. LDWRKL = M
  1494. CHUNK = ( LWORK - M*M - 3*M ) / M
  1495. END IF
  1496. ITAU = IL + LDWRKL*CHUNK
  1497. NWORK = ITAU + M
  1498. *
  1499. * Compute A=L*Q
  1500. * CWorkspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
  1501. * CWorkspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
  1502. * RWorkspace: need 0
  1503. *
  1504. CALL CGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1505. $ LWORK-NWORK+1, IERR )
  1506. *
  1507. * Copy L to WORK(IL), zeroing about above it
  1508. *
  1509. CALL CLACPY( 'L', M, M, A, LDA, WORK( IL ), LDWRKL )
  1510. CALL CLASET( 'U', M-1, M-1, CZERO, CZERO,
  1511. $ WORK( IL+LDWRKL ), LDWRKL )
  1512. *
  1513. * Generate Q in A
  1514. * CWorkspace: need M*M [VT] + M*M [L] + M [tau] + M [work]
  1515. * CWorkspace: prefer M*M [VT] + M*M [L] + M [tau] + M*NB [work]
  1516. * RWorkspace: need 0
  1517. *
  1518. CALL CUNGLQ( M, N, M, A, LDA, WORK( ITAU ),
  1519. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1520. IE = 1
  1521. ITAUQ = ITAU
  1522. ITAUP = ITAUQ + M
  1523. NWORK = ITAUP + M
  1524. *
  1525. * Bidiagonalize L in WORK(IL)
  1526. * CWorkspace: need M*M [VT] + M*M [L] + 2*M [tauq, taup] + M [work]
  1527. * CWorkspace: prefer M*M [VT] + M*M [L] + 2*M [tauq, taup] + 2*M*NB [work]
  1528. * RWorkspace: need M [e]
  1529. *
  1530. CALL CGEBRD( M, M, WORK( IL ), LDWRKL, S, RWORK( IE ),
  1531. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  1532. $ LWORK-NWORK+1, IERR )
  1533. *
  1534. * Perform bidiagonal SVD, computing left singular vectors
  1535. * of bidiagonal matrix in RWORK(IRU) and computing right
  1536. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1537. * CWorkspace: need 0
  1538. * RWorkspace: need M [e] + M*M [RU] + M*M [RVT] + BDSPAC
  1539. *
  1540. IRU = IE + M
  1541. IRVT = IRU + M*M
  1542. NRWORK = IRVT + M*M
  1543. CALL SBDSDC( 'U', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  1544. $ M, RWORK( IRVT ), M, DUM, IDUM,
  1545. $ RWORK( NRWORK ), IWORK, INFO )
  1546. *
  1547. * Copy real matrix RWORK(IRU) to complex matrix WORK(IU)
  1548. * Overwrite WORK(IU) by the left singular vectors of L
  1549. * CWorkspace: need M*M [VT] + M*M [L] + 2*M [tauq, taup] + M [work]
  1550. * CWorkspace: prefer M*M [VT] + M*M [L] + 2*M [tauq, taup] + M*NB [work]
  1551. * RWorkspace: need 0
  1552. *
  1553. CALL CLACP2( 'F', M, M, RWORK( IRU ), M, U, LDU )
  1554. CALL CUNMBR( 'Q', 'L', 'N', M, M, M, WORK( IL ), LDWRKL,
  1555. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1556. $ LWORK-NWORK+1, IERR )
  1557. *
  1558. * Copy real matrix RWORK(IRVT) to complex matrix WORK(IVT)
  1559. * Overwrite WORK(IVT) by the right singular vectors of L
  1560. * CWorkspace: need M*M [VT] + M*M [L] + 2*M [tauq, taup] + M [work]
  1561. * CWorkspace: prefer M*M [VT] + M*M [L] + 2*M [tauq, taup] + M*NB [work]
  1562. * RWorkspace: need 0
  1563. *
  1564. CALL CLACP2( 'F', M, M, RWORK( IRVT ), M, WORK( IVT ),
  1565. $ LDWKVT )
  1566. CALL CUNMBR( 'P', 'R', 'C', M, M, M, WORK( IL ), LDWRKL,
  1567. $ WORK( ITAUP ), WORK( IVT ), LDWKVT,
  1568. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1569. *
  1570. * Multiply right singular vectors of L in WORK(IL) by Q
  1571. * in A, storing result in WORK(IL) and copying to A
  1572. * CWorkspace: need M*M [VT] + M*M [L]
  1573. * CWorkspace: prefer M*M [VT] + M*N [L]
  1574. * RWorkspace: need 0
  1575. *
  1576. DO 40 I = 1, N, CHUNK
  1577. BLK = MIN( N-I+1, CHUNK )
  1578. CALL CGEMM( 'N', 'N', M, BLK, M, CONE, WORK( IVT ), M,
  1579. $ A( 1, I ), LDA, CZERO, WORK( IL ),
  1580. $ LDWRKL )
  1581. CALL CLACPY( 'F', M, BLK, WORK( IL ), LDWRKL,
  1582. $ A( 1, I ), LDA )
  1583. 40 CONTINUE
  1584. *
  1585. ELSE IF( WNTQS ) THEN
  1586. *
  1587. * Path 3t (N >> M, JOBZ='S')
  1588. * M right singular vectors to be computed in VT and
  1589. * M left singular vectors to be computed in U
  1590. *
  1591. IL = 1
  1592. *
  1593. * WORK(IL) is M by M
  1594. *
  1595. LDWRKL = M
  1596. ITAU = IL + LDWRKL*M
  1597. NWORK = ITAU + M
  1598. *
  1599. * Compute A=L*Q
  1600. * CWorkspace: need M*M [L] + M [tau] + M [work]
  1601. * CWorkspace: prefer M*M [L] + M [tau] + M*NB [work]
  1602. * RWorkspace: need 0
  1603. *
  1604. CALL CGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1605. $ LWORK-NWORK+1, IERR )
  1606. *
  1607. * Copy L to WORK(IL), zeroing out above it
  1608. *
  1609. CALL CLACPY( 'L', M, M, A, LDA, WORK( IL ), LDWRKL )
  1610. CALL CLASET( 'U', M-1, M-1, CZERO, CZERO,
  1611. $ WORK( IL+LDWRKL ), LDWRKL )
  1612. *
  1613. * Generate Q in A
  1614. * CWorkspace: need M*M [L] + M [tau] + M [work]
  1615. * CWorkspace: prefer M*M [L] + M [tau] + M*NB [work]
  1616. * RWorkspace: need 0
  1617. *
  1618. CALL CUNGLQ( M, N, M, A, LDA, WORK( ITAU ),
  1619. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1620. IE = 1
  1621. ITAUQ = ITAU
  1622. ITAUP = ITAUQ + M
  1623. NWORK = ITAUP + M
  1624. *
  1625. * Bidiagonalize L in WORK(IL)
  1626. * CWorkspace: need M*M [L] + 2*M [tauq, taup] + M [work]
  1627. * CWorkspace: prefer M*M [L] + 2*M [tauq, taup] + 2*M*NB [work]
  1628. * RWorkspace: need M [e]
  1629. *
  1630. CALL CGEBRD( M, M, WORK( IL ), LDWRKL, S, RWORK( IE ),
  1631. $ WORK( ITAUQ ), WORK( ITAUP ), WORK( NWORK ),
  1632. $ LWORK-NWORK+1, IERR )
  1633. *
  1634. * Perform bidiagonal SVD, computing left singular vectors
  1635. * of bidiagonal matrix in RWORK(IRU) and computing right
  1636. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1637. * CWorkspace: need 0
  1638. * RWorkspace: need M [e] + M*M [RU] + M*M [RVT] + BDSPAC
  1639. *
  1640. IRU = IE + M
  1641. IRVT = IRU + M*M
  1642. NRWORK = IRVT + M*M
  1643. CALL SBDSDC( 'U', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  1644. $ M, RWORK( IRVT ), M, DUM, IDUM,
  1645. $ RWORK( NRWORK ), IWORK, INFO )
  1646. *
  1647. * Copy real matrix RWORK(IRU) to complex matrix U
  1648. * Overwrite U by left singular vectors of L
  1649. * CWorkspace: need M*M [L] + 2*M [tauq, taup] + M [work]
  1650. * CWorkspace: prefer M*M [L] + 2*M [tauq, taup] + M*NB [work]
  1651. * RWorkspace: need 0
  1652. *
  1653. CALL CLACP2( 'F', M, M, RWORK( IRU ), M, U, LDU )
  1654. CALL CUNMBR( 'Q', 'L', 'N', M, M, M, WORK( IL ), LDWRKL,
  1655. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1656. $ LWORK-NWORK+1, IERR )
  1657. *
  1658. * Copy real matrix RWORK(IRVT) to complex matrix VT
  1659. * Overwrite VT by left singular vectors of L
  1660. * CWorkspace: need M*M [L] + 2*M [tauq, taup] + M [work]
  1661. * CWorkspace: prefer M*M [L] + 2*M [tauq, taup] + M*NB [work]
  1662. * RWorkspace: need 0
  1663. *
  1664. CALL CLACP2( 'F', M, M, RWORK( IRVT ), M, VT, LDVT )
  1665. CALL CUNMBR( 'P', 'R', 'C', M, M, M, WORK( IL ), LDWRKL,
  1666. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  1667. $ LWORK-NWORK+1, IERR )
  1668. *
  1669. * Copy VT to WORK(IL), multiply right singular vectors of L
  1670. * in WORK(IL) by Q in A, storing result in VT
  1671. * CWorkspace: need M*M [L]
  1672. * RWorkspace: need 0
  1673. *
  1674. CALL CLACPY( 'F', M, M, VT, LDVT, WORK( IL ), LDWRKL )
  1675. CALL CGEMM( 'N', 'N', M, N, M, CONE, WORK( IL ), LDWRKL,
  1676. $ A, LDA, CZERO, VT, LDVT )
  1677. *
  1678. ELSE IF( WNTQA ) THEN
  1679. *
  1680. * Path 4t (N >> M, JOBZ='A')
  1681. * N right singular vectors to be computed in VT and
  1682. * M left singular vectors to be computed in U
  1683. *
  1684. IVT = 1
  1685. *
  1686. * WORK(IVT) is M by M
  1687. *
  1688. LDWKVT = M
  1689. ITAU = IVT + LDWKVT*M
  1690. NWORK = ITAU + M
  1691. *
  1692. * Compute A=L*Q, copying result to VT
  1693. * CWorkspace: need M*M [VT] + M [tau] + M [work]
  1694. * CWorkspace: prefer M*M [VT] + M [tau] + M*NB [work]
  1695. * RWorkspace: need 0
  1696. *
  1697. CALL CGELQF( M, N, A, LDA, WORK( ITAU ), WORK( NWORK ),
  1698. $ LWORK-NWORK+1, IERR )
  1699. CALL CLACPY( 'U', M, N, A, LDA, VT, LDVT )
  1700. *
  1701. * Generate Q in VT
  1702. * CWorkspace: need M*M [VT] + M [tau] + N [work]
  1703. * CWorkspace: prefer M*M [VT] + M [tau] + N*NB [work]
  1704. * RWorkspace: need 0
  1705. *
  1706. CALL CUNGLQ( N, N, M, VT, LDVT, WORK( ITAU ),
  1707. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1708. *
  1709. * Produce L in A, zeroing out above it
  1710. *
  1711. CALL CLASET( 'U', M-1, M-1, CZERO, CZERO, A( 1, 2 ),
  1712. $ LDA )
  1713. IE = 1
  1714. ITAUQ = ITAU
  1715. ITAUP = ITAUQ + M
  1716. NWORK = ITAUP + M
  1717. *
  1718. * Bidiagonalize L in A
  1719. * CWorkspace: need M*M [VT] + 2*M [tauq, taup] + M [work]
  1720. * CWorkspace: prefer M*M [VT] + 2*M [tauq, taup] + 2*M*NB [work]
  1721. * RWorkspace: need M [e]
  1722. *
  1723. CALL CGEBRD( M, M, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  1724. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1725. $ IERR )
  1726. *
  1727. * Perform bidiagonal SVD, computing left singular vectors
  1728. * of bidiagonal matrix in RWORK(IRU) and computing right
  1729. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1730. * CWorkspace: need 0
  1731. * RWorkspace: need M [e] + M*M [RU] + M*M [RVT] + BDSPAC
  1732. *
  1733. IRU = IE + M
  1734. IRVT = IRU + M*M
  1735. NRWORK = IRVT + M*M
  1736. CALL SBDSDC( 'U', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  1737. $ M, RWORK( IRVT ), M, DUM, IDUM,
  1738. $ RWORK( NRWORK ), IWORK, INFO )
  1739. *
  1740. * Copy real matrix RWORK(IRU) to complex matrix U
  1741. * Overwrite U by left singular vectors of L
  1742. * CWorkspace: need M*M [VT] + 2*M [tauq, taup] + M [work]
  1743. * CWorkspace: prefer M*M [VT] + 2*M [tauq, taup] + M*NB [work]
  1744. * RWorkspace: need 0
  1745. *
  1746. CALL CLACP2( 'F', M, M, RWORK( IRU ), M, U, LDU )
  1747. CALL CUNMBR( 'Q', 'L', 'N', M, M, M, A, LDA,
  1748. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  1749. $ LWORK-NWORK+1, IERR )
  1750. *
  1751. * Copy real matrix RWORK(IRVT) to complex matrix WORK(IVT)
  1752. * Overwrite WORK(IVT) by right singular vectors of L
  1753. * CWorkspace: need M*M [VT] + 2*M [tauq, taup] + M [work]
  1754. * CWorkspace: prefer M*M [VT] + 2*M [tauq, taup] + M*NB [work]
  1755. * RWorkspace: need 0
  1756. *
  1757. CALL CLACP2( 'F', M, M, RWORK( IRVT ), M, WORK( IVT ),
  1758. $ LDWKVT )
  1759. CALL CUNMBR( 'P', 'R', 'C', M, M, M, A, LDA,
  1760. $ WORK( ITAUP ), WORK( IVT ), LDWKVT,
  1761. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1762. *
  1763. * Multiply right singular vectors of L in WORK(IVT) by
  1764. * Q in VT, storing result in A
  1765. * CWorkspace: need M*M [VT]
  1766. * RWorkspace: need 0
  1767. *
  1768. CALL CGEMM( 'N', 'N', M, N, M, CONE, WORK( IVT ), LDWKVT,
  1769. $ VT, LDVT, CZERO, A, LDA )
  1770. *
  1771. * Copy right singular vectors of A from A to VT
  1772. *
  1773. CALL CLACPY( 'F', M, N, A, LDA, VT, LDVT )
  1774. *
  1775. END IF
  1776. *
  1777. ELSE IF( N.GE.MNTHR2 ) THEN
  1778. *
  1779. * MNTHR2 <= N < MNTHR1
  1780. *
  1781. * Path 5t (N >> M, but not as much as MNTHR1)
  1782. * Reduce to bidiagonal form without QR decomposition, use
  1783. * CUNGBR and matrix multiplication to compute singular vectors
  1784. *
  1785. IE = 1
  1786. NRWORK = IE + M
  1787. ITAUQ = 1
  1788. ITAUP = ITAUQ + M
  1789. NWORK = ITAUP + M
  1790. *
  1791. * Bidiagonalize A
  1792. * CWorkspace: need 2*M [tauq, taup] + N [work]
  1793. * CWorkspace: prefer 2*M [tauq, taup] + (M+N)*NB [work]
  1794. * RWorkspace: need M [e]
  1795. *
  1796. CALL CGEBRD( M, N, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  1797. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  1798. $ IERR )
  1799. *
  1800. IF( WNTQN ) THEN
  1801. *
  1802. * Path 5tn (N >> M, JOBZ='N')
  1803. * Compute singular values only
  1804. * CWorkspace: need 0
  1805. * RWorkspace: need M [e] + BDSPAC
  1806. *
  1807. CALL SBDSDC( 'L', 'N', M, S, RWORK( IE ), DUM,1,DUM,1,
  1808. $ DUM, IDUM, RWORK( NRWORK ), IWORK, INFO )
  1809. ELSE IF( WNTQO ) THEN
  1810. IRVT = NRWORK
  1811. IRU = IRVT + M*M
  1812. NRWORK = IRU + M*M
  1813. IVT = NWORK
  1814. *
  1815. * Path 5to (N >> M, JOBZ='O')
  1816. * Copy A to U, generate Q
  1817. * CWorkspace: need 2*M [tauq, taup] + M [work]
  1818. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  1819. * RWorkspace: need 0
  1820. *
  1821. CALL CLACPY( 'L', M, M, A, LDA, U, LDU )
  1822. CALL CUNGBR( 'Q', M, M, N, U, LDU, WORK( ITAUQ ),
  1823. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1824. *
  1825. * Generate P**H in A
  1826. * CWorkspace: need 2*M [tauq, taup] + M [work]
  1827. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  1828. * RWorkspace: need 0
  1829. *
  1830. CALL CUNGBR( 'P', M, N, M, A, LDA, WORK( ITAUP ),
  1831. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1832. *
  1833. LDWKVT = M
  1834. IF( LWORK .GE. M*N + 3*M ) THEN
  1835. *
  1836. * WORK( IVT ) is M by N
  1837. *
  1838. NWORK = IVT + LDWKVT*N
  1839. CHUNK = N
  1840. ELSE
  1841. *
  1842. * WORK( IVT ) is M by CHUNK
  1843. *
  1844. CHUNK = ( LWORK - 3*M ) / M
  1845. NWORK = IVT + LDWKVT*CHUNK
  1846. END IF
  1847. *
  1848. * Perform bidiagonal SVD, computing left singular vectors
  1849. * of bidiagonal matrix in RWORK(IRU) and computing right
  1850. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1851. * CWorkspace: need 0
  1852. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
  1853. *
  1854. CALL SBDSDC( 'L', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  1855. $ M, RWORK( IRVT ), M, DUM, IDUM,
  1856. $ RWORK( NRWORK ), IWORK, INFO )
  1857. *
  1858. * Multiply Q in U by real matrix RWORK(IRVT)
  1859. * storing the result in WORK(IVT), copying to U
  1860. * CWorkspace: need 2*M [tauq, taup] + M*M [VT]
  1861. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + 2*M*M [rwork]
  1862. *
  1863. CALL CLACRM( M, M, U, LDU, RWORK( IRU ), M, WORK( IVT ),
  1864. $ LDWKVT, RWORK( NRWORK ) )
  1865. CALL CLACPY( 'F', M, M, WORK( IVT ), LDWKVT, U, LDU )
  1866. *
  1867. * Multiply RWORK(IRVT) by P**H in A, storing the
  1868. * result in WORK(IVT), copying to A
  1869. * CWorkspace: need 2*M [tauq, taup] + M*M [VT]
  1870. * CWorkspace: prefer 2*M [tauq, taup] + M*N [VT]
  1871. * RWorkspace: need M [e] + M*M [RVT] + 2*M*M [rwork]
  1872. * RWorkspace: prefer M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
  1873. *
  1874. NRWORK = IRU
  1875. DO 50 I = 1, N, CHUNK
  1876. BLK = MIN( N-I+1, CHUNK )
  1877. CALL CLARCM( M, BLK, RWORK( IRVT ), M, A( 1, I ), LDA,
  1878. $ WORK( IVT ), LDWKVT, RWORK( NRWORK ) )
  1879. CALL CLACPY( 'F', M, BLK, WORK( IVT ), LDWKVT,
  1880. $ A( 1, I ), LDA )
  1881. 50 CONTINUE
  1882. ELSE IF( WNTQS ) THEN
  1883. *
  1884. * Path 5ts (N >> M, JOBZ='S')
  1885. * Copy A to U, generate Q
  1886. * CWorkspace: need 2*M [tauq, taup] + M [work]
  1887. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  1888. * RWorkspace: need 0
  1889. *
  1890. CALL CLACPY( 'L', M, M, A, LDA, U, LDU )
  1891. CALL CUNGBR( 'Q', M, M, N, U, LDU, WORK( ITAUQ ),
  1892. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1893. *
  1894. * Copy A to VT, generate P**H
  1895. * CWorkspace: need 2*M [tauq, taup] + M [work]
  1896. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  1897. * RWorkspace: need 0
  1898. *
  1899. CALL CLACPY( 'U', M, N, A, LDA, VT, LDVT )
  1900. CALL CUNGBR( 'P', M, N, M, VT, LDVT, WORK( ITAUP ),
  1901. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1902. *
  1903. * Perform bidiagonal SVD, computing left singular vectors
  1904. * of bidiagonal matrix in RWORK(IRU) and computing right
  1905. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1906. * CWorkspace: need 0
  1907. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
  1908. *
  1909. IRVT = NRWORK
  1910. IRU = IRVT + M*M
  1911. NRWORK = IRU + M*M
  1912. CALL SBDSDC( 'L', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  1913. $ M, RWORK( IRVT ), M, DUM, IDUM,
  1914. $ RWORK( NRWORK ), IWORK, INFO )
  1915. *
  1916. * Multiply Q in U by real matrix RWORK(IRU), storing the
  1917. * result in A, copying to U
  1918. * CWorkspace: need 0
  1919. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + 2*M*M [rwork]
  1920. *
  1921. CALL CLACRM( M, M, U, LDU, RWORK( IRU ), M, A, LDA,
  1922. $ RWORK( NRWORK ) )
  1923. CALL CLACPY( 'F', M, M, A, LDA, U, LDU )
  1924. *
  1925. * Multiply real matrix RWORK(IRVT) by P**H in VT,
  1926. * storing the result in A, copying to VT
  1927. * CWorkspace: need 0
  1928. * RWorkspace: need M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
  1929. *
  1930. NRWORK = IRU
  1931. CALL CLARCM( M, N, RWORK( IRVT ), M, VT, LDVT, A, LDA,
  1932. $ RWORK( NRWORK ) )
  1933. CALL CLACPY( 'F', M, N, A, LDA, VT, LDVT )
  1934. ELSE
  1935. *
  1936. * Path 5ta (N >> M, JOBZ='A')
  1937. * Copy A to U, generate Q
  1938. * CWorkspace: need 2*M [tauq, taup] + M [work]
  1939. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  1940. * RWorkspace: need 0
  1941. *
  1942. CALL CLACPY( 'L', M, M, A, LDA, U, LDU )
  1943. CALL CUNGBR( 'Q', M, M, N, U, LDU, WORK( ITAUQ ),
  1944. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1945. *
  1946. * Copy A to VT, generate P**H
  1947. * CWorkspace: need 2*M [tauq, taup] + N [work]
  1948. * CWorkspace: prefer 2*M [tauq, taup] + N*NB [work]
  1949. * RWorkspace: need 0
  1950. *
  1951. CALL CLACPY( 'U', M, N, A, LDA, VT, LDVT )
  1952. CALL CUNGBR( 'P', N, N, M, VT, LDVT, WORK( ITAUP ),
  1953. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  1954. *
  1955. * Perform bidiagonal SVD, computing left singular vectors
  1956. * of bidiagonal matrix in RWORK(IRU) and computing right
  1957. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  1958. * CWorkspace: need 0
  1959. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
  1960. *
  1961. IRVT = NRWORK
  1962. IRU = IRVT + M*M
  1963. NRWORK = IRU + M*M
  1964. CALL SBDSDC( 'L', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  1965. $ M, RWORK( IRVT ), M, DUM, IDUM,
  1966. $ RWORK( NRWORK ), IWORK, INFO )
  1967. *
  1968. * Multiply Q in U by real matrix RWORK(IRU), storing the
  1969. * result in A, copying to U
  1970. * CWorkspace: need 0
  1971. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + 2*M*M [rwork]
  1972. *
  1973. CALL CLACRM( M, M, U, LDU, RWORK( IRU ), M, A, LDA,
  1974. $ RWORK( NRWORK ) )
  1975. CALL CLACPY( 'F', M, M, A, LDA, U, LDU )
  1976. *
  1977. * Multiply real matrix RWORK(IRVT) by P**H in VT,
  1978. * storing the result in A, copying to VT
  1979. * CWorkspace: need 0
  1980. * RWorkspace: need M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
  1981. *
  1982. NRWORK = IRU
  1983. CALL CLARCM( M, N, RWORK( IRVT ), M, VT, LDVT, A, LDA,
  1984. $ RWORK( NRWORK ) )
  1985. CALL CLACPY( 'F', M, N, A, LDA, VT, LDVT )
  1986. END IF
  1987. *
  1988. ELSE
  1989. *
  1990. * N .LT. MNTHR2
  1991. *
  1992. * Path 6t (N > M, but not much larger)
  1993. * Reduce to bidiagonal form without LQ decomposition
  1994. * Use CUNMBR to compute singular vectors
  1995. *
  1996. IE = 1
  1997. NRWORK = IE + M
  1998. ITAUQ = 1
  1999. ITAUP = ITAUQ + M
  2000. NWORK = ITAUP + M
  2001. *
  2002. * Bidiagonalize A
  2003. * CWorkspace: need 2*M [tauq, taup] + N [work]
  2004. * CWorkspace: prefer 2*M [tauq, taup] + (M+N)*NB [work]
  2005. * RWorkspace: need M [e]
  2006. *
  2007. CALL CGEBRD( M, N, A, LDA, S, RWORK( IE ), WORK( ITAUQ ),
  2008. $ WORK( ITAUP ), WORK( NWORK ), LWORK-NWORK+1,
  2009. $ IERR )
  2010. IF( WNTQN ) THEN
  2011. *
  2012. * Path 6tn (N > M, JOBZ='N')
  2013. * Compute singular values only
  2014. * CWorkspace: need 0
  2015. * RWorkspace: need M [e] + BDSPAC
  2016. *
  2017. CALL SBDSDC( 'L', 'N', M, S, RWORK( IE ), DUM,1,DUM,1,
  2018. $ DUM, IDUM, RWORK( NRWORK ), IWORK, INFO )
  2019. ELSE IF( WNTQO ) THEN
  2020. * Path 6to (N > M, JOBZ='O')
  2021. LDWKVT = M
  2022. IVT = NWORK
  2023. IF( LWORK .GE. M*N + 3*M ) THEN
  2024. *
  2025. * WORK( IVT ) is M by N
  2026. *
  2027. CALL CLASET( 'F', M, N, CZERO, CZERO, WORK( IVT ),
  2028. $ LDWKVT )
  2029. NWORK = IVT + LDWKVT*N
  2030. ELSE
  2031. *
  2032. * WORK( IVT ) is M by CHUNK
  2033. *
  2034. CHUNK = ( LWORK - 3*M ) / M
  2035. NWORK = IVT + LDWKVT*CHUNK
  2036. END IF
  2037. *
  2038. * Perform bidiagonal SVD, computing left singular vectors
  2039. * of bidiagonal matrix in RWORK(IRU) and computing right
  2040. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  2041. * CWorkspace: need 0
  2042. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
  2043. *
  2044. IRVT = NRWORK
  2045. IRU = IRVT + M*M
  2046. NRWORK = IRU + M*M
  2047. CALL SBDSDC( 'L', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  2048. $ M, RWORK( IRVT ), M, DUM, IDUM,
  2049. $ RWORK( NRWORK ), IWORK, INFO )
  2050. *
  2051. * Copy real matrix RWORK(IRU) to complex matrix U
  2052. * Overwrite U by left singular vectors of A
  2053. * CWorkspace: need 2*M [tauq, taup] + M*M [VT] + M [work]
  2054. * CWorkspace: prefer 2*M [tauq, taup] + M*M [VT] + M*NB [work]
  2055. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU]
  2056. *
  2057. CALL CLACP2( 'F', M, M, RWORK( IRU ), M, U, LDU )
  2058. CALL CUNMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  2059. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  2060. $ LWORK-NWORK+1, IERR )
  2061. *
  2062. IF( LWORK .GE. M*N + 3*M ) THEN
  2063. *
  2064. * Path 6to-fast
  2065. * Copy real matrix RWORK(IRVT) to complex matrix WORK(IVT)
  2066. * Overwrite WORK(IVT) by right singular vectors of A,
  2067. * copying to A
  2068. * CWorkspace: need 2*M [tauq, taup] + M*N [VT] + M [work]
  2069. * CWorkspace: prefer 2*M [tauq, taup] + M*N [VT] + M*NB [work]
  2070. * RWorkspace: need M [e] + M*M [RVT]
  2071. *
  2072. CALL CLACP2( 'F', M, M, RWORK( IRVT ), M, WORK( IVT ),
  2073. $ LDWKVT )
  2074. CALL CUNMBR( 'P', 'R', 'C', M, N, M, A, LDA,
  2075. $ WORK( ITAUP ), WORK( IVT ), LDWKVT,
  2076. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  2077. CALL CLACPY( 'F', M, N, WORK( IVT ), LDWKVT, A, LDA )
  2078. ELSE
  2079. *
  2080. * Path 6to-slow
  2081. * Generate P**H in A
  2082. * CWorkspace: need 2*M [tauq, taup] + M*M [VT] + M [work]
  2083. * CWorkspace: prefer 2*M [tauq, taup] + M*M [VT] + M*NB [work]
  2084. * RWorkspace: need 0
  2085. *
  2086. CALL CUNGBR( 'P', M, N, M, A, LDA, WORK( ITAUP ),
  2087. $ WORK( NWORK ), LWORK-NWORK+1, IERR )
  2088. *
  2089. * Multiply Q in A by real matrix RWORK(IRU), storing the
  2090. * result in WORK(IU), copying to A
  2091. * CWorkspace: need 2*M [tauq, taup] + M*M [VT]
  2092. * CWorkspace: prefer 2*M [tauq, taup] + M*N [VT]
  2093. * RWorkspace: need M [e] + M*M [RVT] + 2*M*M [rwork]
  2094. * RWorkspace: prefer M [e] + M*M [RVT] + 2*M*N [rwork] < M + 5*M*M since N < 2*M here
  2095. *
  2096. NRWORK = IRU
  2097. DO 60 I = 1, N, CHUNK
  2098. BLK = MIN( N-I+1, CHUNK )
  2099. CALL CLARCM( M, BLK, RWORK( IRVT ), M, A( 1, I ),
  2100. $ LDA, WORK( IVT ), LDWKVT,
  2101. $ RWORK( NRWORK ) )
  2102. CALL CLACPY( 'F', M, BLK, WORK( IVT ), LDWKVT,
  2103. $ A( 1, I ), LDA )
  2104. 60 CONTINUE
  2105. END IF
  2106. ELSE IF( WNTQS ) THEN
  2107. *
  2108. * Path 6ts (N > M, JOBZ='S')
  2109. * Perform bidiagonal SVD, computing left singular vectors
  2110. * of bidiagonal matrix in RWORK(IRU) and computing right
  2111. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  2112. * CWorkspace: need 0
  2113. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
  2114. *
  2115. IRVT = NRWORK
  2116. IRU = IRVT + M*M
  2117. NRWORK = IRU + M*M
  2118. CALL SBDSDC( 'L', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  2119. $ M, RWORK( IRVT ), M, DUM, IDUM,
  2120. $ RWORK( NRWORK ), IWORK, INFO )
  2121. *
  2122. * Copy real matrix RWORK(IRU) to complex matrix U
  2123. * Overwrite U by left singular vectors of A
  2124. * CWorkspace: need 2*M [tauq, taup] + M [work]
  2125. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  2126. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU]
  2127. *
  2128. CALL CLACP2( 'F', M, M, RWORK( IRU ), M, U, LDU )
  2129. CALL CUNMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  2130. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  2131. $ LWORK-NWORK+1, IERR )
  2132. *
  2133. * Copy real matrix RWORK(IRVT) to complex matrix VT
  2134. * Overwrite VT by right singular vectors of A
  2135. * CWorkspace: need 2*M [tauq, taup] + M [work]
  2136. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  2137. * RWorkspace: need M [e] + M*M [RVT]
  2138. *
  2139. CALL CLASET( 'F', M, N, CZERO, CZERO, VT, LDVT )
  2140. CALL CLACP2( 'F', M, M, RWORK( IRVT ), M, VT, LDVT )
  2141. CALL CUNMBR( 'P', 'R', 'C', M, N, M, A, LDA,
  2142. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  2143. $ LWORK-NWORK+1, IERR )
  2144. ELSE
  2145. *
  2146. * Path 6ta (N > M, JOBZ='A')
  2147. * Perform bidiagonal SVD, computing left singular vectors
  2148. * of bidiagonal matrix in RWORK(IRU) and computing right
  2149. * singular vectors of bidiagonal matrix in RWORK(IRVT)
  2150. * CWorkspace: need 0
  2151. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU] + BDSPAC
  2152. *
  2153. IRVT = NRWORK
  2154. IRU = IRVT + M*M
  2155. NRWORK = IRU + M*M
  2156. *
  2157. CALL SBDSDC( 'L', 'I', M, S, RWORK( IE ), RWORK( IRU ),
  2158. $ M, RWORK( IRVT ), M, DUM, IDUM,
  2159. $ RWORK( NRWORK ), IWORK, INFO )
  2160. *
  2161. * Copy real matrix RWORK(IRU) to complex matrix U
  2162. * Overwrite U by left singular vectors of A
  2163. * CWorkspace: need 2*M [tauq, taup] + M [work]
  2164. * CWorkspace: prefer 2*M [tauq, taup] + M*NB [work]
  2165. * RWorkspace: need M [e] + M*M [RVT] + M*M [RU]
  2166. *
  2167. CALL CLACP2( 'F', M, M, RWORK( IRU ), M, U, LDU )
  2168. CALL CUNMBR( 'Q', 'L', 'N', M, M, N, A, LDA,
  2169. $ WORK( ITAUQ ), U, LDU, WORK( NWORK ),
  2170. $ LWORK-NWORK+1, IERR )
  2171. *
  2172. * Set all of VT to identity matrix
  2173. *
  2174. CALL CLASET( 'F', N, N, CZERO, CONE, VT, LDVT )
  2175. *
  2176. * Copy real matrix RWORK(IRVT) to complex matrix VT
  2177. * Overwrite VT by right singular vectors of A
  2178. * CWorkspace: need 2*M [tauq, taup] + N [work]
  2179. * CWorkspace: prefer 2*M [tauq, taup] + N*NB [work]
  2180. * RWorkspace: need M [e] + M*M [RVT]
  2181. *
  2182. CALL CLACP2( 'F', M, M, RWORK( IRVT ), M, VT, LDVT )
  2183. CALL CUNMBR( 'P', 'R', 'C', N, N, M, A, LDA,
  2184. $ WORK( ITAUP ), VT, LDVT, WORK( NWORK ),
  2185. $ LWORK-NWORK+1, IERR )
  2186. END IF
  2187. *
  2188. END IF
  2189. *
  2190. END IF
  2191. *
  2192. * Undo scaling if necessary
  2193. *
  2194. IF( ISCL.EQ.1 ) THEN
  2195. IF( ANRM.GT.BIGNUM )
  2196. $ CALL SLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN, 1, S, MINMN,
  2197. $ IERR )
  2198. IF( INFO.NE.0 .AND. ANRM.GT.BIGNUM )
  2199. $ CALL SLASCL( 'G', 0, 0, BIGNUM, ANRM, MINMN-1, 1,
  2200. $ RWORK( IE ), MINMN, IERR )
  2201. IF( ANRM.LT.SMLNUM )
  2202. $ CALL SLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN, 1, S, MINMN,
  2203. $ IERR )
  2204. IF( INFO.NE.0 .AND. ANRM.LT.SMLNUM )
  2205. $ CALL SLASCL( 'G', 0, 0, SMLNUM, ANRM, MINMN-1, 1,
  2206. $ RWORK( IE ), MINMN, IERR )
  2207. END IF
  2208. *
  2209. * Return optimal workspace in WORK(1)
  2210. *
  2211. WORK( 1 ) = SROUNDUP_LWORK( MAXWRK )
  2212. *
  2213. RETURN
  2214. *
  2215. * End of CGESDD
  2216. *
  2217. END