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.

zgesdd.f 90 kB

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