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.

cgejsv.f 95 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217221822192220222122222223222422252226222722282229223022312232
  1. *> \brief \b CGEJSV
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download CGEJSV + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgejsv.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgejsv.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgejsv.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE CGEJSV( JOBA, JOBU, JOBV, JOBR, JOBT, JOBP,
  22. * M, N, A, LDA, SVA, U, LDU, V, LDV,
  23. * CWORK, LWORK, RWORK, LRWORK, IWORK, INFO )
  24. *
  25. * .. Scalar Arguments ..
  26. * IMPLICIT NONE
  27. * INTEGER INFO, LDA, LDU, LDV, LWORK, M, N
  28. * ..
  29. * .. Array Arguments ..
  30. * COMPLEX A( LDA, * ), U( LDU, * ), V( LDV, * ), CWORK( LWORK )
  31. * REAL SVA( N ), RWORK( LRWORK )
  32. * INTEGER IWORK( * )
  33. * CHARACTER*1 JOBA, JOBP, JOBR, JOBT, JOBU, JOBV
  34. * ..
  35. *
  36. *
  37. *> \par Purpose:
  38. * =============
  39. *>
  40. *> \verbatim
  41. *>
  42. *> CGEJSV computes the singular value decomposition (SVD) of a complex M-by-N
  43. *> matrix [A], where M >= N. The SVD of [A] is written as
  44. *>
  45. *> [A] = [U] * [SIGMA] * [V]^*,
  46. *>
  47. *> where [SIGMA] is an N-by-N (M-by-N) matrix which is zero except for its N
  48. *> diagonal elements, [U] is an M-by-N (or M-by-M) unitary matrix, and
  49. *> [V] is an N-by-N unitary matrix. The diagonal elements of [SIGMA] are
  50. *> the singular values of [A]. The columns of [U] and [V] are the left and
  51. *> the right singular vectors of [A], respectively. The matrices [U] and [V]
  52. *> are computed and stored in the arrays U and V, respectively. The diagonal
  53. *> of [SIGMA] is computed and stored in the array SVA.
  54. *> \endverbatim
  55. *>
  56. *> Arguments:
  57. *> ==========
  58. *>
  59. *> \param[in] JOBA
  60. *> \verbatim
  61. *> JOBA is CHARACTER*1
  62. *> Specifies the level of accuracy:
  63. *> = 'C': This option works well (high relative accuracy) if A = B * D,
  64. *> with well-conditioned B and arbitrary diagonal matrix D.
  65. *> The accuracy cannot be spoiled by COLUMN scaling. The
  66. *> accuracy of the computed output depends on the condition of
  67. *> B, and the procedure aims at the best theoretical accuracy.
  68. *> The relative error max_{i=1:N}|d sigma_i| / sigma_i is
  69. *> bounded by f(M,N)*epsilon* cond(B), independent of D.
  70. *> The input matrix is preprocessed with the QRF with column
  71. *> pivoting. This initial preprocessing and preconditioning by
  72. *> a rank revealing QR factorization is common for all values of
  73. *> JOBA. Additional actions are specified as follows:
  74. *> = 'E': Computation as with 'C' with an additional estimate of the
  75. *> condition number of B. It provides a realistic error bound.
  76. *> = 'F': If A = D1 * C * D2 with ill-conditioned diagonal scalings
  77. *> D1, D2, and well-conditioned matrix C, this option gives
  78. *> higher accuracy than the 'C' option. If the structure of the
  79. *> input matrix is not known, and relative accuracy is
  80. *> desirable, then this option is advisable. The input matrix A
  81. *> is preprocessed with QR factorization with FULL (row and
  82. *> column) pivoting.
  83. *> = 'G': Computation as with 'F' with an additional estimate of the
  84. *> condition number of B, where A=B*D. If A has heavily weighted
  85. *> rows, then using this condition number gives too pessimistic
  86. *> error bound.
  87. *> = 'A': Small singular values are not well determined by the data
  88. *> and are considered as noisy; the matrix is treated as
  89. *> numerically rank deficient. The error in the computed
  90. *> singular values is bounded by f(m,n)*epsilon*||A||.
  91. *> The computed SVD A = U * S * V^* restores A up to
  92. *> f(m,n)*epsilon*||A||.
  93. *> This gives the procedure the licence to discard (set to zero)
  94. *> all singular values below N*epsilon*||A||.
  95. *> = 'R': Similar as in 'A'. Rank revealing property of the initial
  96. *> QR factorization is used do reveal (using triangular factor)
  97. *> a gap sigma_{r+1} < epsilon * sigma_r in which case the
  98. *> numerical RANK is declared to be r. The SVD is computed with
  99. *> absolute error bounds, but more accurately than with 'A'.
  100. *> \endverbatim
  101. *>
  102. *> \param[in] JOBU
  103. *> \verbatim
  104. *> JOBU is CHARACTER*1
  105. *> Specifies whether to compute the columns of U:
  106. *> = 'U': N columns of U are returned in the array U.
  107. *> = 'F': full set of M left sing. vectors is returned in the array U.
  108. *> = 'W': U may be used as workspace of length M*N. See the description
  109. *> of U.
  110. *> = 'N': U is not computed.
  111. *> \endverbatim
  112. *>
  113. *> \param[in] JOBV
  114. *> \verbatim
  115. *> JOBV is CHARACTER*1
  116. *> Specifies whether to compute the matrix V:
  117. *> = 'V': N columns of V are returned in the array V; Jacobi rotations
  118. *> are not explicitly accumulated.
  119. *> = 'J': N columns of V are returned in the array V, but they are
  120. *> computed as the product of Jacobi rotations, if JOBT = 'N'.
  121. *> = 'W': V may be used as workspace of length N*N. See the description
  122. *> of V.
  123. *> = 'N': V is not computed.
  124. *> \endverbatim
  125. *>
  126. *> \param[in] JOBR
  127. *> \verbatim
  128. *> JOBR is CHARACTER*1
  129. *> Specifies the RANGE for the singular values. Issues the licence to
  130. *> set to zero small positive singular values if they are outside
  131. *> specified range. If A .NE. 0 is scaled so that the largest singular
  132. *> value of c*A is around SQRT(BIG), BIG=SLAMCH('O'), then JOBR issues
  133. *> the licence to kill columns of A whose norm in c*A is less than
  134. *> SQRT(SFMIN) (for JOBR = 'R'), or less than SMALL=SFMIN/EPSLN,
  135. *> where SFMIN=SLAMCH('S'), EPSLN=SLAMCH('E').
  136. *> = 'N': Do not kill small columns of c*A. This option assumes that
  137. *> BLAS and QR factorizations and triangular solvers are
  138. *> implemented to work in that range. If the condition of A
  139. *> is greater than BIG, use CGESVJ.
  140. *> = 'R': RESTRICTED range for sigma(c*A) is [SQRT(SFMIN), SQRT(BIG)]
  141. *> (roughly, as described above). This option is recommended.
  142. *> ===========================
  143. *> For computing the singular values in the FULL range [SFMIN,BIG]
  144. *> use CGESVJ.
  145. *> \endverbatim
  146. *>
  147. *> \param[in] JOBT
  148. *> \verbatim
  149. *> JOBT is CHARACTER*1
  150. *> If the matrix is square then the procedure may determine to use
  151. *> transposed A if A^* seems to be better with respect to convergence.
  152. *> If the matrix is not square, JOBT is ignored.
  153. *> The decision is based on two values of entropy over the adjoint
  154. *> orbit of A^* * A. See the descriptions of RWORK(6) and RWORK(7).
  155. *> = 'T': transpose if entropy test indicates possibly faster
  156. *> convergence of Jacobi process if A^* is taken as input. If A is
  157. *> replaced with A^*, then the row pivoting is included automatically.
  158. *> = 'N': do not speculate.
  159. *> The option 'T' can be used to compute only the singular values, or
  160. *> the full SVD (U, SIGMA and V). For only one set of singular vectors
  161. *> (U or V), the caller should provide both U and V, as one of the
  162. *> matrices is used as workspace if the matrix A is transposed.
  163. *> The implementer can easily remove this constraint and make the
  164. *> code more complicated. See the descriptions of U and V.
  165. *> In general, this option is considered experimental, and 'N'; should
  166. *> be preferred. This is subject to changes in the future.
  167. *> \endverbatim
  168. *>
  169. *> \param[in] JOBP
  170. *> \verbatim
  171. *> JOBP is CHARACTER*1
  172. *> Issues the licence to introduce structured perturbations to drown
  173. *> denormalized numbers. This licence should be active if the
  174. *> denormals are poorly implemented, causing slow computation,
  175. *> especially in cases of fast convergence (!). For details see [1,2].
  176. *> For the sake of simplicity, this perturbations are included only
  177. *> when the full SVD or only the singular values are requested. The
  178. *> implementer/user can easily add the perturbation for the cases of
  179. *> computing one set of singular vectors.
  180. *> = 'P': introduce perturbation
  181. *> = 'N': do not perturb
  182. *> \endverbatim
  183. *>
  184. *> \param[in] M
  185. *> \verbatim
  186. *> M is INTEGER
  187. *> The number of rows of the input matrix A. M >= 0.
  188. *> \endverbatim
  189. *>
  190. *> \param[in] N
  191. *> \verbatim
  192. *> N is INTEGER
  193. *> The number of columns of the input matrix A. M >= N >= 0.
  194. *> \endverbatim
  195. *>
  196. *> \param[in,out] A
  197. *> \verbatim
  198. *> A is COMPLEX array, dimension (LDA,N)
  199. *> On entry, the M-by-N matrix A.
  200. *> \endverbatim
  201. *>
  202. *> \param[in] LDA
  203. *> \verbatim
  204. *> LDA is INTEGER
  205. *> The leading dimension of the array A. LDA >= max(1,M).
  206. *> \endverbatim
  207. *>
  208. *> \param[out] SVA
  209. *> \verbatim
  210. *> SVA is REAL array, dimension (N)
  211. *> On exit,
  212. *> - For RWORK(1)/RWORK(2) = ONE: The singular values of A. During
  213. *> the computation SVA contains Euclidean column norms of the
  214. *> iterated matrices in the array A.
  215. *> - For RWORK(1) .NE. RWORK(2): The singular values of A are
  216. *> (RWORK(1)/RWORK(2)) * SVA(1:N). This factored form is used if
  217. *> sigma_max(A) overflows or if small singular values have been
  218. *> saved from underflow by scaling the input matrix A.
  219. *> - If JOBR='R' then some of the singular values may be returned
  220. *> as exact zeros obtained by "set to zero" because they are
  221. *> below the numerical rank threshold or are denormalized numbers.
  222. *> \endverbatim
  223. *>
  224. *> \param[out] U
  225. *> \verbatim
  226. *> U is COMPLEX array, dimension ( LDU, N ) or ( LDU, M )
  227. *> If JOBU = 'U', then U contains on exit the M-by-N matrix of
  228. *> the left singular vectors.
  229. *> If JOBU = 'F', then U contains on exit the M-by-M matrix of
  230. *> the left singular vectors, including an ONB
  231. *> of the orthogonal complement of the Range(A).
  232. *> If JOBU = 'W' .AND. (JOBV = 'V' .AND. JOBT = 'T' .AND. M = N),
  233. *> then U is used as workspace if the procedure
  234. *> replaces A with A^*. In that case, [V] is computed
  235. *> in U as left singular vectors of A^* and then
  236. *> copied back to the V array. This 'W' option is just
  237. *> a reminder to the caller that in this case U is
  238. *> reserved as workspace of length N*N.
  239. *> If JOBU = 'N' U is not referenced, unless JOBT='T'.
  240. *> \endverbatim
  241. *>
  242. *> \param[in] LDU
  243. *> \verbatim
  244. *> LDU is INTEGER
  245. *> The leading dimension of the array U, LDU >= 1.
  246. *> IF JOBU = 'U' or 'F' or 'W', then LDU >= M.
  247. *> \endverbatim
  248. *>
  249. *> \param[out] V
  250. *> \verbatim
  251. *> V is COMPLEX array, dimension ( LDV, N )
  252. *> If JOBV = 'V', 'J' then V contains on exit the N-by-N matrix of
  253. *> the right singular vectors;
  254. *> If JOBV = 'W', AND (JOBU = 'U' AND JOBT = 'T' AND M = N),
  255. *> then V is used as workspace if the pprocedure
  256. *> replaces A with A^*. In that case, [U] is computed
  257. *> in V as right singular vectors of A^* and then
  258. *> copied back to the U array. This 'W' option is just
  259. *> a reminder to the caller that in this case V is
  260. *> reserved as workspace of length N*N.
  261. *> If JOBV = 'N' V is not referenced, unless JOBT='T'.
  262. *> \endverbatim
  263. *>
  264. *> \param[in] LDV
  265. *> \verbatim
  266. *> LDV is INTEGER
  267. *> The leading dimension of the array V, LDV >= 1.
  268. *> If JOBV = 'V' or 'J' or 'W', then LDV >= N.
  269. *> \endverbatim
  270. *>
  271. *> \param[out] CWORK
  272. *> \verbatim
  273. *> CWORK is COMPLEX array, dimension (MAX(2,LWORK))
  274. *> If the call to CGEJSV is a workspace query (indicated by LWORK=-1 or
  275. *> LRWORK=-1), then on exit CWORK(1) contains the required length of
  276. *> CWORK for the job parameters used in the call.
  277. *> \endverbatim
  278. *>
  279. *> \param[in] LWORK
  280. *> \verbatim
  281. *> LWORK is INTEGER
  282. *> Length of CWORK to confirm proper allocation of workspace.
  283. *> LWORK depends on the job:
  284. *>
  285. *> 1. If only SIGMA is needed ( JOBU = 'N', JOBV = 'N' ) and
  286. *> 1.1 .. no scaled condition estimate required (JOBA.NE.'E'.AND.JOBA.NE.'G'):
  287. *> LWORK >= 2*N+1. This is the minimal requirement.
  288. *> ->> For optimal performance (blocked code) the optimal value
  289. *> is LWORK >= N + (N+1)*NB. Here NB is the optimal
  290. *> block size for CGEQP3 and CGEQRF.
  291. *> In general, optimal LWORK is computed as
  292. *> LWORK >= max(N+LWORK(CGEQP3),N+LWORK(CGEQRF), LWORK(CGESVJ)).
  293. *> 1.2. .. an estimate of the scaled condition number of A is
  294. *> required (JOBA='E', or 'G'). In this case, LWORK the minimal
  295. *> requirement is LWORK >= N*N + 2*N.
  296. *> ->> For optimal performance (blocked code) the optimal value
  297. *> is LWORK >= max(N+(N+1)*NB, N*N+2*N)=N**2+2*N.
  298. *> In general, the optimal length LWORK is computed as
  299. *> LWORK >= max(N+LWORK(CGEQP3),N+LWORK(CGEQRF), LWORK(CGESVJ),
  300. *> N*N+LWORK(CPOCON)).
  301. *> 2. If SIGMA and the right singular vectors are needed (JOBV = 'V'),
  302. *> (JOBU = 'N')
  303. *> 2.1 .. no scaled condition estimate requested (JOBE = 'N'):
  304. *> -> the minimal requirement is LWORK >= 3*N.
  305. *> -> For optimal performance,
  306. *> LWORK >= max(N+(N+1)*NB, 2*N+N*NB)=2*N+N*NB,
  307. *> where NB is the optimal block size for CGEQP3, CGEQRF, CGELQF,
  308. *> CUNMLQ. In general, the optimal length LWORK is computed as
  309. *> LWORK >= max(N+LWORK(CGEQP3), N+LWORK(CGESVJ),
  310. *> N+LWORK(CGELQF), 2*N+LWORK(CGEQRF), N+LWORK(CUNMLQ)).
  311. *> 2.2 .. an estimate of the scaled condition number of A is
  312. *> required (JOBA='E', or 'G').
  313. *> -> the minimal requirement is LWORK >= 3*N.
  314. *> -> For optimal performance,
  315. *> LWORK >= max(N+(N+1)*NB, 2*N,2*N+N*NB)=2*N+N*NB,
  316. *> where NB is the optimal block size for CGEQP3, CGEQRF, CGELQF,
  317. *> CUNMLQ. In general, the optimal length LWORK is computed as
  318. *> LWORK >= max(N+LWORK(CGEQP3), LWORK(CPOCON), N+LWORK(CGESVJ),
  319. *> N+LWORK(CGELQF), 2*N+LWORK(CGEQRF), N+LWORK(CUNMLQ)).
  320. *> 3. If SIGMA and the left singular vectors are needed
  321. *> 3.1 .. no scaled condition estimate requested (JOBE = 'N'):
  322. *> -> the minimal requirement is LWORK >= 3*N.
  323. *> -> For optimal performance:
  324. *> if JOBU = 'U' :: LWORK >= max(3*N, N+(N+1)*NB, 2*N+N*NB)=2*N+N*NB,
  325. *> where NB is the optimal block size for CGEQP3, CGEQRF, CUNMQR.
  326. *> In general, the optimal length LWORK is computed as
  327. *> LWORK >= max(N+LWORK(CGEQP3), 2*N+LWORK(CGEQRF), N+LWORK(CUNMQR)).
  328. *> 3.2 .. an estimate of the scaled condition number of A is
  329. *> required (JOBA='E', or 'G').
  330. *> -> the minimal requirement is LWORK >= 3*N.
  331. *> -> For optimal performance:
  332. *> if JOBU = 'U' :: LWORK >= max(3*N, N+(N+1)*NB, 2*N+N*NB)=2*N+N*NB,
  333. *> where NB is the optimal block size for CGEQP3, CGEQRF, CUNMQR.
  334. *> In general, the optimal length LWORK is computed as
  335. *> LWORK >= max(N+LWORK(CGEQP3),N+LWORK(CPOCON),
  336. *> 2*N+LWORK(CGEQRF), N+LWORK(CUNMQR)).
  337. *>
  338. *> 4. If the full SVD is needed: (JOBU = 'U' or JOBU = 'F') and
  339. *> 4.1. if JOBV = 'V'
  340. *> the minimal requirement is LWORK >= 5*N+2*N*N.
  341. *> 4.2. if JOBV = 'J' the minimal requirement is
  342. *> LWORK >= 4*N+N*N.
  343. *> In both cases, the allocated CWORK can accommodate blocked runs
  344. *> of CGEQP3, CGEQRF, CGELQF, CUNMQR, CUNMLQ.
  345. *>
  346. *> If the call to CGEJSV is a workspace query (indicated by LWORK=-1 or
  347. *> LRWORK=-1), then on exit CWORK(1) contains the optimal and CWORK(2) contains the
  348. *> minimal length of CWORK for the job parameters used in the call.
  349. *> \endverbatim
  350. *>
  351. *> \param[out] RWORK
  352. *> \verbatim
  353. *> RWORK is REAL array, dimension (MAX(7,LRWORK))
  354. *> On exit,
  355. *> RWORK(1) = Determines the scaling factor SCALE = RWORK(2) / RWORK(1)
  356. *> such that SCALE*SVA(1:N) are the computed singular values
  357. *> of A. (See the description of SVA().)
  358. *> RWORK(2) = See the description of RWORK(1).
  359. *> RWORK(3) = SCONDA is an estimate for the condition number of
  360. *> column equilibrated A. (If JOBA = 'E' or 'G')
  361. *> SCONDA is an estimate of SQRT(||(R^* * R)^(-1)||_1).
  362. *> It is computed using CPOCON. It holds
  363. *> N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
  364. *> where R is the triangular factor from the QRF of A.
  365. *> However, if R is truncated and the numerical rank is
  366. *> determined to be strictly smaller than N, SCONDA is
  367. *> returned as -1, thus indicating that the smallest
  368. *> singular values might be lost.
  369. *>
  370. *> If full SVD is needed, the following two condition numbers are
  371. *> useful for the analysis of the algorithm. They are provided for
  372. *> a developer/implementer who is familiar with the details of
  373. *> the method.
  374. *>
  375. *> RWORK(4) = an estimate of the scaled condition number of the
  376. *> triangular factor in the first QR factorization.
  377. *> RWORK(5) = an estimate of the scaled condition number of the
  378. *> triangular factor in the second QR factorization.
  379. *> The following two parameters are computed if JOBT = 'T'.
  380. *> They are provided for a developer/implementer who is familiar
  381. *> with the details of the method.
  382. *> RWORK(6) = the entropy of A^* * A :: this is the Shannon entropy
  383. *> of diag(A^* * A) / Trace(A^* * A) taken as point in the
  384. *> probability simplex.
  385. *> RWORK(7) = the entropy of A * A^*. (See the description of RWORK(6).)
  386. *> If the call to CGEJSV is a workspace query (indicated by LWORK=-1 or
  387. *> LRWORK=-1), then on exit RWORK(1) contains the required length of
  388. *> RWORK for the job parameters used in the call.
  389. *> \endverbatim
  390. *>
  391. *> \param[in] LRWORK
  392. *> \verbatim
  393. *> LRWORK is INTEGER
  394. *> Length of RWORK to confirm proper allocation of workspace.
  395. *> LRWORK depends on the job:
  396. *>
  397. *> 1. If only the singular values are requested i.e. if
  398. *> LSAME(JOBU,'N') .AND. LSAME(JOBV,'N')
  399. *> then:
  400. *> 1.1. If LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G'),
  401. *> then: LRWORK = max( 7, 2 * M ).
  402. *> 1.2. Otherwise, LRWORK = max( 7, N ).
  403. *> 2. If singular values with the right singular vectors are requested
  404. *> i.e. if
  405. *> (LSAME(JOBV,'V').OR.LSAME(JOBV,'J')) .AND.
  406. *> .NOT.(LSAME(JOBU,'U').OR.LSAME(JOBU,'F'))
  407. *> then:
  408. *> 2.1. If LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G'),
  409. *> then LRWORK = max( 7, 2 * M ).
  410. *> 2.2. Otherwise, LRWORK = max( 7, N ).
  411. *> 3. If singular values with the left singular vectors are requested, i.e. if
  412. *> (LSAME(JOBU,'U').OR.LSAME(JOBU,'F')) .AND.
  413. *> .NOT.(LSAME(JOBV,'V').OR.LSAME(JOBV,'J'))
  414. *> then:
  415. *> 3.1. If LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G'),
  416. *> then LRWORK = max( 7, 2 * M ).
  417. *> 3.2. Otherwise, LRWORK = max( 7, N ).
  418. *> 4. If singular values with both the left and the right singular vectors
  419. *> are requested, i.e. if
  420. *> (LSAME(JOBU,'U').OR.LSAME(JOBU,'F')) .AND.
  421. *> (LSAME(JOBV,'V').OR.LSAME(JOBV,'J'))
  422. *> then:
  423. *> 4.1. If LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G'),
  424. *> then LRWORK = max( 7, 2 * M ).
  425. *> 4.2. Otherwise, LRWORK = max( 7, N ).
  426. *>
  427. *> If, on entry, LRWORK = -1 or LWORK=-1, a workspace query is assumed and
  428. *> the length of RWORK is returned in RWORK(1).
  429. *> \endverbatim
  430. *>
  431. *> \param[out] IWORK
  432. *> \verbatim
  433. *> IWORK is INTEGER array, of dimension at least 4, that further depends
  434. *> on the job:
  435. *>
  436. *> 1. If only the singular values are requested then:
  437. *> If ( LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G') )
  438. *> then the length of IWORK is N+M; otherwise the length of IWORK is N.
  439. *> 2. If the singular values and the right singular vectors are requested then:
  440. *> If ( LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G') )
  441. *> then the length of IWORK is N+M; otherwise the length of IWORK is N.
  442. *> 3. If the singular values and the left singular vectors are requested then:
  443. *> If ( LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G') )
  444. *> then the length of IWORK is N+M; otherwise the length of IWORK is N.
  445. *> 4. If the singular values with both the left and the right singular vectors
  446. *> are requested, then:
  447. *> 4.1. If LSAME(JOBV,'J') the length of IWORK is determined as follows:
  448. *> If ( LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G') )
  449. *> then the length of IWORK is N+M; otherwise the length of IWORK is N.
  450. *> 4.2. If LSAME(JOBV,'V') the length of IWORK is determined as follows:
  451. *> If ( LSAME(JOBT,'T') .OR. LSAME(JOBA,'F') .OR. LSAME(JOBA,'G') )
  452. *> then the length of IWORK is 2*N+M; otherwise the length of IWORK is 2*N.
  453. *>
  454. *> On exit,
  455. *> IWORK(1) = the numerical rank determined after the initial
  456. *> QR factorization with pivoting. See the descriptions
  457. *> of JOBA and JOBR.
  458. *> IWORK(2) = the number of the computed nonzero singular values
  459. *> IWORK(3) = if nonzero, a warning message:
  460. *> If IWORK(3) = 1 then some of the column norms of A
  461. *> were denormalized floats. The requested high accuracy
  462. *> is not warranted by the data.
  463. *> IWORK(4) = 1 or -1. If IWORK(4) = 1, then the procedure used A^* to
  464. *> do the job as specified by the JOB parameters.
  465. *> If the call to CGEJSV is a workspace query (indicated by LWORK = -1 and
  466. *> LRWORK = -1), then on exit IWORK(1) contains the required length of
  467. *> IWORK for the job parameters used in the call.
  468. *> \endverbatim
  469. *>
  470. *> \param[out] INFO
  471. *> \verbatim
  472. *> INFO is INTEGER
  473. *> < 0: if INFO = -i, then the i-th argument had an illegal value.
  474. *> = 0: successful exit;
  475. *> > 0: CGEJSV did not converge in the maximal allowed number
  476. *> of sweeps. The computed values may be inaccurate.
  477. *> \endverbatim
  478. *
  479. * Authors:
  480. * ========
  481. *
  482. *> \author Univ. of Tennessee
  483. *> \author Univ. of California Berkeley
  484. *> \author Univ. of Colorado Denver
  485. *> \author NAG Ltd.
  486. *
  487. *> \ingroup complexGEsing
  488. *
  489. *> \par Further Details:
  490. * =====================
  491. *>
  492. *> \verbatim
  493. *> CGEJSV implements a preconditioned Jacobi SVD algorithm. It uses CGEQP3,
  494. *> CGEQRF, and CGELQF as preprocessors and preconditioners. Optionally, an
  495. *> additional row pivoting can be used as a preprocessor, which in some
  496. *> cases results in much higher accuracy. An example is matrix A with the
  497. *> structure A = D1 * C * D2, where D1, D2 are arbitrarily ill-conditioned
  498. *> diagonal matrices and C is well-conditioned matrix. In that case, complete
  499. *> pivoting in the first QR factorizations provides accuracy dependent on the
  500. *> condition number of C, and independent of D1, D2. Such higher accuracy is
  501. *> not completely understood theoretically, but it works well in practice.
  502. *> Further, if A can be written as A = B*D, with well-conditioned B and some
  503. *> diagonal D, then the high accuracy is guaranteed, both theoretically and
  504. *> in software, independent of D. For more details see [1], [2].
  505. *> The computational range for the singular values can be the full range
  506. *> ( UNDERFLOW,OVERFLOW ), provided that the machine arithmetic and the BLAS
  507. *> & LAPACK routines called by CGEJSV are implemented to work in that range.
  508. *> If that is not the case, then the restriction for safe computation with
  509. *> the singular values in the range of normalized IEEE numbers is that the
  510. *> spectral condition number kappa(A)=sigma_max(A)/sigma_min(A) does not
  511. *> overflow. This code (CGEJSV) is best used in this restricted range,
  512. *> meaning that singular values of magnitude below ||A||_2 / SLAMCH('O') are
  513. *> returned as zeros. See JOBR for details on this.
  514. *> Further, this implementation is somewhat slower than the one described
  515. *> in [1,2] due to replacement of some non-LAPACK components, and because
  516. *> the choice of some tuning parameters in the iterative part (CGESVJ) is
  517. *> left to the implementer on a particular machine.
  518. *> The rank revealing QR factorization (in this code: CGEQP3) should be
  519. *> implemented as in [3]. We have a new version of CGEQP3 under development
  520. *> that is more robust than the current one in LAPACK, with a cleaner cut in
  521. *> rank deficient cases. It will be available in the SIGMA library [4].
  522. *> If M is much larger than N, it is obvious that the initial QRF with
  523. *> column pivoting can be preprocessed by the QRF without pivoting. That
  524. *> well known trick is not used in CGEJSV because in some cases heavy row
  525. *> weighting can be treated with complete pivoting. The overhead in cases
  526. *> M much larger than N is then only due to pivoting, but the benefits in
  527. *> terms of accuracy have prevailed. The implementer/user can incorporate
  528. *> this extra QRF step easily. The implementer can also improve data movement
  529. *> (matrix transpose, matrix copy, matrix transposed copy) - this
  530. *> implementation of CGEJSV uses only the simplest, naive data movement.
  531. *> \endverbatim
  532. *
  533. *> \par Contributor:
  534. * ==================
  535. *>
  536. *> Zlatko Drmac (Zagreb, Croatia)
  537. *
  538. *> \par References:
  539. * ================
  540. *>
  541. *> \verbatim
  542. *>
  543. *> [1] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm I.
  544. *> SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1322-1342.
  545. *> LAPACK Working note 169.
  546. *> [2] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm II.
  547. *> SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1343-1362.
  548. *> LAPACK Working note 170.
  549. *> [3] Z. Drmac and Z. Bujanovic: On the failure of rank-revealing QR
  550. *> factorization software - a case study.
  551. *> ACM Trans. Math. Softw. Vol. 35, No 2 (2008), pp. 1-28.
  552. *> LAPACK Working note 176.
  553. *> [4] Z. Drmac: SIGMA - mathematical software library for accurate SVD, PSV,
  554. *> QSVD, (H,K)-SVD computations.
  555. *> Department of Mathematics, University of Zagreb, 2008, 2016.
  556. *> \endverbatim
  557. *
  558. *> \par Bugs, examples and comments:
  559. * =================================
  560. *>
  561. *> Please report all bugs and send interesting examples and/or comments to
  562. *> drmac@math.hr. Thank you.
  563. *>
  564. * =====================================================================
  565. SUBROUTINE CGEJSV( JOBA, JOBU, JOBV, JOBR, JOBT, JOBP,
  566. $ M, N, A, LDA, SVA, U, LDU, V, LDV,
  567. $ CWORK, LWORK, RWORK, LRWORK, IWORK, INFO )
  568. *
  569. * -- LAPACK computational routine --
  570. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  571. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  572. *
  573. * .. Scalar Arguments ..
  574. IMPLICIT NONE
  575. INTEGER INFO, LDA, LDU, LDV, LWORK, LRWORK, M, N
  576. * ..
  577. * .. Array Arguments ..
  578. COMPLEX A( LDA, * ), U( LDU, * ), V( LDV, * ), CWORK( LWORK )
  579. REAL SVA( N ), RWORK( LRWORK )
  580. INTEGER IWORK( * )
  581. CHARACTER*1 JOBA, JOBP, JOBR, JOBT, JOBU, JOBV
  582. * ..
  583. *
  584. * ===========================================================================
  585. *
  586. * .. Local Parameters ..
  587. REAL ZERO, ONE
  588. PARAMETER ( ZERO = 0.0E0, ONE = 1.0E0 )
  589. COMPLEX CZERO, CONE
  590. PARAMETER ( CZERO = ( 0.0E0, 0.0E0 ), CONE = ( 1.0E0, 0.0E0 ) )
  591. * ..
  592. * .. Local Scalars ..
  593. COMPLEX CTEMP
  594. REAL AAPP, AAQQ, AATMAX, AATMIN, BIG, BIG1, COND_OK,
  595. $ CONDR1, CONDR2, ENTRA, ENTRAT, EPSLN, MAXPRJ, SCALEM,
  596. $ SCONDA, SFMIN, SMALL, TEMP1, USCAL1, USCAL2, XSC
  597. INTEGER IERR, N1, NR, NUMRANK, p, q, WARNING
  598. LOGICAL ALMORT, DEFR, ERREST, GOSCAL, JRACC, KILL, LQUERY,
  599. $ LSVEC, L2ABER, L2KILL, L2PERT, L2RANK, L2TRAN, NOSCAL,
  600. $ ROWPIV, RSVEC, TRANSP
  601. *
  602. INTEGER OPTWRK, MINWRK, MINRWRK, MINIWRK
  603. INTEGER LWCON, LWLQF, LWQP3, LWQRF, LWUNMLQ, LWUNMQR, LWUNMQRM,
  604. $ LWSVDJ, LWSVDJV, LRWQP3, LRWCON, LRWSVDJ, IWOFF
  605. INTEGER LWRK_CGELQF, LWRK_CGEQP3, LWRK_CGEQP3N, LWRK_CGEQRF,
  606. $ LWRK_CGESVJ, LWRK_CGESVJV, LWRK_CGESVJU, LWRK_CUNMLQ,
  607. $ LWRK_CUNMQR, LWRK_CUNMQRM
  608. * ..
  609. * .. Local Arrays
  610. COMPLEX CDUMMY(1)
  611. REAL RDUMMY(1)
  612. *
  613. * .. Intrinsic Functions ..
  614. INTRINSIC ABS, CMPLX, CONJG, ALOG, MAX, MIN, REAL, NINT, SQRT
  615. * ..
  616. * .. External Functions ..
  617. REAL SLAMCH, SCNRM2
  618. INTEGER ISAMAX, ICAMAX
  619. LOGICAL LSAME
  620. EXTERNAL ISAMAX, ICAMAX, LSAME, SLAMCH, SCNRM2
  621. * ..
  622. * .. External Subroutines ..
  623. EXTERNAL SLASSQ, CCOPY, CGELQF, CGEQP3, CGEQRF, CLACPY, CLAPMR,
  624. $ CLASCL, SLASCL, CLASET, CLASSQ, CLASWP, CUNGQR, CUNMLQ,
  625. $ CUNMQR, CPOCON, SSCAL, CSSCAL, CSWAP, CTRSM, CLACGV,
  626. $ XERBLA
  627. *
  628. EXTERNAL CGESVJ
  629. * ..
  630. *
  631. * Test the input arguments
  632. *
  633. LSVEC = LSAME( JOBU, 'U' ) .OR. LSAME( JOBU, 'F' )
  634. JRACC = LSAME( JOBV, 'J' )
  635. RSVEC = LSAME( JOBV, 'V' ) .OR. JRACC
  636. ROWPIV = LSAME( JOBA, 'F' ) .OR. LSAME( JOBA, 'G' )
  637. L2RANK = LSAME( JOBA, 'R' )
  638. L2ABER = LSAME( JOBA, 'A' )
  639. ERREST = LSAME( JOBA, 'E' ) .OR. LSAME( JOBA, 'G' )
  640. L2TRAN = LSAME( JOBT, 'T' ) .AND. ( M .EQ. N )
  641. L2KILL = LSAME( JOBR, 'R' )
  642. DEFR = LSAME( JOBR, 'N' )
  643. L2PERT = LSAME( JOBP, 'P' )
  644. *
  645. LQUERY = ( LWORK .EQ. -1 ) .OR. ( LRWORK .EQ. -1 )
  646. *
  647. IF ( .NOT.(ROWPIV .OR. L2RANK .OR. L2ABER .OR.
  648. $ ERREST .OR. LSAME( JOBA, 'C' ) )) THEN
  649. INFO = - 1
  650. ELSE IF ( .NOT.( LSVEC .OR. LSAME( JOBU, 'N' ) .OR.
  651. $ ( LSAME( JOBU, 'W' ) .AND. RSVEC .AND. L2TRAN ) ) ) THEN
  652. INFO = - 2
  653. ELSE IF ( .NOT.( RSVEC .OR. LSAME( JOBV, 'N' ) .OR.
  654. $ ( LSAME( JOBV, 'W' ) .AND. LSVEC .AND. L2TRAN ) ) ) THEN
  655. INFO = - 3
  656. ELSE IF ( .NOT. ( L2KILL .OR. DEFR ) ) THEN
  657. INFO = - 4
  658. ELSE IF ( .NOT. ( LSAME(JOBT,'T') .OR. LSAME(JOBT,'N') ) ) THEN
  659. INFO = - 5
  660. ELSE IF ( .NOT. ( L2PERT .OR. LSAME( JOBP, 'N' ) ) ) THEN
  661. INFO = - 6
  662. ELSE IF ( M .LT. 0 ) THEN
  663. INFO = - 7
  664. ELSE IF ( ( N .LT. 0 ) .OR. ( N .GT. M ) ) THEN
  665. INFO = - 8
  666. ELSE IF ( LDA .LT. M ) THEN
  667. INFO = - 10
  668. ELSE IF ( LSVEC .AND. ( LDU .LT. M ) ) THEN
  669. INFO = - 13
  670. ELSE IF ( RSVEC .AND. ( LDV .LT. N ) ) THEN
  671. INFO = - 15
  672. ELSE
  673. * #:)
  674. INFO = 0
  675. END IF
  676. *
  677. IF ( INFO .EQ. 0 ) THEN
  678. * .. compute the minimal and the optimal workspace lengths
  679. * [[The expressions for computing the minimal and the optimal
  680. * values of LCWORK, LRWORK are written with a lot of redundancy and
  681. * can be simplified. However, this verbose form is useful for
  682. * maintenance and modifications of the code.]]
  683. *
  684. * .. minimal workspace length for CGEQP3 of an M x N matrix,
  685. * CGEQRF of an N x N matrix, CGELQF of an N x N matrix,
  686. * CUNMLQ for computing N x N matrix, CUNMQR for computing N x N
  687. * matrix, CUNMQR for computing M x N matrix, respectively.
  688. LWQP3 = N+1
  689. LWQRF = MAX( 1, N )
  690. LWLQF = MAX( 1, N )
  691. LWUNMLQ = MAX( 1, N )
  692. LWUNMQR = MAX( 1, N )
  693. LWUNMQRM = MAX( 1, M )
  694. * .. minimal workspace length for CPOCON of an N x N matrix
  695. LWCON = 2 * N
  696. * .. minimal workspace length for CGESVJ of an N x N matrix,
  697. * without and with explicit accumulation of Jacobi rotations
  698. LWSVDJ = MAX( 2 * N, 1 )
  699. LWSVDJV = MAX( 2 * N, 1 )
  700. * .. minimal REAL workspace length for CGEQP3, CPOCON, CGESVJ
  701. LRWQP3 = 2 * N
  702. LRWCON = N
  703. LRWSVDJ = N
  704. IF ( LQUERY ) THEN
  705. CALL CGEQP3( M, N, A, LDA, IWORK, CDUMMY, CDUMMY, -1,
  706. $ RDUMMY, IERR )
  707. LWRK_CGEQP3 = INT( CDUMMY(1) )
  708. CALL CGEQRF( N, N, A, LDA, CDUMMY, CDUMMY,-1, IERR )
  709. LWRK_CGEQRF = INT( CDUMMY(1) )
  710. CALL CGELQF( N, N, A, LDA, CDUMMY, CDUMMY,-1, IERR )
  711. LWRK_CGELQF = INT( CDUMMY(1) )
  712. END IF
  713. MINWRK = 2
  714. OPTWRK = 2
  715. MINIWRK = N
  716. IF ( .NOT. (LSVEC .OR. RSVEC ) ) THEN
  717. * .. minimal and optimal sizes of the complex workspace if
  718. * only the singular values are requested
  719. IF ( ERREST ) THEN
  720. MINWRK = MAX( N+LWQP3, N**2+LWCON, N+LWQRF, LWSVDJ )
  721. ELSE
  722. MINWRK = MAX( N+LWQP3, N+LWQRF, LWSVDJ )
  723. END IF
  724. IF ( LQUERY ) THEN
  725. CALL CGESVJ( 'L', 'N', 'N', N, N, A, LDA, SVA, N, V,
  726. $ LDV, CDUMMY, -1, RDUMMY, -1, IERR )
  727. LWRK_CGESVJ = INT( CDUMMY(1) )
  728. IF ( ERREST ) THEN
  729. OPTWRK = MAX( N+LWRK_CGEQP3, N**2+LWCON,
  730. $ N+LWRK_CGEQRF, LWRK_CGESVJ )
  731. ELSE
  732. OPTWRK = MAX( N+LWRK_CGEQP3, N+LWRK_CGEQRF,
  733. $ LWRK_CGESVJ )
  734. END IF
  735. END IF
  736. IF ( L2TRAN .OR. ROWPIV ) THEN
  737. IF ( ERREST ) THEN
  738. MINRWRK = MAX( 7, 2*M, LRWQP3, LRWCON, LRWSVDJ )
  739. ELSE
  740. MINRWRK = MAX( 7, 2*M, LRWQP3, LRWSVDJ )
  741. END IF
  742. ELSE
  743. IF ( ERREST ) THEN
  744. MINRWRK = MAX( 7, LRWQP3, LRWCON, LRWSVDJ )
  745. ELSE
  746. MINRWRK = MAX( 7, LRWQP3, LRWSVDJ )
  747. END IF
  748. END IF
  749. IF ( ROWPIV .OR. L2TRAN ) MINIWRK = MINIWRK + M
  750. ELSE IF ( RSVEC .AND. (.NOT.LSVEC) ) THEN
  751. * .. minimal and optimal sizes of the complex workspace if the
  752. * singular values and the right singular vectors are requested
  753. IF ( ERREST ) THEN
  754. MINWRK = MAX( N+LWQP3, LWCON, LWSVDJ, N+LWLQF,
  755. $ 2*N+LWQRF, N+LWSVDJ, N+LWUNMLQ )
  756. ELSE
  757. MINWRK = MAX( N+LWQP3, LWSVDJ, N+LWLQF, 2*N+LWQRF,
  758. $ N+LWSVDJ, N+LWUNMLQ )
  759. END IF
  760. IF ( LQUERY ) THEN
  761. CALL CGESVJ( 'L', 'U', 'N', N,N, U, LDU, SVA, N, A,
  762. $ LDA, CDUMMY, -1, RDUMMY, -1, IERR )
  763. LWRK_CGESVJ = INT( CDUMMY(1) )
  764. CALL CUNMLQ( 'L', 'C', N, N, N, A, LDA, CDUMMY,
  765. $ V, LDV, CDUMMY, -1, IERR )
  766. LWRK_CUNMLQ = INT( CDUMMY(1) )
  767. IF ( ERREST ) THEN
  768. OPTWRK = MAX( N+LWRK_CGEQP3, LWCON, LWRK_CGESVJ,
  769. $ N+LWRK_CGELQF, 2*N+LWRK_CGEQRF,
  770. $ N+LWRK_CGESVJ, N+LWRK_CUNMLQ )
  771. ELSE
  772. OPTWRK = MAX( N+LWRK_CGEQP3, LWRK_CGESVJ,N+LWRK_CGELQF,
  773. $ 2*N+LWRK_CGEQRF, N+LWRK_CGESVJ,
  774. $ N+LWRK_CUNMLQ )
  775. END IF
  776. END IF
  777. IF ( L2TRAN .OR. ROWPIV ) THEN
  778. IF ( ERREST ) THEN
  779. MINRWRK = MAX( 7, 2*M, LRWQP3, LRWSVDJ, LRWCON )
  780. ELSE
  781. MINRWRK = MAX( 7, 2*M, LRWQP3, LRWSVDJ )
  782. END IF
  783. ELSE
  784. IF ( ERREST ) THEN
  785. MINRWRK = MAX( 7, LRWQP3, LRWSVDJ, LRWCON )
  786. ELSE
  787. MINRWRK = MAX( 7, LRWQP3, LRWSVDJ )
  788. END IF
  789. END IF
  790. IF ( ROWPIV .OR. L2TRAN ) MINIWRK = MINIWRK + M
  791. ELSE IF ( LSVEC .AND. (.NOT.RSVEC) ) THEN
  792. * .. minimal and optimal sizes of the complex workspace if the
  793. * singular values and the left singular vectors are requested
  794. IF ( ERREST ) THEN
  795. MINWRK = N + MAX( LWQP3,LWCON,N+LWQRF,LWSVDJ,LWUNMQRM )
  796. ELSE
  797. MINWRK = N + MAX( LWQP3, N+LWQRF, LWSVDJ, LWUNMQRM )
  798. END IF
  799. IF ( LQUERY ) THEN
  800. CALL CGESVJ( 'L', 'U', 'N', N,N, U, LDU, SVA, N, A,
  801. $ LDA, CDUMMY, -1, RDUMMY, -1, IERR )
  802. LWRK_CGESVJ = INT( CDUMMY(1) )
  803. CALL CUNMQR( 'L', 'N', M, N, N, A, LDA, CDUMMY, U,
  804. $ LDU, CDUMMY, -1, IERR )
  805. LWRK_CUNMQRM = INT( CDUMMY(1) )
  806. IF ( ERREST ) THEN
  807. OPTWRK = N + MAX( LWRK_CGEQP3, LWCON, N+LWRK_CGEQRF,
  808. $ LWRK_CGESVJ, LWRK_CUNMQRM )
  809. ELSE
  810. OPTWRK = N + MAX( LWRK_CGEQP3, N+LWRK_CGEQRF,
  811. $ LWRK_CGESVJ, LWRK_CUNMQRM )
  812. END IF
  813. END IF
  814. IF ( L2TRAN .OR. ROWPIV ) THEN
  815. IF ( ERREST ) THEN
  816. MINRWRK = MAX( 7, 2*M, LRWQP3, LRWSVDJ, LRWCON )
  817. ELSE
  818. MINRWRK = MAX( 7, 2*M, LRWQP3, LRWSVDJ )
  819. END IF
  820. ELSE
  821. IF ( ERREST ) THEN
  822. MINRWRK = MAX( 7, LRWQP3, LRWSVDJ, LRWCON )
  823. ELSE
  824. MINRWRK = MAX( 7, LRWQP3, LRWSVDJ )
  825. END IF
  826. END IF
  827. IF ( ROWPIV .OR. L2TRAN ) MINIWRK = MINIWRK + M
  828. ELSE
  829. * .. minimal and optimal sizes of the complex workspace if the
  830. * full SVD is requested
  831. IF ( .NOT. JRACC ) THEN
  832. IF ( ERREST ) THEN
  833. MINWRK = MAX( N+LWQP3, N+LWCON, 2*N+N**2+LWCON,
  834. $ 2*N+LWQRF, 2*N+LWQP3,
  835. $ 2*N+N**2+N+LWLQF, 2*N+N**2+N+N**2+LWCON,
  836. $ 2*N+N**2+N+LWSVDJ, 2*N+N**2+N+LWSVDJV,
  837. $ 2*N+N**2+N+LWUNMQR,2*N+N**2+N+LWUNMLQ,
  838. $ N+N**2+LWSVDJ, N+LWUNMQRM )
  839. ELSE
  840. MINWRK = MAX( N+LWQP3, 2*N+N**2+LWCON,
  841. $ 2*N+LWQRF, 2*N+LWQP3,
  842. $ 2*N+N**2+N+LWLQF, 2*N+N**2+N+N**2+LWCON,
  843. $ 2*N+N**2+N+LWSVDJ, 2*N+N**2+N+LWSVDJV,
  844. $ 2*N+N**2+N+LWUNMQR,2*N+N**2+N+LWUNMLQ,
  845. $ N+N**2+LWSVDJ, N+LWUNMQRM )
  846. END IF
  847. MINIWRK = MINIWRK + N
  848. IF ( ROWPIV .OR. L2TRAN ) MINIWRK = MINIWRK + M
  849. ELSE
  850. IF ( ERREST ) THEN
  851. MINWRK = MAX( N+LWQP3, N+LWCON, 2*N+LWQRF,
  852. $ 2*N+N**2+LWSVDJV, 2*N+N**2+N+LWUNMQR,
  853. $ N+LWUNMQRM )
  854. ELSE
  855. MINWRK = MAX( N+LWQP3, 2*N+LWQRF,
  856. $ 2*N+N**2+LWSVDJV, 2*N+N**2+N+LWUNMQR,
  857. $ N+LWUNMQRM )
  858. END IF
  859. IF ( ROWPIV .OR. L2TRAN ) MINIWRK = MINIWRK + M
  860. END IF
  861. IF ( LQUERY ) THEN
  862. CALL CUNMQR( 'L', 'N', M, N, N, A, LDA, CDUMMY, U,
  863. $ LDU, CDUMMY, -1, IERR )
  864. LWRK_CUNMQRM = INT( CDUMMY(1) )
  865. CALL CUNMQR( 'L', 'N', N, N, N, A, LDA, CDUMMY, U,
  866. $ LDU, CDUMMY, -1, IERR )
  867. LWRK_CUNMQR = INT( CDUMMY(1) )
  868. IF ( .NOT. JRACC ) THEN
  869. CALL CGEQP3( N,N, A, LDA, IWORK, CDUMMY,CDUMMY, -1,
  870. $ RDUMMY, IERR )
  871. LWRK_CGEQP3N = INT( CDUMMY(1) )
  872. CALL CGESVJ( 'L', 'U', 'N', N, N, U, LDU, SVA,
  873. $ N, V, LDV, CDUMMY, -1, RDUMMY, -1, IERR )
  874. LWRK_CGESVJ = INT( CDUMMY(1) )
  875. CALL CGESVJ( 'U', 'U', 'N', N, N, U, LDU, SVA,
  876. $ N, V, LDV, CDUMMY, -1, RDUMMY, -1, IERR )
  877. LWRK_CGESVJU = INT( CDUMMY(1) )
  878. CALL CGESVJ( 'L', 'U', 'V', N, N, U, LDU, SVA,
  879. $ N, V, LDV, CDUMMY, -1, RDUMMY, -1, IERR )
  880. LWRK_CGESVJV = INT( CDUMMY(1) )
  881. CALL CUNMLQ( 'L', 'C', N, N, N, A, LDA, CDUMMY,
  882. $ V, LDV, CDUMMY, -1, IERR )
  883. LWRK_CUNMLQ = INT( CDUMMY(1) )
  884. IF ( ERREST ) THEN
  885. OPTWRK = MAX( N+LWRK_CGEQP3, N+LWCON,
  886. $ 2*N+N**2+LWCON, 2*N+LWRK_CGEQRF,
  887. $ 2*N+LWRK_CGEQP3N,
  888. $ 2*N+N**2+N+LWRK_CGELQF,
  889. $ 2*N+N**2+N+N**2+LWCON,
  890. $ 2*N+N**2+N+LWRK_CGESVJ,
  891. $ 2*N+N**2+N+LWRK_CGESVJV,
  892. $ 2*N+N**2+N+LWRK_CUNMQR,
  893. $ 2*N+N**2+N+LWRK_CUNMLQ,
  894. $ N+N**2+LWRK_CGESVJU,
  895. $ N+LWRK_CUNMQRM )
  896. ELSE
  897. OPTWRK = MAX( N+LWRK_CGEQP3,
  898. $ 2*N+N**2+LWCON, 2*N+LWRK_CGEQRF,
  899. $ 2*N+LWRK_CGEQP3N,
  900. $ 2*N+N**2+N+LWRK_CGELQF,
  901. $ 2*N+N**2+N+N**2+LWCON,
  902. $ 2*N+N**2+N+LWRK_CGESVJ,
  903. $ 2*N+N**2+N+LWRK_CGESVJV,
  904. $ 2*N+N**2+N+LWRK_CUNMQR,
  905. $ 2*N+N**2+N+LWRK_CUNMLQ,
  906. $ N+N**2+LWRK_CGESVJU,
  907. $ N+LWRK_CUNMQRM )
  908. END IF
  909. ELSE
  910. CALL CGESVJ( 'L', 'U', 'V', N, N, U, LDU, SVA,
  911. $ N, V, LDV, CDUMMY, -1, RDUMMY, -1, IERR )
  912. LWRK_CGESVJV = INT( CDUMMY(1) )
  913. CALL CUNMQR( 'L', 'N', N, N, N, CDUMMY, N, CDUMMY,
  914. $ V, LDV, CDUMMY, -1, IERR )
  915. LWRK_CUNMQR = INT( CDUMMY(1) )
  916. CALL CUNMQR( 'L', 'N', M, N, N, A, LDA, CDUMMY, U,
  917. $ LDU, CDUMMY, -1, IERR )
  918. LWRK_CUNMQRM = INT( CDUMMY(1) )
  919. IF ( ERREST ) THEN
  920. OPTWRK = MAX( N+LWRK_CGEQP3, N+LWCON,
  921. $ 2*N+LWRK_CGEQRF, 2*N+N**2,
  922. $ 2*N+N**2+LWRK_CGESVJV,
  923. $ 2*N+N**2+N+LWRK_CUNMQR,N+LWRK_CUNMQRM )
  924. ELSE
  925. OPTWRK = MAX( N+LWRK_CGEQP3, 2*N+LWRK_CGEQRF,
  926. $ 2*N+N**2, 2*N+N**2+LWRK_CGESVJV,
  927. $ 2*N+N**2+N+LWRK_CUNMQR,
  928. $ N+LWRK_CUNMQRM )
  929. END IF
  930. END IF
  931. END IF
  932. IF ( L2TRAN .OR. ROWPIV ) THEN
  933. MINRWRK = MAX( 7, 2*M, LRWQP3, LRWSVDJ, LRWCON )
  934. ELSE
  935. MINRWRK = MAX( 7, LRWQP3, LRWSVDJ, LRWCON )
  936. END IF
  937. END IF
  938. MINWRK = MAX( 2, MINWRK )
  939. OPTWRK = MAX( OPTWRK, MINWRK )
  940. IF ( LWORK .LT. MINWRK .AND. (.NOT.LQUERY) ) INFO = - 17
  941. IF ( LRWORK .LT. MINRWRK .AND. (.NOT.LQUERY) ) INFO = - 19
  942. END IF
  943. *
  944. IF ( INFO .NE. 0 ) THEN
  945. * #:(
  946. CALL XERBLA( 'CGEJSV', - INFO )
  947. RETURN
  948. ELSE IF ( LQUERY ) THEN
  949. CWORK(1) = OPTWRK
  950. CWORK(2) = MINWRK
  951. RWORK(1) = MINRWRK
  952. IWORK(1) = MAX( 4, MINIWRK )
  953. RETURN
  954. END IF
  955. *
  956. * Quick return for void matrix (Y3K safe)
  957. * #:)
  958. IF ( ( M .EQ. 0 ) .OR. ( N .EQ. 0 ) ) THEN
  959. IWORK(1:4) = 0
  960. RWORK(1:7) = 0
  961. RETURN
  962. ENDIF
  963. *
  964. * Determine whether the matrix U should be M x N or M x M
  965. *
  966. IF ( LSVEC ) THEN
  967. N1 = N
  968. IF ( LSAME( JOBU, 'F' ) ) N1 = M
  969. END IF
  970. *
  971. * Set numerical parameters
  972. *
  973. *! NOTE: Make sure SLAMCH() does not fail on the target architecture.
  974. *
  975. EPSLN = SLAMCH('Epsilon')
  976. SFMIN = SLAMCH('SafeMinimum')
  977. SMALL = SFMIN / EPSLN
  978. BIG = SLAMCH('O')
  979. * BIG = ONE / SFMIN
  980. *
  981. * Initialize SVA(1:N) = diag( ||A e_i||_2 )_1^N
  982. *
  983. *(!) If necessary, scale SVA() to protect the largest norm from
  984. * overflow. It is possible that this scaling pushes the smallest
  985. * column norm left from the underflow threshold (extreme case).
  986. *
  987. SCALEM = ONE / SQRT(REAL(M)*REAL(N))
  988. NOSCAL = .TRUE.
  989. GOSCAL = .TRUE.
  990. DO 1874 p = 1, N
  991. AAPP = ZERO
  992. AAQQ = ONE
  993. CALL CLASSQ( M, A(1,p), 1, AAPP, AAQQ )
  994. IF ( AAPP .GT. BIG ) THEN
  995. INFO = - 9
  996. CALL XERBLA( 'CGEJSV', -INFO )
  997. RETURN
  998. END IF
  999. AAQQ = SQRT(AAQQ)
  1000. IF ( ( AAPP .LT. (BIG / AAQQ) ) .AND. NOSCAL ) THEN
  1001. SVA(p) = AAPP * AAQQ
  1002. ELSE
  1003. NOSCAL = .FALSE.
  1004. SVA(p) = AAPP * ( AAQQ * SCALEM )
  1005. IF ( GOSCAL ) THEN
  1006. GOSCAL = .FALSE.
  1007. CALL SSCAL( p-1, SCALEM, SVA, 1 )
  1008. END IF
  1009. END IF
  1010. 1874 CONTINUE
  1011. *
  1012. IF ( NOSCAL ) SCALEM = ONE
  1013. *
  1014. AAPP = ZERO
  1015. AAQQ = BIG
  1016. DO 4781 p = 1, N
  1017. AAPP = MAX( AAPP, SVA(p) )
  1018. IF ( SVA(p) .NE. ZERO ) AAQQ = MIN( AAQQ, SVA(p) )
  1019. 4781 CONTINUE
  1020. *
  1021. * Quick return for zero M x N matrix
  1022. * #:)
  1023. IF ( AAPP .EQ. ZERO ) THEN
  1024. IF ( LSVEC ) CALL CLASET( 'G', M, N1, CZERO, CONE, U, LDU )
  1025. IF ( RSVEC ) CALL CLASET( 'G', N, N, CZERO, CONE, V, LDV )
  1026. RWORK(1) = ONE
  1027. RWORK(2) = ONE
  1028. IF ( ERREST ) RWORK(3) = ONE
  1029. IF ( LSVEC .AND. RSVEC ) THEN
  1030. RWORK(4) = ONE
  1031. RWORK(5) = ONE
  1032. END IF
  1033. IF ( L2TRAN ) THEN
  1034. RWORK(6) = ZERO
  1035. RWORK(7) = ZERO
  1036. END IF
  1037. IWORK(1) = 0
  1038. IWORK(2) = 0
  1039. IWORK(3) = 0
  1040. IWORK(4) = -1
  1041. RETURN
  1042. END IF
  1043. *
  1044. * Issue warning if denormalized column norms detected. Override the
  1045. * high relative accuracy request. Issue licence to kill nonzero columns
  1046. * (set them to zero) whose norm is less than sigma_max / BIG (roughly).
  1047. * #:(
  1048. WARNING = 0
  1049. IF ( AAQQ .LE. SFMIN ) THEN
  1050. L2RANK = .TRUE.
  1051. L2KILL = .TRUE.
  1052. WARNING = 1
  1053. END IF
  1054. *
  1055. * Quick return for one-column matrix
  1056. * #:)
  1057. IF ( N .EQ. 1 ) THEN
  1058. *
  1059. IF ( LSVEC ) THEN
  1060. CALL CLASCL( 'G',0,0,SVA(1),SCALEM, M,1,A(1,1),LDA,IERR )
  1061. CALL CLACPY( 'A', M, 1, A, LDA, U, LDU )
  1062. * computing all M left singular vectors of the M x 1 matrix
  1063. IF ( N1 .NE. N ) THEN
  1064. CALL CGEQRF( M, N, U,LDU, CWORK, CWORK(N+1),LWORK-N,IERR )
  1065. CALL CUNGQR( M,N1,1, U,LDU,CWORK,CWORK(N+1),LWORK-N,IERR )
  1066. CALL CCOPY( M, A(1,1), 1, U(1,1), 1 )
  1067. END IF
  1068. END IF
  1069. IF ( RSVEC ) THEN
  1070. V(1,1) = CONE
  1071. END IF
  1072. IF ( SVA(1) .LT. (BIG*SCALEM) ) THEN
  1073. SVA(1) = SVA(1) / SCALEM
  1074. SCALEM = ONE
  1075. END IF
  1076. RWORK(1) = ONE / SCALEM
  1077. RWORK(2) = ONE
  1078. IF ( SVA(1) .NE. ZERO ) THEN
  1079. IWORK(1) = 1
  1080. IF ( ( SVA(1) / SCALEM) .GE. SFMIN ) THEN
  1081. IWORK(2) = 1
  1082. ELSE
  1083. IWORK(2) = 0
  1084. END IF
  1085. ELSE
  1086. IWORK(1) = 0
  1087. IWORK(2) = 0
  1088. END IF
  1089. IWORK(3) = 0
  1090. IWORK(4) = -1
  1091. IF ( ERREST ) RWORK(3) = ONE
  1092. IF ( LSVEC .AND. RSVEC ) THEN
  1093. RWORK(4) = ONE
  1094. RWORK(5) = ONE
  1095. END IF
  1096. IF ( L2TRAN ) THEN
  1097. RWORK(6) = ZERO
  1098. RWORK(7) = ZERO
  1099. END IF
  1100. RETURN
  1101. *
  1102. END IF
  1103. *
  1104. TRANSP = .FALSE.
  1105. *
  1106. AATMAX = -ONE
  1107. AATMIN = BIG
  1108. IF ( ROWPIV .OR. L2TRAN ) THEN
  1109. *
  1110. * Compute the row norms, needed to determine row pivoting sequence
  1111. * (in the case of heavily row weighted A, row pivoting is strongly
  1112. * advised) and to collect information needed to compare the
  1113. * structures of A * A^* and A^* * A (in the case L2TRAN.EQ..TRUE.).
  1114. *
  1115. IF ( L2TRAN ) THEN
  1116. DO 1950 p = 1, M
  1117. XSC = ZERO
  1118. TEMP1 = ONE
  1119. CALL CLASSQ( N, A(p,1), LDA, XSC, TEMP1 )
  1120. * CLASSQ gets both the ell_2 and the ell_infinity norm
  1121. * in one pass through the vector
  1122. RWORK(M+p) = XSC * SCALEM
  1123. RWORK(p) = XSC * (SCALEM*SQRT(TEMP1))
  1124. AATMAX = MAX( AATMAX, RWORK(p) )
  1125. IF (RWORK(p) .NE. ZERO)
  1126. $ AATMIN = MIN(AATMIN,RWORK(p))
  1127. 1950 CONTINUE
  1128. ELSE
  1129. DO 1904 p = 1, M
  1130. RWORK(M+p) = SCALEM*ABS( A(p,ICAMAX(N,A(p,1),LDA)) )
  1131. AATMAX = MAX( AATMAX, RWORK(M+p) )
  1132. AATMIN = MIN( AATMIN, RWORK(M+p) )
  1133. 1904 CONTINUE
  1134. END IF
  1135. *
  1136. END IF
  1137. *
  1138. * For square matrix A try to determine whether A^* would be better
  1139. * input for the preconditioned Jacobi SVD, with faster convergence.
  1140. * The decision is based on an O(N) function of the vector of column
  1141. * and row norms of A, based on the Shannon entropy. This should give
  1142. * the right choice in most cases when the difference actually matters.
  1143. * It may fail and pick the slower converging side.
  1144. *
  1145. ENTRA = ZERO
  1146. ENTRAT = ZERO
  1147. IF ( L2TRAN ) THEN
  1148. *
  1149. XSC = ZERO
  1150. TEMP1 = ONE
  1151. CALL SLASSQ( N, SVA, 1, XSC, TEMP1 )
  1152. TEMP1 = ONE / TEMP1
  1153. *
  1154. ENTRA = ZERO
  1155. DO 1113 p = 1, N
  1156. BIG1 = ( ( SVA(p) / XSC )**2 ) * TEMP1
  1157. IF ( BIG1 .NE. ZERO ) ENTRA = ENTRA + BIG1 * ALOG(BIG1)
  1158. 1113 CONTINUE
  1159. ENTRA = - ENTRA / ALOG(REAL(N))
  1160. *
  1161. * Now, SVA().^2/Trace(A^* * A) is a point in the probability simplex.
  1162. * It is derived from the diagonal of A^* * A. Do the same with the
  1163. * diagonal of A * A^*, compute the entropy of the corresponding
  1164. * probability distribution. Note that A * A^* and A^* * A have the
  1165. * same trace.
  1166. *
  1167. ENTRAT = ZERO
  1168. DO 1114 p = 1, M
  1169. BIG1 = ( ( RWORK(p) / XSC )**2 ) * TEMP1
  1170. IF ( BIG1 .NE. ZERO ) ENTRAT = ENTRAT + BIG1 * ALOG(BIG1)
  1171. 1114 CONTINUE
  1172. ENTRAT = - ENTRAT / ALOG(REAL(M))
  1173. *
  1174. * Analyze the entropies and decide A or A^*. Smaller entropy
  1175. * usually means better input for the algorithm.
  1176. *
  1177. TRANSP = ( ENTRAT .LT. ENTRA )
  1178. *
  1179. * If A^* is better than A, take the adjoint of A. This is allowed
  1180. * only for square matrices, M=N.
  1181. IF ( TRANSP ) THEN
  1182. * In an optimal implementation, this trivial transpose
  1183. * should be replaced with faster transpose.
  1184. DO 1115 p = 1, N - 1
  1185. A(p,p) = CONJG(A(p,p))
  1186. DO 1116 q = p + 1, N
  1187. CTEMP = CONJG(A(q,p))
  1188. A(q,p) = CONJG(A(p,q))
  1189. A(p,q) = CTEMP
  1190. 1116 CONTINUE
  1191. 1115 CONTINUE
  1192. A(N,N) = CONJG(A(N,N))
  1193. DO 1117 p = 1, N
  1194. RWORK(M+p) = SVA(p)
  1195. SVA(p) = RWORK(p)
  1196. * previously computed row 2-norms are now column 2-norms
  1197. * of the transposed matrix
  1198. 1117 CONTINUE
  1199. TEMP1 = AAPP
  1200. AAPP = AATMAX
  1201. AATMAX = TEMP1
  1202. TEMP1 = AAQQ
  1203. AAQQ = AATMIN
  1204. AATMIN = TEMP1
  1205. KILL = LSVEC
  1206. LSVEC = RSVEC
  1207. RSVEC = KILL
  1208. IF ( LSVEC ) N1 = N
  1209. *
  1210. ROWPIV = .TRUE.
  1211. END IF
  1212. *
  1213. END IF
  1214. * END IF L2TRAN
  1215. *
  1216. * Scale the matrix so that its maximal singular value remains less
  1217. * than SQRT(BIG) -- the matrix is scaled so that its maximal column
  1218. * has Euclidean norm equal to SQRT(BIG/N). The only reason to keep
  1219. * SQRT(BIG) instead of BIG is the fact that CGEJSV uses LAPACK and
  1220. * BLAS routines that, in some implementations, are not capable of
  1221. * working in the full interval [SFMIN,BIG] and that they may provoke
  1222. * overflows in the intermediate results. If the singular values spread
  1223. * from SFMIN to BIG, then CGESVJ will compute them. So, in that case,
  1224. * one should use CGESVJ instead of CGEJSV.
  1225. BIG1 = SQRT( BIG )
  1226. TEMP1 = SQRT( BIG / REAL(N) )
  1227. * >> for future updates: allow bigger range, i.e. the largest column
  1228. * will be allowed up to BIG/N and CGESVJ will do the rest. However, for
  1229. * this all other (LAPACK) components must allow such a range.
  1230. * TEMP1 = BIG/REAL(N)
  1231. * TEMP1 = BIG * EPSLN this should 'almost' work with current LAPACK components
  1232. CALL SLASCL( 'G', 0, 0, AAPP, TEMP1, N, 1, SVA, N, IERR )
  1233. IF ( AAQQ .GT. (AAPP * SFMIN) ) THEN
  1234. AAQQ = ( AAQQ / AAPP ) * TEMP1
  1235. ELSE
  1236. AAQQ = ( AAQQ * TEMP1 ) / AAPP
  1237. END IF
  1238. TEMP1 = TEMP1 * SCALEM
  1239. CALL CLASCL( 'G', 0, 0, AAPP, TEMP1, M, N, A, LDA, IERR )
  1240. *
  1241. * To undo scaling at the end of this procedure, multiply the
  1242. * computed singular values with USCAL2 / USCAL1.
  1243. *
  1244. USCAL1 = TEMP1
  1245. USCAL2 = AAPP
  1246. *
  1247. IF ( L2KILL ) THEN
  1248. * L2KILL enforces computation of nonzero singular values in
  1249. * the restricted range of condition number of the initial A,
  1250. * sigma_max(A) / sigma_min(A) approx. SQRT(BIG)/SQRT(SFMIN).
  1251. XSC = SQRT( SFMIN )
  1252. ELSE
  1253. XSC = SMALL
  1254. *
  1255. * Now, if the condition number of A is too big,
  1256. * sigma_max(A) / sigma_min(A) .GT. SQRT(BIG/N) * EPSLN / SFMIN,
  1257. * as a precaution measure, the full SVD is computed using CGESVJ
  1258. * with accumulated Jacobi rotations. This provides numerically
  1259. * more robust computation, at the cost of slightly increased run
  1260. * time. Depending on the concrete implementation of BLAS and LAPACK
  1261. * (i.e. how they behave in presence of extreme ill-conditioning) the
  1262. * implementor may decide to remove this switch.
  1263. IF ( ( AAQQ.LT.SQRT(SFMIN) ) .AND. LSVEC .AND. RSVEC ) THEN
  1264. JRACC = .TRUE.
  1265. END IF
  1266. *
  1267. END IF
  1268. IF ( AAQQ .LT. XSC ) THEN
  1269. DO 700 p = 1, N
  1270. IF ( SVA(p) .LT. XSC ) THEN
  1271. CALL CLASET( 'A', M, 1, CZERO, CZERO, A(1,p), LDA )
  1272. SVA(p) = ZERO
  1273. END IF
  1274. 700 CONTINUE
  1275. END IF
  1276. *
  1277. * Preconditioning using QR factorization with pivoting
  1278. *
  1279. IF ( ROWPIV ) THEN
  1280. * Optional row permutation (Bjoerck row pivoting):
  1281. * A result by Cox and Higham shows that the Bjoerck's
  1282. * row pivoting combined with standard column pivoting
  1283. * has similar effect as Powell-Reid complete pivoting.
  1284. * The ell-infinity norms of A are made nonincreasing.
  1285. IF ( ( LSVEC .AND. RSVEC ) .AND. .NOT.( JRACC ) ) THEN
  1286. IWOFF = 2*N
  1287. ELSE
  1288. IWOFF = N
  1289. END IF
  1290. DO 1952 p = 1, M - 1
  1291. q = ISAMAX( M-p+1, RWORK(M+p), 1 ) + p - 1
  1292. IWORK(IWOFF+p) = q
  1293. IF ( p .NE. q ) THEN
  1294. TEMP1 = RWORK(M+p)
  1295. RWORK(M+p) = RWORK(M+q)
  1296. RWORK(M+q) = TEMP1
  1297. END IF
  1298. 1952 CONTINUE
  1299. CALL CLASWP( N, A, LDA, 1, M-1, IWORK(IWOFF+1), 1 )
  1300. END IF
  1301. *
  1302. * End of the preparation phase (scaling, optional sorting and
  1303. * transposing, optional flushing of small columns).
  1304. *
  1305. * Preconditioning
  1306. *
  1307. * If the full SVD is needed, the right singular vectors are computed
  1308. * from a matrix equation, and for that we need theoretical analysis
  1309. * of the Businger-Golub pivoting. So we use CGEQP3 as the first RR QRF.
  1310. * In all other cases the first RR QRF can be chosen by other criteria
  1311. * (eg speed by replacing global with restricted window pivoting, such
  1312. * as in xGEQPX from TOMS # 782). Good results will be obtained using
  1313. * xGEQPX with properly (!) chosen numerical parameters.
  1314. * Any improvement of CGEQP3 improves overall performance of CGEJSV.
  1315. *
  1316. * A * P1 = Q1 * [ R1^* 0]^*:
  1317. DO 1963 p = 1, N
  1318. * .. all columns are free columns
  1319. IWORK(p) = 0
  1320. 1963 CONTINUE
  1321. CALL CGEQP3( M, N, A, LDA, IWORK, CWORK, CWORK(N+1), LWORK-N,
  1322. $ RWORK, IERR )
  1323. *
  1324. * The upper triangular matrix R1 from the first QRF is inspected for
  1325. * rank deficiency and possibilities for deflation, or possible
  1326. * ill-conditioning. Depending on the user specified flag L2RANK,
  1327. * the procedure explores possibilities to reduce the numerical
  1328. * rank by inspecting the computed upper triangular factor. If
  1329. * L2RANK or L2ABER are up, then CGEJSV will compute the SVD of
  1330. * A + dA, where ||dA|| <= f(M,N)*EPSLN.
  1331. *
  1332. NR = 1
  1333. IF ( L2ABER ) THEN
  1334. * Standard absolute error bound suffices. All sigma_i with
  1335. * sigma_i < N*EPSLN*||A|| are flushed to zero. This is an
  1336. * aggressive enforcement of lower numerical rank by introducing a
  1337. * backward error of the order of N*EPSLN*||A||.
  1338. TEMP1 = SQRT(REAL(N))*EPSLN
  1339. DO 3001 p = 2, N
  1340. IF ( ABS(A(p,p)) .GE. (TEMP1*ABS(A(1,1))) ) THEN
  1341. NR = NR + 1
  1342. ELSE
  1343. GO TO 3002
  1344. END IF
  1345. 3001 CONTINUE
  1346. 3002 CONTINUE
  1347. ELSE IF ( L2RANK ) THEN
  1348. * .. similarly as above, only slightly more gentle (less aggressive).
  1349. * Sudden drop on the diagonal of R1 is used as the criterion for
  1350. * close-to-rank-deficient.
  1351. TEMP1 = SQRT(SFMIN)
  1352. DO 3401 p = 2, N
  1353. IF ( ( ABS(A(p,p)) .LT. (EPSLN*ABS(A(p-1,p-1))) ) .OR.
  1354. $ ( ABS(A(p,p)) .LT. SMALL ) .OR.
  1355. $ ( L2KILL .AND. (ABS(A(p,p)) .LT. TEMP1) ) ) GO TO 3402
  1356. NR = NR + 1
  1357. 3401 CONTINUE
  1358. 3402 CONTINUE
  1359. *
  1360. ELSE
  1361. * The goal is high relative accuracy. However, if the matrix
  1362. * has high scaled condition number the relative accuracy is in
  1363. * general not feasible. Later on, a condition number estimator
  1364. * will be deployed to estimate the scaled condition number.
  1365. * Here we just remove the underflowed part of the triangular
  1366. * factor. This prevents the situation in which the code is
  1367. * working hard to get the accuracy not warranted by the data.
  1368. TEMP1 = SQRT(SFMIN)
  1369. DO 3301 p = 2, N
  1370. IF ( ( ABS(A(p,p)) .LT. SMALL ) .OR.
  1371. $ ( L2KILL .AND. (ABS(A(p,p)) .LT. TEMP1) ) ) GO TO 3302
  1372. NR = NR + 1
  1373. 3301 CONTINUE
  1374. 3302 CONTINUE
  1375. *
  1376. END IF
  1377. *
  1378. ALMORT = .FALSE.
  1379. IF ( NR .EQ. N ) THEN
  1380. MAXPRJ = ONE
  1381. DO 3051 p = 2, N
  1382. TEMP1 = ABS(A(p,p)) / SVA(IWORK(p))
  1383. MAXPRJ = MIN( MAXPRJ, TEMP1 )
  1384. 3051 CONTINUE
  1385. IF ( MAXPRJ**2 .GE. ONE - REAL(N)*EPSLN ) ALMORT = .TRUE.
  1386. END IF
  1387. *
  1388. *
  1389. SCONDA = - ONE
  1390. CONDR1 = - ONE
  1391. CONDR2 = - ONE
  1392. *
  1393. IF ( ERREST ) THEN
  1394. IF ( N .EQ. NR ) THEN
  1395. IF ( RSVEC ) THEN
  1396. * .. V is available as workspace
  1397. CALL CLACPY( 'U', N, N, A, LDA, V, LDV )
  1398. DO 3053 p = 1, N
  1399. TEMP1 = SVA(IWORK(p))
  1400. CALL CSSCAL( p, ONE/TEMP1, V(1,p), 1 )
  1401. 3053 CONTINUE
  1402. IF ( LSVEC )THEN
  1403. CALL CPOCON( 'U', N, V, LDV, ONE, TEMP1,
  1404. $ CWORK(N+1), RWORK, IERR )
  1405. ELSE
  1406. CALL CPOCON( 'U', N, V, LDV, ONE, TEMP1,
  1407. $ CWORK, RWORK, IERR )
  1408. END IF
  1409. *
  1410. ELSE IF ( LSVEC ) THEN
  1411. * .. U is available as workspace
  1412. CALL CLACPY( 'U', N, N, A, LDA, U, LDU )
  1413. DO 3054 p = 1, N
  1414. TEMP1 = SVA(IWORK(p))
  1415. CALL CSSCAL( p, ONE/TEMP1, U(1,p), 1 )
  1416. 3054 CONTINUE
  1417. CALL CPOCON( 'U', N, U, LDU, ONE, TEMP1,
  1418. $ CWORK(N+1), RWORK, IERR )
  1419. ELSE
  1420. CALL CLACPY( 'U', N, N, A, LDA, CWORK, N )
  1421. *[] CALL CLACPY( 'U', N, N, A, LDA, CWORK(N+1), N )
  1422. * Change: here index shifted by N to the left, CWORK(1:N)
  1423. * not needed for SIGMA only computation
  1424. DO 3052 p = 1, N
  1425. TEMP1 = SVA(IWORK(p))
  1426. *[] CALL CSSCAL( p, ONE/TEMP1, CWORK(N+(p-1)*N+1), 1 )
  1427. CALL CSSCAL( p, ONE/TEMP1, CWORK((p-1)*N+1), 1 )
  1428. 3052 CONTINUE
  1429. * .. the columns of R are scaled to have unit Euclidean lengths.
  1430. *[] CALL CPOCON( 'U', N, CWORK(N+1), N, ONE, TEMP1,
  1431. *[] $ CWORK(N+N*N+1), RWORK, IERR )
  1432. CALL CPOCON( 'U', N, CWORK, N, ONE, TEMP1,
  1433. $ CWORK(N*N+1), RWORK, IERR )
  1434. *
  1435. END IF
  1436. IF ( TEMP1 .NE. ZERO ) THEN
  1437. SCONDA = ONE / SQRT(TEMP1)
  1438. ELSE
  1439. SCONDA = - ONE
  1440. END IF
  1441. * SCONDA is an estimate of SQRT(||(R^* * R)^(-1)||_1).
  1442. * N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
  1443. ELSE
  1444. SCONDA = - ONE
  1445. END IF
  1446. END IF
  1447. *
  1448. L2PERT = L2PERT .AND. ( ABS( A(1,1)/A(NR,NR) ) .GT. SQRT(BIG1) )
  1449. * If there is no violent scaling, artificial perturbation is not needed.
  1450. *
  1451. * Phase 3:
  1452. *
  1453. IF ( .NOT. ( RSVEC .OR. LSVEC ) ) THEN
  1454. *
  1455. * Singular Values only
  1456. *
  1457. * .. transpose A(1:NR,1:N)
  1458. DO 1946 p = 1, MIN( N-1, NR )
  1459. CALL CCOPY( N-p, A(p,p+1), LDA, A(p+1,p), 1 )
  1460. CALL CLACGV( N-p+1, A(p,p), 1 )
  1461. 1946 CONTINUE
  1462. IF ( NR .EQ. N ) A(N,N) = CONJG(A(N,N))
  1463. *
  1464. * The following two DO-loops introduce small relative perturbation
  1465. * into the strict upper triangle of the lower triangular matrix.
  1466. * Small entries below the main diagonal are also changed.
  1467. * This modification is useful if the computing environment does not
  1468. * provide/allow FLUSH TO ZERO underflow, for it prevents many
  1469. * annoying denormalized numbers in case of strongly scaled matrices.
  1470. * The perturbation is structured so that it does not introduce any
  1471. * new perturbation of the singular values, and it does not destroy
  1472. * the job done by the preconditioner.
  1473. * The licence for this perturbation is in the variable L2PERT, which
  1474. * should be .FALSE. if FLUSH TO ZERO underflow is active.
  1475. *
  1476. IF ( .NOT. ALMORT ) THEN
  1477. *
  1478. IF ( L2PERT ) THEN
  1479. * XSC = SQRT(SMALL)
  1480. XSC = EPSLN / REAL(N)
  1481. DO 4947 q = 1, NR
  1482. CTEMP = CMPLX(XSC*ABS(A(q,q)),ZERO)
  1483. DO 4949 p = 1, N
  1484. IF ( ( (p.GT.q) .AND. (ABS(A(p,q)).LE.TEMP1) )
  1485. $ .OR. ( p .LT. q ) )
  1486. * $ A(p,q) = TEMP1 * ( A(p,q) / ABS(A(p,q)) )
  1487. $ A(p,q) = CTEMP
  1488. 4949 CONTINUE
  1489. 4947 CONTINUE
  1490. ELSE
  1491. CALL CLASET( 'U', NR-1,NR-1, CZERO,CZERO, A(1,2),LDA )
  1492. END IF
  1493. *
  1494. * .. second preconditioning using the QR factorization
  1495. *
  1496. CALL CGEQRF( N,NR, A,LDA, CWORK, CWORK(N+1),LWORK-N, IERR )
  1497. *
  1498. * .. and transpose upper to lower triangular
  1499. DO 1948 p = 1, NR - 1
  1500. CALL CCOPY( NR-p, A(p,p+1), LDA, A(p+1,p), 1 )
  1501. CALL CLACGV( NR-p+1, A(p,p), 1 )
  1502. 1948 CONTINUE
  1503. *
  1504. END IF
  1505. *
  1506. * Row-cyclic Jacobi SVD algorithm with column pivoting
  1507. *
  1508. * .. again some perturbation (a "background noise") is added
  1509. * to drown denormals
  1510. IF ( L2PERT ) THEN
  1511. * XSC = SQRT(SMALL)
  1512. XSC = EPSLN / REAL(N)
  1513. DO 1947 q = 1, NR
  1514. CTEMP = CMPLX(XSC*ABS(A(q,q)),ZERO)
  1515. DO 1949 p = 1, NR
  1516. IF ( ( (p.GT.q) .AND. (ABS(A(p,q)).LE.TEMP1) )
  1517. $ .OR. ( p .LT. q ) )
  1518. * $ A(p,q) = TEMP1 * ( A(p,q) / ABS(A(p,q)) )
  1519. $ A(p,q) = CTEMP
  1520. 1949 CONTINUE
  1521. 1947 CONTINUE
  1522. ELSE
  1523. CALL CLASET( 'U', NR-1, NR-1, CZERO, CZERO, A(1,2), LDA )
  1524. END IF
  1525. *
  1526. * .. and one-sided Jacobi rotations are started on a lower
  1527. * triangular matrix (plus perturbation which is ignored in
  1528. * the part which destroys triangular form (confusing?!))
  1529. *
  1530. CALL CGESVJ( 'L', 'N', 'N', NR, NR, A, LDA, SVA,
  1531. $ N, V, LDV, CWORK, LWORK, RWORK, LRWORK, INFO )
  1532. *
  1533. SCALEM = RWORK(1)
  1534. NUMRANK = NINT(RWORK(2))
  1535. *
  1536. *
  1537. ELSE IF ( ( RSVEC .AND. ( .NOT. LSVEC ) .AND. ( .NOT. JRACC ) )
  1538. $ .OR.
  1539. $ ( JRACC .AND. ( .NOT. LSVEC ) .AND. ( NR .NE. N ) ) ) THEN
  1540. *
  1541. * -> Singular Values and Right Singular Vectors <-
  1542. *
  1543. IF ( ALMORT ) THEN
  1544. *
  1545. * .. in this case NR equals N
  1546. DO 1998 p = 1, NR
  1547. CALL CCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
  1548. CALL CLACGV( N-p+1, V(p,p), 1 )
  1549. 1998 CONTINUE
  1550. CALL CLASET( 'U', NR-1,NR-1, CZERO, CZERO, V(1,2), LDV )
  1551. *
  1552. CALL CGESVJ( 'L','U','N', N, NR, V, LDV, SVA, NR, A, LDA,
  1553. $ CWORK, LWORK, RWORK, LRWORK, INFO )
  1554. SCALEM = RWORK(1)
  1555. NUMRANK = NINT(RWORK(2))
  1556. ELSE
  1557. *
  1558. * .. two more QR factorizations ( one QRF is not enough, two require
  1559. * accumulated product of Jacobi rotations, three are perfect )
  1560. *
  1561. CALL CLASET( 'L', NR-1,NR-1, CZERO, CZERO, A(2,1), LDA )
  1562. CALL CGELQF( NR,N, A, LDA, CWORK, CWORK(N+1), LWORK-N, IERR)
  1563. CALL CLACPY( 'L', NR, NR, A, LDA, V, LDV )
  1564. CALL CLASET( 'U', NR-1,NR-1, CZERO, CZERO, V(1,2), LDV )
  1565. CALL CGEQRF( NR, NR, V, LDV, CWORK(N+1), CWORK(2*N+1),
  1566. $ LWORK-2*N, IERR )
  1567. DO 8998 p = 1, NR
  1568. CALL CCOPY( NR-p+1, V(p,p), LDV, V(p,p), 1 )
  1569. CALL CLACGV( NR-p+1, V(p,p), 1 )
  1570. 8998 CONTINUE
  1571. CALL CLASET('U', NR-1, NR-1, CZERO, CZERO, V(1,2), LDV)
  1572. *
  1573. CALL CGESVJ( 'L', 'U','N', NR, NR, V,LDV, SVA, NR, U,
  1574. $ LDU, CWORK(N+1), LWORK-N, RWORK, LRWORK, INFO )
  1575. SCALEM = RWORK(1)
  1576. NUMRANK = NINT(RWORK(2))
  1577. IF ( NR .LT. N ) THEN
  1578. CALL CLASET( 'A',N-NR, NR, CZERO,CZERO, V(NR+1,1), LDV )
  1579. CALL CLASET( 'A',NR, N-NR, CZERO,CZERO, V(1,NR+1), LDV )
  1580. CALL CLASET( 'A',N-NR,N-NR,CZERO,CONE, V(NR+1,NR+1),LDV )
  1581. END IF
  1582. *
  1583. CALL CUNMLQ( 'L', 'C', N, N, NR, A, LDA, CWORK,
  1584. $ V, LDV, CWORK(N+1), LWORK-N, IERR )
  1585. *
  1586. END IF
  1587. * .. permute the rows of V
  1588. * DO 8991 p = 1, N
  1589. * CALL CCOPY( N, V(p,1), LDV, A(IWORK(p),1), LDA )
  1590. * 8991 CONTINUE
  1591. * CALL CLACPY( 'All', N, N, A, LDA, V, LDV )
  1592. CALL CLAPMR( .FALSE., N, N, V, LDV, IWORK )
  1593. *
  1594. IF ( TRANSP ) THEN
  1595. CALL CLACPY( 'A', N, N, V, LDV, U, LDU )
  1596. END IF
  1597. *
  1598. ELSE IF ( JRACC .AND. (.NOT. LSVEC) .AND. ( NR.EQ. N ) ) THEN
  1599. *
  1600. CALL CLASET( 'L', N-1,N-1, CZERO, CZERO, A(2,1), LDA )
  1601. *
  1602. CALL CGESVJ( 'U','N','V', N, N, A, LDA, SVA, N, V, LDV,
  1603. $ CWORK, LWORK, RWORK, LRWORK, INFO )
  1604. SCALEM = RWORK(1)
  1605. NUMRANK = NINT(RWORK(2))
  1606. CALL CLAPMR( .FALSE., N, N, V, LDV, IWORK )
  1607. *
  1608. ELSE IF ( LSVEC .AND. ( .NOT. RSVEC ) ) THEN
  1609. *
  1610. * .. Singular Values and Left Singular Vectors ..
  1611. *
  1612. * .. second preconditioning step to avoid need to accumulate
  1613. * Jacobi rotations in the Jacobi iterations.
  1614. DO 1965 p = 1, NR
  1615. CALL CCOPY( N-p+1, A(p,p), LDA, U(p,p), 1 )
  1616. CALL CLACGV( N-p+1, U(p,p), 1 )
  1617. 1965 CONTINUE
  1618. CALL CLASET( 'U', NR-1, NR-1, CZERO, CZERO, U(1,2), LDU )
  1619. *
  1620. CALL CGEQRF( N, NR, U, LDU, CWORK(N+1), CWORK(2*N+1),
  1621. $ LWORK-2*N, IERR )
  1622. *
  1623. DO 1967 p = 1, NR - 1
  1624. CALL CCOPY( NR-p, U(p,p+1), LDU, U(p+1,p), 1 )
  1625. CALL CLACGV( N-p+1, U(p,p), 1 )
  1626. 1967 CONTINUE
  1627. CALL CLASET( 'U', NR-1, NR-1, CZERO, CZERO, U(1,2), LDU )
  1628. *
  1629. CALL CGESVJ( 'L', 'U', 'N', NR,NR, U, LDU, SVA, NR, A,
  1630. $ LDA, CWORK(N+1), LWORK-N, RWORK, LRWORK, INFO )
  1631. SCALEM = RWORK(1)
  1632. NUMRANK = NINT(RWORK(2))
  1633. *
  1634. IF ( NR .LT. M ) THEN
  1635. CALL CLASET( 'A', M-NR, NR,CZERO, CZERO, U(NR+1,1), LDU )
  1636. IF ( NR .LT. N1 ) THEN
  1637. CALL CLASET( 'A',NR, N1-NR, CZERO, CZERO, U(1,NR+1),LDU )
  1638. CALL CLASET( 'A',M-NR,N1-NR,CZERO,CONE,U(NR+1,NR+1),LDU )
  1639. END IF
  1640. END IF
  1641. *
  1642. CALL CUNMQR( 'L', 'N', M, N1, N, A, LDA, CWORK, U,
  1643. $ LDU, CWORK(N+1), LWORK-N, IERR )
  1644. *
  1645. IF ( ROWPIV )
  1646. $ CALL CLASWP( N1, U, LDU, 1, M-1, IWORK(IWOFF+1), -1 )
  1647. *
  1648. DO 1974 p = 1, N1
  1649. XSC = ONE / SCNRM2( M, U(1,p), 1 )
  1650. CALL CSSCAL( M, XSC, U(1,p), 1 )
  1651. 1974 CONTINUE
  1652. *
  1653. IF ( TRANSP ) THEN
  1654. CALL CLACPY( 'A', N, N, U, LDU, V, LDV )
  1655. END IF
  1656. *
  1657. ELSE
  1658. *
  1659. * .. Full SVD ..
  1660. *
  1661. IF ( .NOT. JRACC ) THEN
  1662. *
  1663. IF ( .NOT. ALMORT ) THEN
  1664. *
  1665. * Second Preconditioning Step (QRF [with pivoting])
  1666. * Note that the composition of TRANSPOSE, QRF and TRANSPOSE is
  1667. * equivalent to an LQF CALL. Since in many libraries the QRF
  1668. * seems to be better optimized than the LQF, we do explicit
  1669. * transpose and use the QRF. This is subject to changes in an
  1670. * optimized implementation of CGEJSV.
  1671. *
  1672. DO 1968 p = 1, NR
  1673. CALL CCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
  1674. CALL CLACGV( N-p+1, V(p,p), 1 )
  1675. 1968 CONTINUE
  1676. *
  1677. * .. the following two loops perturb small entries to avoid
  1678. * denormals in the second QR factorization, where they are
  1679. * as good as zeros. This is done to avoid painfully slow
  1680. * computation with denormals. The relative size of the perturbation
  1681. * is a parameter that can be changed by the implementer.
  1682. * This perturbation device will be obsolete on machines with
  1683. * properly implemented arithmetic.
  1684. * To switch it off, set L2PERT=.FALSE. To remove it from the
  1685. * code, remove the action under L2PERT=.TRUE., leave the ELSE part.
  1686. * The following two loops should be blocked and fused with the
  1687. * transposed copy above.
  1688. *
  1689. IF ( L2PERT ) THEN
  1690. XSC = SQRT(SMALL)
  1691. DO 2969 q = 1, NR
  1692. CTEMP = CMPLX(XSC*ABS( V(q,q) ),ZERO)
  1693. DO 2968 p = 1, N
  1694. IF ( ( p .GT. q ) .AND. ( ABS(V(p,q)) .LE. TEMP1 )
  1695. $ .OR. ( p .LT. q ) )
  1696. * $ V(p,q) = TEMP1 * ( V(p,q) / ABS(V(p,q)) )
  1697. $ V(p,q) = CTEMP
  1698. IF ( p .LT. q ) V(p,q) = - V(p,q)
  1699. 2968 CONTINUE
  1700. 2969 CONTINUE
  1701. ELSE
  1702. CALL CLASET( 'U', NR-1, NR-1, CZERO, CZERO, V(1,2), LDV )
  1703. END IF
  1704. *
  1705. * Estimate the row scaled condition number of R1
  1706. * (If R1 is rectangular, N > NR, then the condition number
  1707. * of the leading NR x NR submatrix is estimated.)
  1708. *
  1709. CALL CLACPY( 'L', NR, NR, V, LDV, CWORK(2*N+1), NR )
  1710. DO 3950 p = 1, NR
  1711. TEMP1 = SCNRM2(NR-p+1,CWORK(2*N+(p-1)*NR+p),1)
  1712. CALL CSSCAL(NR-p+1,ONE/TEMP1,CWORK(2*N+(p-1)*NR+p),1)
  1713. 3950 CONTINUE
  1714. CALL CPOCON('L',NR,CWORK(2*N+1),NR,ONE,TEMP1,
  1715. $ CWORK(2*N+NR*NR+1),RWORK,IERR)
  1716. CONDR1 = ONE / SQRT(TEMP1)
  1717. * .. here need a second opinion on the condition number
  1718. * .. then assume worst case scenario
  1719. * R1 is OK for inverse <=> CONDR1 .LT. REAL(N)
  1720. * more conservative <=> CONDR1 .LT. SQRT(REAL(N))
  1721. *
  1722. COND_OK = SQRT(SQRT(REAL(NR)))
  1723. *[TP] COND_OK is a tuning parameter.
  1724. *
  1725. IF ( CONDR1 .LT. COND_OK ) THEN
  1726. * .. the second QRF without pivoting. Note: in an optimized
  1727. * implementation, this QRF should be implemented as the QRF
  1728. * of a lower triangular matrix.
  1729. * R1^* = Q2 * R2
  1730. CALL CGEQRF( N, NR, V, LDV, CWORK(N+1), CWORK(2*N+1),
  1731. $ LWORK-2*N, IERR )
  1732. *
  1733. IF ( L2PERT ) THEN
  1734. XSC = SQRT(SMALL)/EPSLN
  1735. DO 3959 p = 2, NR
  1736. DO 3958 q = 1, p - 1
  1737. CTEMP=CMPLX(XSC*MIN(ABS(V(p,p)),ABS(V(q,q))),
  1738. $ ZERO)
  1739. IF ( ABS(V(q,p)) .LE. TEMP1 )
  1740. * $ V(q,p) = TEMP1 * ( V(q,p) / ABS(V(q,p)) )
  1741. $ V(q,p) = CTEMP
  1742. 3958 CONTINUE
  1743. 3959 CONTINUE
  1744. END IF
  1745. *
  1746. IF ( NR .NE. N )
  1747. $ CALL CLACPY( 'A', N, NR, V, LDV, CWORK(2*N+1), N )
  1748. * .. save ...
  1749. *
  1750. * .. this transposed copy should be better than naive
  1751. DO 1969 p = 1, NR - 1
  1752. CALL CCOPY( NR-p, V(p,p+1), LDV, V(p+1,p), 1 )
  1753. CALL CLACGV(NR-p+1, V(p,p), 1 )
  1754. 1969 CONTINUE
  1755. V(NR,NR)=CONJG(V(NR,NR))
  1756. *
  1757. CONDR2 = CONDR1
  1758. *
  1759. ELSE
  1760. *
  1761. * .. ill-conditioned case: second QRF with pivoting
  1762. * Note that windowed pivoting would be equally good
  1763. * numerically, and more run-time efficient. So, in
  1764. * an optimal implementation, the next call to CGEQP3
  1765. * should be replaced with eg. CALL CGEQPX (ACM TOMS #782)
  1766. * with properly (carefully) chosen parameters.
  1767. *
  1768. * R1^* * P2 = Q2 * R2
  1769. DO 3003 p = 1, NR
  1770. IWORK(N+p) = 0
  1771. 3003 CONTINUE
  1772. CALL CGEQP3( N, NR, V, LDV, IWORK(N+1), CWORK(N+1),
  1773. $ CWORK(2*N+1), LWORK-2*N, RWORK, IERR )
  1774. ** CALL CGEQRF( N, NR, V, LDV, CWORK(N+1), CWORK(2*N+1),
  1775. ** $ LWORK-2*N, IERR )
  1776. IF ( L2PERT ) THEN
  1777. XSC = SQRT(SMALL)
  1778. DO 3969 p = 2, NR
  1779. DO 3968 q = 1, p - 1
  1780. CTEMP=CMPLX(XSC*MIN(ABS(V(p,p)),ABS(V(q,q))),
  1781. $ ZERO)
  1782. IF ( ABS(V(q,p)) .LE. TEMP1 )
  1783. * $ V(q,p) = TEMP1 * ( V(q,p) / ABS(V(q,p)) )
  1784. $ V(q,p) = CTEMP
  1785. 3968 CONTINUE
  1786. 3969 CONTINUE
  1787. END IF
  1788. *
  1789. CALL CLACPY( 'A', N, NR, V, LDV, CWORK(2*N+1), N )
  1790. *
  1791. IF ( L2PERT ) THEN
  1792. XSC = SQRT(SMALL)
  1793. DO 8970 p = 2, NR
  1794. DO 8971 q = 1, p - 1
  1795. CTEMP=CMPLX(XSC*MIN(ABS(V(p,p)),ABS(V(q,q))),
  1796. $ ZERO)
  1797. * V(p,q) = - TEMP1*( V(q,p) / ABS(V(q,p)) )
  1798. V(p,q) = - CTEMP
  1799. 8971 CONTINUE
  1800. 8970 CONTINUE
  1801. ELSE
  1802. CALL CLASET( 'L',NR-1,NR-1,CZERO,CZERO,V(2,1),LDV )
  1803. END IF
  1804. * Now, compute R2 = L3 * Q3, the LQ factorization.
  1805. CALL CGELQF( NR, NR, V, LDV, CWORK(2*N+N*NR+1),
  1806. $ CWORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR, IERR )
  1807. * .. and estimate the condition number
  1808. CALL CLACPY( 'L',NR,NR,V,LDV,CWORK(2*N+N*NR+NR+1),NR )
  1809. DO 4950 p = 1, NR
  1810. TEMP1 = SCNRM2( p, CWORK(2*N+N*NR+NR+p), NR )
  1811. CALL CSSCAL( p, ONE/TEMP1, CWORK(2*N+N*NR+NR+p), NR )
  1812. 4950 CONTINUE
  1813. CALL CPOCON( 'L',NR,CWORK(2*N+N*NR+NR+1),NR,ONE,TEMP1,
  1814. $ CWORK(2*N+N*NR+NR+NR*NR+1),RWORK,IERR )
  1815. CONDR2 = ONE / SQRT(TEMP1)
  1816. *
  1817. *
  1818. IF ( CONDR2 .GE. COND_OK ) THEN
  1819. * .. save the Householder vectors used for Q3
  1820. * (this overwrites the copy of R2, as it will not be
  1821. * needed in this branch, but it does not overwritte the
  1822. * Huseholder vectors of Q2.).
  1823. CALL CLACPY( 'U', NR, NR, V, LDV, CWORK(2*N+1), N )
  1824. * .. and the rest of the information on Q3 is in
  1825. * WORK(2*N+N*NR+1:2*N+N*NR+N)
  1826. END IF
  1827. *
  1828. END IF
  1829. *
  1830. IF ( L2PERT ) THEN
  1831. XSC = SQRT(SMALL)
  1832. DO 4968 q = 2, NR
  1833. CTEMP = XSC * V(q,q)
  1834. DO 4969 p = 1, q - 1
  1835. * V(p,q) = - TEMP1*( V(p,q) / ABS(V(p,q)) )
  1836. V(p,q) = - CTEMP
  1837. 4969 CONTINUE
  1838. 4968 CONTINUE
  1839. ELSE
  1840. CALL CLASET( 'U', NR-1,NR-1, CZERO,CZERO, V(1,2), LDV )
  1841. END IF
  1842. *
  1843. * Second preconditioning finished; continue with Jacobi SVD
  1844. * The input matrix is lower trinagular.
  1845. *
  1846. * Recover the right singular vectors as solution of a well
  1847. * conditioned triangular matrix equation.
  1848. *
  1849. IF ( CONDR1 .LT. COND_OK ) THEN
  1850. *
  1851. CALL CGESVJ( 'L','U','N',NR,NR,V,LDV,SVA,NR,U, LDU,
  1852. $ CWORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,RWORK,
  1853. $ LRWORK, INFO )
  1854. SCALEM = RWORK(1)
  1855. NUMRANK = NINT(RWORK(2))
  1856. DO 3970 p = 1, NR
  1857. CALL CCOPY( NR, V(1,p), 1, U(1,p), 1 )
  1858. CALL CSSCAL( NR, SVA(p), V(1,p), 1 )
  1859. 3970 CONTINUE
  1860. * .. pick the right matrix equation and solve it
  1861. *
  1862. IF ( NR .EQ. N ) THEN
  1863. * :)) .. best case, R1 is inverted. The solution of this matrix
  1864. * equation is Q2*V2 = the product of the Jacobi rotations
  1865. * used in CGESVJ, premultiplied with the orthogonal matrix
  1866. * from the second QR factorization.
  1867. CALL CTRSM('L','U','N','N', NR,NR,CONE, A,LDA, V,LDV)
  1868. ELSE
  1869. * .. R1 is well conditioned, but non-square. Adjoint of R2
  1870. * is inverted to get the product of the Jacobi rotations
  1871. * used in CGESVJ. The Q-factor from the second QR
  1872. * factorization is then built in explicitly.
  1873. CALL CTRSM('L','U','C','N',NR,NR,CONE,CWORK(2*N+1),
  1874. $ N,V,LDV)
  1875. IF ( NR .LT. N ) THEN
  1876. CALL CLASET('A',N-NR,NR,CZERO,CZERO,V(NR+1,1),LDV)
  1877. CALL CLASET('A',NR,N-NR,CZERO,CZERO,V(1,NR+1),LDV)
  1878. CALL CLASET('A',N-NR,N-NR,CZERO,CONE,V(NR+1,NR+1),LDV)
  1879. END IF
  1880. CALL CUNMQR('L','N',N,N,NR,CWORK(2*N+1),N,CWORK(N+1),
  1881. $ V,LDV,CWORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR)
  1882. END IF
  1883. *
  1884. ELSE IF ( CONDR2 .LT. COND_OK ) THEN
  1885. *
  1886. * The matrix R2 is inverted. The solution of the matrix equation
  1887. * is Q3^* * V3 = the product of the Jacobi rotations (appplied to
  1888. * the lower triangular L3 from the LQ factorization of
  1889. * R2=L3*Q3), pre-multiplied with the transposed Q3.
  1890. CALL CGESVJ( 'L', 'U', 'N', NR, NR, V, LDV, SVA, NR, U,
  1891. $ LDU, CWORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR,
  1892. $ RWORK, LRWORK, INFO )
  1893. SCALEM = RWORK(1)
  1894. NUMRANK = NINT(RWORK(2))
  1895. DO 3870 p = 1, NR
  1896. CALL CCOPY( NR, V(1,p), 1, U(1,p), 1 )
  1897. CALL CSSCAL( NR, SVA(p), U(1,p), 1 )
  1898. 3870 CONTINUE
  1899. CALL CTRSM('L','U','N','N',NR,NR,CONE,CWORK(2*N+1),N,
  1900. $ U,LDU)
  1901. * .. apply the permutation from the second QR factorization
  1902. DO 873 q = 1, NR
  1903. DO 872 p = 1, NR
  1904. CWORK(2*N+N*NR+NR+IWORK(N+p)) = U(p,q)
  1905. 872 CONTINUE
  1906. DO 874 p = 1, NR
  1907. U(p,q) = CWORK(2*N+N*NR+NR+p)
  1908. 874 CONTINUE
  1909. 873 CONTINUE
  1910. IF ( NR .LT. N ) THEN
  1911. CALL CLASET( 'A',N-NR,NR,CZERO,CZERO,V(NR+1,1),LDV )
  1912. CALL CLASET( 'A',NR,N-NR,CZERO,CZERO,V(1,NR+1),LDV )
  1913. CALL CLASET('A',N-NR,N-NR,CZERO,CONE,V(NR+1,NR+1),LDV)
  1914. END IF
  1915. CALL CUNMQR( 'L','N',N,N,NR,CWORK(2*N+1),N,CWORK(N+1),
  1916. $ V,LDV,CWORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
  1917. ELSE
  1918. * Last line of defense.
  1919. * #:( This is a rather pathological case: no scaled condition
  1920. * improvement after two pivoted QR factorizations. Other
  1921. * possibility is that the rank revealing QR factorization
  1922. * or the condition estimator has failed, or the COND_OK
  1923. * is set very close to ONE (which is unnecessary). Normally,
  1924. * this branch should never be executed, but in rare cases of
  1925. * failure of the RRQR or condition estimator, the last line of
  1926. * defense ensures that CGEJSV completes the task.
  1927. * Compute the full SVD of L3 using CGESVJ with explicit
  1928. * accumulation of Jacobi rotations.
  1929. CALL CGESVJ( 'L', 'U', 'V', NR, NR, V, LDV, SVA, NR, U,
  1930. $ LDU, CWORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR,
  1931. $ RWORK, LRWORK, INFO )
  1932. SCALEM = RWORK(1)
  1933. NUMRANK = NINT(RWORK(2))
  1934. IF ( NR .LT. N ) THEN
  1935. CALL CLASET( 'A',N-NR,NR,CZERO,CZERO,V(NR+1,1),LDV )
  1936. CALL CLASET( 'A',NR,N-NR,CZERO,CZERO,V(1,NR+1),LDV )
  1937. CALL CLASET('A',N-NR,N-NR,CZERO,CONE,V(NR+1,NR+1),LDV)
  1938. END IF
  1939. CALL CUNMQR( 'L','N',N,N,NR,CWORK(2*N+1),N,CWORK(N+1),
  1940. $ V,LDV,CWORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
  1941. *
  1942. CALL CUNMLQ( 'L', 'C', NR, NR, NR, CWORK(2*N+1), N,
  1943. $ CWORK(2*N+N*NR+1), U, LDU, CWORK(2*N+N*NR+NR+1),
  1944. $ LWORK-2*N-N*NR-NR, IERR )
  1945. DO 773 q = 1, NR
  1946. DO 772 p = 1, NR
  1947. CWORK(2*N+N*NR+NR+IWORK(N+p)) = U(p,q)
  1948. 772 CONTINUE
  1949. DO 774 p = 1, NR
  1950. U(p,q) = CWORK(2*N+N*NR+NR+p)
  1951. 774 CONTINUE
  1952. 773 CONTINUE
  1953. *
  1954. END IF
  1955. *
  1956. * Permute the rows of V using the (column) permutation from the
  1957. * first QRF. Also, scale the columns to make them unit in
  1958. * Euclidean norm. This applies to all cases.
  1959. *
  1960. TEMP1 = SQRT(REAL(N)) * EPSLN
  1961. DO 1972 q = 1, N
  1962. DO 972 p = 1, N
  1963. CWORK(2*N+N*NR+NR+IWORK(p)) = V(p,q)
  1964. 972 CONTINUE
  1965. DO 973 p = 1, N
  1966. V(p,q) = CWORK(2*N+N*NR+NR+p)
  1967. 973 CONTINUE
  1968. XSC = ONE / SCNRM2( N, V(1,q), 1 )
  1969. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  1970. $ CALL CSSCAL( N, XSC, V(1,q), 1 )
  1971. 1972 CONTINUE
  1972. * At this moment, V contains the right singular vectors of A.
  1973. * Next, assemble the left singular vector matrix U (M x N).
  1974. IF ( NR .LT. M ) THEN
  1975. CALL CLASET('A', M-NR, NR, CZERO, CZERO, U(NR+1,1), LDU)
  1976. IF ( NR .LT. N1 ) THEN
  1977. CALL CLASET('A',NR,N1-NR,CZERO,CZERO,U(1,NR+1),LDU)
  1978. CALL CLASET('A',M-NR,N1-NR,CZERO,CONE,
  1979. $ U(NR+1,NR+1),LDU)
  1980. END IF
  1981. END IF
  1982. *
  1983. * The Q matrix from the first QRF is built into the left singular
  1984. * matrix U. This applies to all cases.
  1985. *
  1986. CALL CUNMQR( 'L', 'N', M, N1, N, A, LDA, CWORK, U,
  1987. $ LDU, CWORK(N+1), LWORK-N, IERR )
  1988. * The columns of U are normalized. The cost is O(M*N) flops.
  1989. TEMP1 = SQRT(REAL(M)) * EPSLN
  1990. DO 1973 p = 1, NR
  1991. XSC = ONE / SCNRM2( M, U(1,p), 1 )
  1992. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  1993. $ CALL CSSCAL( M, XSC, U(1,p), 1 )
  1994. 1973 CONTINUE
  1995. *
  1996. * If the initial QRF is computed with row pivoting, the left
  1997. * singular vectors must be adjusted.
  1998. *
  1999. IF ( ROWPIV )
  2000. $ CALL CLASWP( N1, U, LDU, 1, M-1, IWORK(IWOFF+1), -1 )
  2001. *
  2002. ELSE
  2003. *
  2004. * .. the initial matrix A has almost orthogonal columns and
  2005. * the second QRF is not needed
  2006. *
  2007. CALL CLACPY( 'U', N, N, A, LDA, CWORK(N+1), N )
  2008. IF ( L2PERT ) THEN
  2009. XSC = SQRT(SMALL)
  2010. DO 5970 p = 2, N
  2011. CTEMP = XSC * CWORK( N + (p-1)*N + p )
  2012. DO 5971 q = 1, p - 1
  2013. * CWORK(N+(q-1)*N+p)=-TEMP1 * ( CWORK(N+(p-1)*N+q) /
  2014. * $ ABS(CWORK(N+(p-1)*N+q)) )
  2015. CWORK(N+(q-1)*N+p)=-CTEMP
  2016. 5971 CONTINUE
  2017. 5970 CONTINUE
  2018. ELSE
  2019. CALL CLASET( 'L',N-1,N-1,CZERO,CZERO,CWORK(N+2),N )
  2020. END IF
  2021. *
  2022. CALL CGESVJ( 'U', 'U', 'N', N, N, CWORK(N+1), N, SVA,
  2023. $ N, U, LDU, CWORK(N+N*N+1), LWORK-N-N*N, RWORK, LRWORK,
  2024. $ INFO )
  2025. *
  2026. SCALEM = RWORK(1)
  2027. NUMRANK = NINT(RWORK(2))
  2028. DO 6970 p = 1, N
  2029. CALL CCOPY( N, CWORK(N+(p-1)*N+1), 1, U(1,p), 1 )
  2030. CALL CSSCAL( N, SVA(p), CWORK(N+(p-1)*N+1), 1 )
  2031. 6970 CONTINUE
  2032. *
  2033. CALL CTRSM( 'L', 'U', 'N', 'N', N, N,
  2034. $ CONE, A, LDA, CWORK(N+1), N )
  2035. DO 6972 p = 1, N
  2036. CALL CCOPY( N, CWORK(N+p), N, V(IWORK(p),1), LDV )
  2037. 6972 CONTINUE
  2038. TEMP1 = SQRT(REAL(N))*EPSLN
  2039. DO 6971 p = 1, N
  2040. XSC = ONE / SCNRM2( N, V(1,p), 1 )
  2041. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  2042. $ CALL CSSCAL( N, XSC, V(1,p), 1 )
  2043. 6971 CONTINUE
  2044. *
  2045. * Assemble the left singular vector matrix U (M x N).
  2046. *
  2047. IF ( N .LT. M ) THEN
  2048. CALL CLASET( 'A', M-N, N, CZERO, CZERO, U(N+1,1), LDU )
  2049. IF ( N .LT. N1 ) THEN
  2050. CALL CLASET('A',N, N1-N, CZERO, CZERO, U(1,N+1),LDU)
  2051. CALL CLASET( 'A',M-N,N1-N, CZERO, CONE,U(N+1,N+1),LDU)
  2052. END IF
  2053. END IF
  2054. CALL CUNMQR( 'L', 'N', M, N1, N, A, LDA, CWORK, U,
  2055. $ LDU, CWORK(N+1), LWORK-N, IERR )
  2056. TEMP1 = SQRT(REAL(M))*EPSLN
  2057. DO 6973 p = 1, N1
  2058. XSC = ONE / SCNRM2( M, U(1,p), 1 )
  2059. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  2060. $ CALL CSSCAL( M, XSC, U(1,p), 1 )
  2061. 6973 CONTINUE
  2062. *
  2063. IF ( ROWPIV )
  2064. $ CALL CLASWP( N1, U, LDU, 1, M-1, IWORK(IWOFF+1), -1 )
  2065. *
  2066. END IF
  2067. *
  2068. * end of the >> almost orthogonal case << in the full SVD
  2069. *
  2070. ELSE
  2071. *
  2072. * This branch deploys a preconditioned Jacobi SVD with explicitly
  2073. * accumulated rotations. It is included as optional, mainly for
  2074. * experimental purposes. It does perform well, and can also be used.
  2075. * In this implementation, this branch will be automatically activated
  2076. * if the condition number sigma_max(A) / sigma_min(A) is predicted
  2077. * to be greater than the overflow threshold. This is because the
  2078. * a posteriori computation of the singular vectors assumes robust
  2079. * implementation of BLAS and some LAPACK procedures, capable of working
  2080. * in presence of extreme values, e.g. when the singular values spread from
  2081. * the underflow to the overflow threshold.
  2082. *
  2083. DO 7968 p = 1, NR
  2084. CALL CCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
  2085. CALL CLACGV( N-p+1, V(p,p), 1 )
  2086. 7968 CONTINUE
  2087. *
  2088. IF ( L2PERT ) THEN
  2089. XSC = SQRT(SMALL/EPSLN)
  2090. DO 5969 q = 1, NR
  2091. CTEMP = CMPLX(XSC*ABS( V(q,q) ),ZERO)
  2092. DO 5968 p = 1, N
  2093. IF ( ( p .GT. q ) .AND. ( ABS(V(p,q)) .LE. TEMP1 )
  2094. $ .OR. ( p .LT. q ) )
  2095. * $ V(p,q) = TEMP1 * ( V(p,q) / ABS(V(p,q)) )
  2096. $ V(p,q) = CTEMP
  2097. IF ( p .LT. q ) V(p,q) = - V(p,q)
  2098. 5968 CONTINUE
  2099. 5969 CONTINUE
  2100. ELSE
  2101. CALL CLASET( 'U', NR-1, NR-1, CZERO, CZERO, V(1,2), LDV )
  2102. END IF
  2103. CALL CGEQRF( N, NR, V, LDV, CWORK(N+1), CWORK(2*N+1),
  2104. $ LWORK-2*N, IERR )
  2105. CALL CLACPY( 'L', N, NR, V, LDV, CWORK(2*N+1), N )
  2106. *
  2107. DO 7969 p = 1, NR
  2108. CALL CCOPY( NR-p+1, V(p,p), LDV, U(p,p), 1 )
  2109. CALL CLACGV( NR-p+1, U(p,p), 1 )
  2110. 7969 CONTINUE
  2111. IF ( L2PERT ) THEN
  2112. XSC = SQRT(SMALL/EPSLN)
  2113. DO 9970 q = 2, NR
  2114. DO 9971 p = 1, q - 1
  2115. CTEMP = CMPLX(XSC * MIN(ABS(U(p,p)),ABS(U(q,q))),
  2116. $ ZERO)
  2117. * U(p,q) = - TEMP1 * ( U(q,p) / ABS(U(q,p)) )
  2118. U(p,q) = - CTEMP
  2119. 9971 CONTINUE
  2120. 9970 CONTINUE
  2121. ELSE
  2122. CALL CLASET('U', NR-1, NR-1, CZERO, CZERO, U(1,2), LDU )
  2123. END IF
  2124. CALL CGESVJ( 'L', 'U', 'V', NR, NR, U, LDU, SVA,
  2125. $ N, V, LDV, CWORK(2*N+N*NR+1), LWORK-2*N-N*NR,
  2126. $ RWORK, LRWORK, INFO )
  2127. SCALEM = RWORK(1)
  2128. NUMRANK = NINT(RWORK(2))
  2129. IF ( NR .LT. N ) THEN
  2130. CALL CLASET( 'A',N-NR,NR,CZERO,CZERO,V(NR+1,1),LDV )
  2131. CALL CLASET( 'A',NR,N-NR,CZERO,CZERO,V(1,NR+1),LDV )
  2132. CALL CLASET( 'A',N-NR,N-NR,CZERO,CONE,V(NR+1,NR+1),LDV )
  2133. END IF
  2134. CALL CUNMQR( 'L','N',N,N,NR,CWORK(2*N+1),N,CWORK(N+1),
  2135. $ V,LDV,CWORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
  2136. *
  2137. * Permute the rows of V using the (column) permutation from the
  2138. * first QRF. Also, scale the columns to make them unit in
  2139. * Euclidean norm. This applies to all cases.
  2140. *
  2141. TEMP1 = SQRT(REAL(N)) * EPSLN
  2142. DO 7972 q = 1, N
  2143. DO 8972 p = 1, N
  2144. CWORK(2*N+N*NR+NR+IWORK(p)) = V(p,q)
  2145. 8972 CONTINUE
  2146. DO 8973 p = 1, N
  2147. V(p,q) = CWORK(2*N+N*NR+NR+p)
  2148. 8973 CONTINUE
  2149. XSC = ONE / SCNRM2( N, V(1,q), 1 )
  2150. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  2151. $ CALL CSSCAL( N, XSC, V(1,q), 1 )
  2152. 7972 CONTINUE
  2153. *
  2154. * At this moment, V contains the right singular vectors of A.
  2155. * Next, assemble the left singular vector matrix U (M x N).
  2156. *
  2157. IF ( NR .LT. M ) THEN
  2158. CALL CLASET( 'A', M-NR, NR, CZERO, CZERO, U(NR+1,1), LDU )
  2159. IF ( NR .LT. N1 ) THEN
  2160. CALL CLASET('A',NR, N1-NR, CZERO, CZERO, U(1,NR+1),LDU)
  2161. CALL CLASET('A',M-NR,N1-NR, CZERO, CONE,U(NR+1,NR+1),LDU)
  2162. END IF
  2163. END IF
  2164. *
  2165. CALL CUNMQR( 'L', 'N', M, N1, N, A, LDA, CWORK, U,
  2166. $ LDU, CWORK(N+1), LWORK-N, IERR )
  2167. *
  2168. IF ( ROWPIV )
  2169. $ CALL CLASWP( N1, U, LDU, 1, M-1, IWORK(IWOFF+1), -1 )
  2170. *
  2171. *
  2172. END IF
  2173. IF ( TRANSP ) THEN
  2174. * .. swap U and V because the procedure worked on A^*
  2175. DO 6974 p = 1, N
  2176. CALL CSWAP( N, U(1,p), 1, V(1,p), 1 )
  2177. 6974 CONTINUE
  2178. END IF
  2179. *
  2180. END IF
  2181. * end of the full SVD
  2182. *
  2183. * Undo scaling, if necessary (and possible)
  2184. *
  2185. IF ( USCAL2 .LE. (BIG/SVA(1))*USCAL1 ) THEN
  2186. CALL SLASCL( 'G', 0, 0, USCAL1, USCAL2, NR, 1, SVA, N, IERR )
  2187. USCAL1 = ONE
  2188. USCAL2 = ONE
  2189. END IF
  2190. *
  2191. IF ( NR .LT. N ) THEN
  2192. DO 3004 p = NR+1, N
  2193. SVA(p) = ZERO
  2194. 3004 CONTINUE
  2195. END IF
  2196. *
  2197. RWORK(1) = USCAL2 * SCALEM
  2198. RWORK(2) = USCAL1
  2199. IF ( ERREST ) RWORK(3) = SCONDA
  2200. IF ( LSVEC .AND. RSVEC ) THEN
  2201. RWORK(4) = CONDR1
  2202. RWORK(5) = CONDR2
  2203. END IF
  2204. IF ( L2TRAN ) THEN
  2205. RWORK(6) = ENTRA
  2206. RWORK(7) = ENTRAT
  2207. END IF
  2208. *
  2209. IWORK(1) = NR
  2210. IWORK(2) = NUMRANK
  2211. IWORK(3) = WARNING
  2212. IF ( TRANSP ) THEN
  2213. IWORK(4) = 1
  2214. ELSE
  2215. IWORK(4) = -1
  2216. END IF
  2217. *
  2218. RETURN
  2219. * ..
  2220. * .. END OF CGEJSV
  2221. * ..
  2222. END
  2223. *