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.

cgesvdq.f 59 kB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391
  1. *> \brief <b> CGESVDQ computes the singular value decomposition (SVD) with a QR-Preconditioned QR SVD Method for GE matrices</b>
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download CGESVDQ + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgesvdq.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgesvdq.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgesvdq.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE CGESVDQ( JOBA, JOBP, JOBR, JOBU, JOBV, M, N, A, LDA,
  22. * S, U, LDU, V, LDV, NUMRANK, IWORK, LIWORK,
  23. * CWORK, LCWORK, RWORK, LRWORK, INFO )
  24. *
  25. * .. Scalar Arguments ..
  26. * IMPLICIT NONE
  27. * CHARACTER JOBA, JOBP, JOBR, JOBU, JOBV
  28. * INTEGER M, N, LDA, LDU, LDV, NUMRANK, LIWORK, LCWORK, LRWORK,
  29. * INFO
  30. * ..
  31. * .. Array Arguments ..
  32. * COMPLEX A( LDA, * ), U( LDU, * ), V( LDV, * ), CWORK( * )
  33. * REAL S( * ), RWORK( * )
  34. * INTEGER IWORK( * )
  35. * ..
  36. *
  37. *
  38. *> \par Purpose:
  39. * =============
  40. *>
  41. *> \verbatim
  42. *>
  43. *> CGESVDQ computes the singular value decomposition (SVD) of a complex
  44. *> M-by-N matrix A, where M >= N. The SVD of A is written as
  45. *> [++] [xx] [x0] [xx]
  46. *> A = U * SIGMA * V^*, [++] = [xx] * [ox] * [xx]
  47. *> [++] [xx]
  48. *> where SIGMA is an N-by-N diagonal matrix, U is an M-by-N orthonormal
  49. *> matrix, and V is an N-by-N unitary matrix. The diagonal elements
  50. *> of SIGMA are the singular values of A. The columns of U and V are the
  51. *> left and the right singular vectors of A, respectively.
  52. *> \endverbatim
  53. *
  54. * Arguments:
  55. * ==========
  56. *
  57. *> \param[in] JOBA
  58. *> \verbatim
  59. *> JOBA is CHARACTER*1
  60. *> Specifies the level of accuracy in the computed SVD
  61. *> = 'A' The requested accuracy corresponds to having the backward
  62. *> error bounded by || delta A ||_F <= f(m,n) * EPS * || A ||_F,
  63. *> where EPS = SLAMCH('Epsilon'). This authorises CGESVDQ to
  64. *> truncate the computed triangular factor in a rank revealing
  65. *> QR factorization whenever the truncated part is below the
  66. *> threshold of the order of EPS * ||A||_F. This is aggressive
  67. *> truncation level.
  68. *> = 'M' Similarly as with 'A', but the truncation is more gentle: it
  69. *> is allowed only when there is a drop on the diagonal of the
  70. *> triangular factor in the QR factorization. This is medium
  71. *> truncation level.
  72. *> = 'H' High accuracy requested. No numerical rank determination based
  73. *> on the rank revealing QR factorization is attempted.
  74. *> = 'E' Same as 'H', and in addition the condition number of column
  75. *> scaled A is estimated and returned in RWORK(1).
  76. *> N^(-1/4)*RWORK(1) <= ||pinv(A_scaled)||_2 <= N^(1/4)*RWORK(1)
  77. *> \endverbatim
  78. *>
  79. *> \param[in] JOBP
  80. *> \verbatim
  81. *> JOBP is CHARACTER*1
  82. *> = 'P' The rows of A are ordered in decreasing order with respect to
  83. *> ||A(i,:)||_\infty. This enhances numerical accuracy at the cost
  84. *> of extra data movement. Recommended for numerical robustness.
  85. *> = 'N' No row pivoting.
  86. *> \endverbatim
  87. *>
  88. *> \param[in] JOBR
  89. *> \verbatim
  90. *> JOBR is CHARACTER*1
  91. *> = 'T' After the initial pivoted QR factorization, CGESVD is applied to
  92. *> the adjoint R**H of the computed triangular factor R. This involves
  93. *> some extra data movement (matrix transpositions). Useful for
  94. *> experiments, research and development.
  95. *> = 'N' The triangular factor R is given as input to CGESVD. This may be
  96. *> preferred as it involves less data movement.
  97. *> \endverbatim
  98. *>
  99. *> \param[in] JOBU
  100. *> \verbatim
  101. *> JOBU is CHARACTER*1
  102. *> = 'A' All M left singular vectors are computed and returned in the
  103. *> matrix U. See the description of U.
  104. *> = 'S' or 'U' N = min(M,N) left singular vectors are computed and returned
  105. *> in the matrix U. See the description of U.
  106. *> = 'R' Numerical rank NUMRANK is determined and only NUMRANK left singular
  107. *> vectors are computed and returned in the matrix U.
  108. *> = 'F' The N left singular vectors are returned in factored form as the
  109. *> product of the Q factor from the initial QR factorization and the
  110. *> N left singular vectors of (R**H , 0)**H. If row pivoting is used,
  111. *> then the necessary information on the row pivoting is stored in
  112. *> IWORK(N+1:N+M-1).
  113. *> = 'N' The left singular vectors are not computed.
  114. *> \endverbatim
  115. *>
  116. *> \param[in] JOBV
  117. *> \verbatim
  118. *> JOBV is CHARACTER*1
  119. *> = 'A', 'V' All N right singular vectors are computed and returned in
  120. *> the matrix V.
  121. *> = 'R' Numerical rank NUMRANK is determined and only NUMRANK right singular
  122. *> vectors are computed and returned in the matrix V. This option is
  123. *> allowed only if JOBU = 'R' or JOBU = 'N'; otherwise it is illegal.
  124. *> = 'N' The right singular vectors are not computed.
  125. *> \endverbatim
  126. *>
  127. *> \param[in] M
  128. *> \verbatim
  129. *> M is INTEGER
  130. *> The number of rows of the input matrix A. M >= 0.
  131. *> \endverbatim
  132. *>
  133. *> \param[in] N
  134. *> \verbatim
  135. *> N is INTEGER
  136. *> The number of columns of the input matrix A. M >= N >= 0.
  137. *> \endverbatim
  138. *>
  139. *> \param[in,out] A
  140. *> \verbatim
  141. *> A is COMPLEX array of dimensions LDA x N
  142. *> On entry, the input matrix A.
  143. *> On exit, if JOBU .NE. 'N' or JOBV .NE. 'N', the lower triangle of A contains
  144. *> the Householder vectors as stored by CGEQP3. If JOBU = 'F', these Householder
  145. *> vectors together with CWORK(1:N) can be used to restore the Q factors from
  146. *> the initial pivoted QR factorization of A. See the description of U.
  147. *> \endverbatim
  148. *>
  149. *> \param[in] LDA
  150. *> \verbatim
  151. *> LDA is INTEGER.
  152. *> The leading dimension of the array A. LDA >= max(1,M).
  153. *> \endverbatim
  154. *>
  155. *> \param[out] S
  156. *> \verbatim
  157. *> S is REAL array of dimension N.
  158. *> The singular values of A, ordered so that S(i) >= S(i+1).
  159. *> \endverbatim
  160. *>
  161. *> \param[out] U
  162. *> \verbatim
  163. *> U is COMPLEX array, dimension
  164. *> LDU x M if JOBU = 'A'; see the description of LDU. In this case,
  165. *> on exit, U contains the M left singular vectors.
  166. *> LDU x N if JOBU = 'S', 'U', 'R' ; see the description of LDU. In this
  167. *> case, U contains the leading N or the leading NUMRANK left singular vectors.
  168. *> LDU x N if JOBU = 'F' ; see the description of LDU. In this case U
  169. *> contains N x N unitary matrix that can be used to form the left
  170. *> singular vectors.
  171. *> If JOBU = 'N', U is not referenced.
  172. *> \endverbatim
  173. *>
  174. *> \param[in] LDU
  175. *> \verbatim
  176. *> LDU is INTEGER.
  177. *> The leading dimension of the array U.
  178. *> If JOBU = 'A', 'S', 'U', 'R', LDU >= max(1,M).
  179. *> If JOBU = 'F', LDU >= max(1,N).
  180. *> Otherwise, LDU >= 1.
  181. *> \endverbatim
  182. *>
  183. *> \param[out] V
  184. *> \verbatim
  185. *> V is COMPLEX array, dimension
  186. *> LDV x N if JOBV = 'A', 'V', 'R' or if JOBA = 'E' .
  187. *> If JOBV = 'A', or 'V', V contains the N-by-N unitary matrix V**H;
  188. *> If JOBV = 'R', V contains the first NUMRANK rows of V**H (the right
  189. *> singular vectors, stored rowwise, of the NUMRANK largest singular values).
  190. *> If JOBV = 'N' and JOBA = 'E', V is used as a workspace.
  191. *> If JOBV = 'N', and JOBA.NE.'E', V is not referenced.
  192. *> \endverbatim
  193. *>
  194. *> \param[in] LDV
  195. *> \verbatim
  196. *> LDV is INTEGER
  197. *> The leading dimension of the array V.
  198. *> If JOBV = 'A', 'V', 'R', or JOBA = 'E', LDV >= max(1,N).
  199. *> Otherwise, LDV >= 1.
  200. *> \endverbatim
  201. *>
  202. *> \param[out] NUMRANK
  203. *> \verbatim
  204. *> NUMRANK is INTEGER
  205. *> NUMRANK is the numerical rank first determined after the rank
  206. *> revealing QR factorization, following the strategy specified by the
  207. *> value of JOBA. If JOBV = 'R' and JOBU = 'R', only NUMRANK
  208. *> leading singular values and vectors are then requested in the call
  209. *> of CGESVD. The final value of NUMRANK might be further reduced if
  210. *> some singular values are computed as zeros.
  211. *> \endverbatim
  212. *>
  213. *> \param[out] IWORK
  214. *> \verbatim
  215. *> IWORK is INTEGER array, dimension (max(1, LIWORK)).
  216. *> On exit, IWORK(1:N) contains column pivoting permutation of the
  217. *> rank revealing QR factorization.
  218. *> If JOBP = 'P', IWORK(N+1:N+M-1) contains the indices of the sequence
  219. *> of row swaps used in row pivoting. These can be used to restore the
  220. *> left singular vectors in the case JOBU = 'F'.
  221. *>
  222. *> If LIWORK, LCWORK, or LRWORK = -1, then on exit, if INFO = 0,
  223. *> LIWORK(1) returns the minimal LIWORK.
  224. *> \endverbatim
  225. *>
  226. *> \param[in] LIWORK
  227. *> \verbatim
  228. *> LIWORK is INTEGER
  229. *> The dimension of the array IWORK.
  230. *> LIWORK >= N + M - 1, if JOBP = 'P';
  231. *> LIWORK >= N if JOBP = 'N'.
  232. *>
  233. *> If LIWORK = -1, then a workspace query is assumed; the routine
  234. *> only calculates and returns the optimal and minimal sizes
  235. *> for the CWORK, IWORK, and RWORK arrays, and no error
  236. *> message related to LCWORK is issued by XERBLA.
  237. *> \endverbatim
  238. *>
  239. *> \param[out] CWORK
  240. *> \verbatim
  241. *> CWORK is COMPLEX array, dimension (max(2, LCWORK)), used as a workspace.
  242. *> On exit, if, on entry, LCWORK.NE.-1, CWORK(1:N) contains parameters
  243. *> needed to recover the Q factor from the QR factorization computed by
  244. *> CGEQP3.
  245. *>
  246. *> If LIWORK, LCWORK, or LRWORK = -1, then on exit, if INFO = 0,
  247. *> CWORK(1) returns the optimal LCWORK, and
  248. *> CWORK(2) returns the minimal LCWORK.
  249. *> \endverbatim
  250. *>
  251. *> \param[in,out] LCWORK
  252. *> \verbatim
  253. *> LCWORK is INTEGER
  254. *> The dimension of the array CWORK. It is determined as follows:
  255. *> Let LWQP3 = N+1, LWCON = 2*N, and let
  256. *> LWUNQ = { MAX( N, 1 ), if JOBU = 'R', 'S', or 'U'
  257. *> { MAX( M, 1 ), if JOBU = 'A'
  258. *> LWSVD = MAX( 3*N, 1 )
  259. *> LWLQF = MAX( N/2, 1 ), LWSVD2 = MAX( 3*(N/2), 1 ), LWUNLQ = MAX( N, 1 ),
  260. *> LWQRF = MAX( N/2, 1 ), LWUNQ2 = MAX( N, 1 )
  261. *> Then the minimal value of LCWORK is:
  262. *> = MAX( N + LWQP3, LWSVD ) if only the singular values are needed;
  263. *> = MAX( N + LWQP3, LWCON, LWSVD ) if only the singular values are needed,
  264. *> and a scaled condition estimate requested;
  265. *>
  266. *> = N + MAX( LWQP3, LWSVD, LWUNQ ) if the singular values and the left
  267. *> singular vectors are requested;
  268. *> = N + MAX( LWQP3, LWCON, LWSVD, LWUNQ ) if the singular values and the left
  269. *> singular vectors are requested, and also
  270. *> a scaled condition estimate requested;
  271. *>
  272. *> = N + MAX( LWQP3, LWSVD ) if the singular values and the right
  273. *> singular vectors are requested;
  274. *> = N + MAX( LWQP3, LWCON, LWSVD ) if the singular values and the right
  275. *> singular vectors are requested, and also
  276. *> a scaled condition etimate requested;
  277. *>
  278. *> = N + MAX( LWQP3, LWSVD, LWUNQ ) if the full SVD is requested with JOBV = 'R';
  279. *> independent of JOBR;
  280. *> = N + MAX( LWQP3, LWCON, LWSVD, LWUNQ ) if the full SVD is requested,
  281. *> JOBV = 'R' and, also a scaled condition
  282. *> estimate requested; independent of JOBR;
  283. *> = MAX( N + MAX( LWQP3, LWSVD, LWUNQ ),
  284. *> N + MAX( LWQP3, N/2+LWLQF, N/2+LWSVD2, N/2+LWUNLQ, LWUNQ) ) if the
  285. *> full SVD is requested with JOBV = 'A' or 'V', and
  286. *> JOBR ='N'
  287. *> = MAX( N + MAX( LWQP3, LWCON, LWSVD, LWUNQ ),
  288. *> N + MAX( LWQP3, LWCON, N/2+LWLQF, N/2+LWSVD2, N/2+LWUNLQ, LWUNQ ) )
  289. *> if the full SVD is requested with JOBV = 'A' or 'V', and
  290. *> JOBR ='N', and also a scaled condition number estimate
  291. *> requested.
  292. *> = MAX( N + MAX( LWQP3, LWSVD, LWUNQ ),
  293. *> N + MAX( LWQP3, N/2+LWQRF, N/2+LWSVD2, N/2+LWUNQ2, LWUNQ ) ) if the
  294. *> full SVD is requested with JOBV = 'A', 'V', and JOBR ='T'
  295. *> = MAX( N + MAX( LWQP3, LWCON, LWSVD, LWUNQ ),
  296. *> N + MAX( LWQP3, LWCON, N/2+LWQRF, N/2+LWSVD2, N/2+LWUNQ2, LWUNQ ) )
  297. *> if the full SVD is requested with JOBV = 'A', 'V' and
  298. *> JOBR ='T', and also a scaled condition number estimate
  299. *> requested.
  300. *> Finally, LCWORK must be at least two: LCWORK = MAX( 2, LCWORK ).
  301. *>
  302. *> If LCWORK = -1, then a workspace query is assumed; the routine
  303. *> only calculates and returns the optimal and minimal sizes
  304. *> for the CWORK, IWORK, and RWORK arrays, and no error
  305. *> message related to LCWORK is issued by XERBLA.
  306. *> \endverbatim
  307. *>
  308. *> \param[out] RWORK
  309. *> \verbatim
  310. *> RWORK is REAL array, dimension (max(1, LRWORK)).
  311. *> On exit,
  312. *> 1. If JOBA = 'E', RWORK(1) contains an estimate of the condition
  313. *> number of column scaled A. If A = C * D where D is diagonal and C
  314. *> has unit columns in the Euclidean norm, then, assuming full column rank,
  315. *> N^(-1/4) * RWORK(1) <= ||pinv(C)||_2 <= N^(1/4) * RWORK(1).
  316. *> Otherwise, RWORK(1) = -1.
  317. *> 2. RWORK(2) contains the number of singular values computed as
  318. *> exact zeros in CGESVD applied to the upper triangular or trapeziodal
  319. *> R (from the initial QR factorization). In case of early exit (no call to
  320. *> CGESVD, such as in the case of zero matrix) RWORK(2) = -1.
  321. *>
  322. *> If LIWORK, LCWORK, or LRWORK = -1, then on exit, if INFO = 0,
  323. *> RWORK(1) returns the minimal LRWORK.
  324. *> \endverbatim
  325. *>
  326. *> \param[in] LRWORK
  327. *> \verbatim
  328. *> LRWORK is INTEGER.
  329. *> The dimension of the array RWORK.
  330. *> If JOBP ='P', then LRWORK >= MAX(2, M, 5*N);
  331. *> Otherwise, LRWORK >= MAX(2, 5*N).
  332. *>
  333. *> If LRWORK = -1, then a workspace query is assumed; the routine
  334. *> only calculates and returns the optimal and minimal sizes
  335. *> for the CWORK, IWORK, and RWORK arrays, and no error
  336. *> message related to LCWORK is issued by XERBLA.
  337. *> \endverbatim
  338. *>
  339. *> \param[out] INFO
  340. *> \verbatim
  341. *> INFO is INTEGER
  342. *> = 0: successful exit.
  343. *> < 0: if INFO = -i, the i-th argument had an illegal value.
  344. *> > 0: if CBDSQR did not converge, INFO specifies how many superdiagonals
  345. *> of an intermediate bidiagonal form B (computed in CGESVD) did not
  346. *> converge to zero.
  347. *> \endverbatim
  348. *
  349. *> \par Further Details:
  350. * ========================
  351. *>
  352. *> \verbatim
  353. *>
  354. *> 1. The data movement (matrix transpose) is coded using simple nested
  355. *> DO-loops because BLAS and LAPACK do not provide corresponding subroutines.
  356. *> Those DO-loops are easily identified in this source code - by the CONTINUE
  357. *> statements labeled with 11**. In an optimized version of this code, the
  358. *> nested DO loops should be replaced with calls to an optimized subroutine.
  359. *> 2. This code scales A by 1/SQRT(M) if the largest ABS(A(i,j)) could cause
  360. *> column norm overflow. This is the minial precaution and it is left to the
  361. *> SVD routine (CGESVD) to do its own preemptive scaling if potential over-
  362. *> or underflows are detected. To avoid repeated scanning of the array A,
  363. *> an optimal implementation would do all necessary scaling before calling
  364. *> CGESVD and the scaling in CGESVD can be switched off.
  365. *> 3. Other comments related to code optimization are given in comments in the
  366. *> code, enlosed in [[double brackets]].
  367. *> \endverbatim
  368. *
  369. *> \par Bugs, examples and comments
  370. * ===========================
  371. *
  372. *> \verbatim
  373. *> Please report all bugs and send interesting examples and/or comments to
  374. *> drmac@math.hr. Thank you.
  375. *> \endverbatim
  376. *
  377. *> \par References
  378. * ===============
  379. *
  380. *> \verbatim
  381. *> [1] Zlatko Drmac, Algorithm 977: A QR-Preconditioned QR SVD Method for
  382. *> Computing the SVD with High Accuracy. ACM Trans. Math. Softw.
  383. *> 44(1): 11:1-11:30 (2017)
  384. *>
  385. *> SIGMA library, xGESVDQ section updated February 2016.
  386. *> Developed and coded by Zlatko Drmac, Department of Mathematics
  387. *> University of Zagreb, Croatia, drmac@math.hr
  388. *> \endverbatim
  389. *
  390. *
  391. *> \par Contributors:
  392. * ==================
  393. *>
  394. *> \verbatim
  395. *> Developed and coded by Zlatko Drmac, Department of Mathematics
  396. *> University of Zagreb, Croatia, drmac@math.hr
  397. *> \endverbatim
  398. *
  399. * Authors:
  400. * ========
  401. *
  402. *> \author Univ. of Tennessee
  403. *> \author Univ. of California Berkeley
  404. *> \author Univ. of Colorado Denver
  405. *> \author NAG Ltd.
  406. *
  407. *> \date November 2018
  408. *
  409. *> \ingroup complexGEsing
  410. *
  411. * =====================================================================
  412. SUBROUTINE CGESVDQ( JOBA, JOBP, JOBR, JOBU, JOBV, M, N, A, LDA,
  413. $ S, U, LDU, V, LDV, NUMRANK, IWORK, LIWORK,
  414. $ CWORK, LCWORK, RWORK, LRWORK, INFO )
  415. * .. Scalar Arguments ..
  416. IMPLICIT NONE
  417. CHARACTER JOBA, JOBP, JOBR, JOBU, JOBV
  418. INTEGER M, N, LDA, LDU, LDV, NUMRANK, LIWORK, LCWORK, LRWORK,
  419. $ INFO
  420. * ..
  421. * .. Array Arguments ..
  422. COMPLEX A( LDA, * ), U( LDU, * ), V( LDV, * ), CWORK( * )
  423. REAL S( * ), RWORK( * )
  424. INTEGER IWORK( * )
  425. *
  426. * =====================================================================
  427. *
  428. * .. Parameters ..
  429. REAL ZERO, ONE
  430. PARAMETER ( ZERO = 0.0E0, ONE = 1.0E0 )
  431. COMPLEX CZERO, CONE
  432. PARAMETER ( CZERO = ( 0.0E0, 0.0E0 ), CONE = ( 1.0E0, 0.0E0 ) )
  433. * ..
  434. * .. Local Scalars ..
  435. INTEGER IERR, NR, N1, OPTRATIO, p, q
  436. INTEGER LWCON, LWQP3, LWRK_CGELQF, LWRK_CGESVD, LWRK_CGESVD2,
  437. $ LWRK_CGEQP3, LWRK_CGEQRF, LWRK_CUNMLQ, LWRK_CUNMQR,
  438. $ LWRK_CUNMQR2, LWLQF, LWQRF, LWSVD, LWSVD2, LWUNQ,
  439. $ LWUNQ2, LWUNLQ, MINWRK, MINWRK2, OPTWRK, OPTWRK2,
  440. $ IMINWRK, RMINWRK
  441. LOGICAL ACCLA, ACCLM, ACCLH, ASCALED, CONDA, DNTWU, DNTWV,
  442. $ LQUERY, LSVC0, LSVEC, ROWPRM, RSVEC, RTRANS, WNTUA,
  443. $ WNTUF, WNTUR, WNTUS, WNTVA, WNTVR
  444. REAL BIG, EPSLN, RTMP, SCONDA, SFMIN
  445. COMPLEX CTMP
  446. * ..
  447. * .. Local Arrays
  448. COMPLEX CDUMMY(1)
  449. REAL RDUMMY(1)
  450. * ..
  451. * .. External Subroutines (BLAS, LAPACK)
  452. EXTERNAL CGELQF, CGEQP3, CGEQRF, CGESVD, CLACPY, CLAPMT,
  453. $ CLASCL, CLASET, CLASWP, CSSCAL, SLASET, SLASCL,
  454. $ CPOCON, CUNMLQ, CUNMQR, XERBLA
  455. * ..
  456. * .. External Functions (BLAS, LAPACK)
  457. LOGICAL LSAME
  458. INTEGER ISAMAX
  459. REAL CLANGE, SCNRM2, SLAMCH
  460. EXTERNAL CLANGE, LSAME, ISAMAX, SCNRM2, SLAMCH
  461. * ..
  462. * .. Intrinsic Functions ..
  463. INTRINSIC ABS, CONJG, MAX, MIN, REAL, SQRT
  464. * ..
  465. * .. Executable Statements ..
  466. *
  467. * Test the input arguments
  468. *
  469. WNTUS = LSAME( JOBU, 'S' ) .OR. LSAME( JOBU, 'U' )
  470. WNTUR = LSAME( JOBU, 'R' )
  471. WNTUA = LSAME( JOBU, 'A' )
  472. WNTUF = LSAME( JOBU, 'F' )
  473. LSVC0 = WNTUS .OR. WNTUR .OR. WNTUA
  474. LSVEC = LSVC0 .OR. WNTUF
  475. DNTWU = LSAME( JOBU, 'N' )
  476. *
  477. WNTVR = LSAME( JOBV, 'R' )
  478. WNTVA = LSAME( JOBV, 'A' ) .OR. LSAME( JOBV, 'V' )
  479. RSVEC = WNTVR .OR. WNTVA
  480. DNTWV = LSAME( JOBV, 'N' )
  481. *
  482. ACCLA = LSAME( JOBA, 'A' )
  483. ACCLM = LSAME( JOBA, 'M' )
  484. CONDA = LSAME( JOBA, 'E' )
  485. ACCLH = LSAME( JOBA, 'H' ) .OR. CONDA
  486. *
  487. ROWPRM = LSAME( JOBP, 'P' )
  488. RTRANS = LSAME( JOBR, 'T' )
  489. *
  490. IF ( ROWPRM ) THEN
  491. IMINWRK = MAX( 1, N + M - 1 )
  492. RMINWRK = MAX( 2, M, 5*N )
  493. ELSE
  494. IMINWRK = MAX( 1, N )
  495. RMINWRK = MAX( 2, 5*N )
  496. END IF
  497. LQUERY = (LIWORK .EQ. -1 .OR. LCWORK .EQ. -1 .OR. LRWORK .EQ. -1)
  498. INFO = 0
  499. IF ( .NOT. ( ACCLA .OR. ACCLM .OR. ACCLH ) ) THEN
  500. INFO = -1
  501. ELSE IF ( .NOT.( ROWPRM .OR. LSAME( JOBP, 'N' ) ) ) THEN
  502. INFO = -2
  503. ELSE IF ( .NOT.( RTRANS .OR. LSAME( JOBR, 'N' ) ) ) THEN
  504. INFO = -3
  505. ELSE IF ( .NOT.( LSVEC .OR. DNTWU ) ) THEN
  506. INFO = -4
  507. ELSE IF ( WNTUR .AND. WNTVA ) THEN
  508. INFO = -5
  509. ELSE IF ( .NOT.( RSVEC .OR. DNTWV )) THEN
  510. INFO = -5
  511. ELSE IF ( M.LT.0 ) THEN
  512. INFO = -6
  513. ELSE IF ( ( N.LT.0 ) .OR. ( N.GT.M ) ) THEN
  514. INFO = -7
  515. ELSE IF ( LDA.LT.MAX( 1, M ) ) THEN
  516. INFO = -9
  517. ELSE IF ( LDU.LT.1 .OR. ( LSVC0 .AND. LDU.LT.M ) .OR.
  518. $ ( WNTUF .AND. LDU.LT.N ) ) THEN
  519. INFO = -12
  520. ELSE IF ( LDV.LT.1 .OR. ( RSVEC .AND. LDV.LT.N ) .OR.
  521. $ ( CONDA .AND. LDV.LT.N ) ) THEN
  522. INFO = -14
  523. ELSE IF ( LIWORK .LT. IMINWRK .AND. .NOT. LQUERY ) THEN
  524. INFO = -17
  525. END IF
  526. *
  527. *
  528. IF ( INFO .EQ. 0 ) THEN
  529. *
  530. * Compute workspace
  531. * .. compute the minimal and the optimal workspace lengths
  532. * [[The expressions for computing the minimal and the optimal
  533. * values of LCWORK are written with a lot of redundancy and
  534. * can be simplified. However, this detailed form is easier for
  535. * maintenance and modifications of the code.]]
  536. *
  537. * .. minimal workspace length for CGEQP3 of an M x N matrix
  538. LWQP3 = N+1
  539. * .. minimal workspace length for CUNMQR to build left singular vectors
  540. IF ( WNTUS .OR. WNTUR ) THEN
  541. LWUNQ = MAX( N , 1 )
  542. ELSE IF ( WNTUA ) THEN
  543. LWUNQ = MAX( M , 1 )
  544. END IF
  545. * .. minimal workspace length for CPOCON of an N x N matrix
  546. LWCON = 2 * N
  547. * .. CGESVD of an N x N matrix
  548. LWSVD = MAX( 3 * N, 1 )
  549. IF ( LQUERY ) THEN
  550. CALL CGEQP3( M, N, A, LDA, IWORK, CDUMMY, CDUMMY, -1,
  551. $ RDUMMY, IERR )
  552. LWRK_CGEQP3 = INT( CDUMMY(1) )
  553. IF ( WNTUS .OR. WNTUR ) THEN
  554. CALL CUNMQR( 'L', 'N', M, N, N, A, LDA, CDUMMY, U,
  555. $ LDU, CDUMMY, -1, IERR )
  556. LWRK_CUNMQR = INT( CDUMMY(1) )
  557. ELSE IF ( WNTUA ) THEN
  558. CALL CUNMQR( 'L', 'N', M, M, N, A, LDA, CDUMMY, U,
  559. $ LDU, CDUMMY, -1, IERR )
  560. LWRK_CUNMQR = INT( CDUMMY(1) )
  561. ELSE
  562. LWRK_CUNMQR = 0
  563. END IF
  564. END IF
  565. MINWRK = 2
  566. OPTWRK = 2
  567. IF ( .NOT. (LSVEC .OR. RSVEC )) THEN
  568. * .. minimal and optimal sizes of the complex workspace if
  569. * only the singular values are requested
  570. IF ( CONDA ) THEN
  571. MINWRK = MAX( N+LWQP3, LWCON, LWSVD )
  572. ELSE
  573. MINWRK = MAX( N+LWQP3, LWSVD )
  574. END IF
  575. IF ( LQUERY ) THEN
  576. CALL CGESVD( 'N', 'N', N, N, A, LDA, S, U, LDU,
  577. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  578. LWRK_CGESVD = INT( CDUMMY(1) )
  579. IF ( CONDA ) THEN
  580. OPTWRK = MAX( N+LWRK_CGEQP3, N+LWCON, LWRK_CGESVD )
  581. ELSE
  582. OPTWRK = MAX( N+LWRK_CGEQP3, LWRK_CGESVD )
  583. END IF
  584. END IF
  585. ELSE IF ( LSVEC .AND. (.NOT.RSVEC) ) THEN
  586. * .. minimal and optimal sizes of the complex workspace if the
  587. * singular values and the left singular vectors are requested
  588. IF ( CONDA ) THEN
  589. MINWRK = N + MAX( LWQP3, LWCON, LWSVD, LWUNQ )
  590. ELSE
  591. MINWRK = N + MAX( LWQP3, LWSVD, LWUNQ )
  592. END IF
  593. IF ( LQUERY ) THEN
  594. IF ( RTRANS ) THEN
  595. CALL CGESVD( 'N', 'O', N, N, A, LDA, S, U, LDU,
  596. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  597. ELSE
  598. CALL CGESVD( 'O', 'N', N, N, A, LDA, S, U, LDU,
  599. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  600. END IF
  601. LWRK_CGESVD = INT( CDUMMY(1) )
  602. IF ( CONDA ) THEN
  603. OPTWRK = N + MAX( LWRK_CGEQP3, LWCON, LWRK_CGESVD,
  604. $ LWRK_CUNMQR )
  605. ELSE
  606. OPTWRK = N + MAX( LWRK_CGEQP3, LWRK_CGESVD,
  607. $ LWRK_CUNMQR )
  608. END IF
  609. END IF
  610. ELSE IF ( RSVEC .AND. (.NOT.LSVEC) ) THEN
  611. * .. minimal and optimal sizes of the complex workspace if the
  612. * singular values and the right singular vectors are requested
  613. IF ( CONDA ) THEN
  614. MINWRK = N + MAX( LWQP3, LWCON, LWSVD )
  615. ELSE
  616. MINWRK = N + MAX( LWQP3, LWSVD )
  617. END IF
  618. IF ( LQUERY ) THEN
  619. IF ( RTRANS ) THEN
  620. CALL CGESVD( 'O', 'N', N, N, A, LDA, S, U, LDU,
  621. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  622. ELSE
  623. CALL CGESVD( 'N', 'O', N, N, A, LDA, S, U, LDU,
  624. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  625. END IF
  626. LWRK_CGESVD = INT( CDUMMY(1) )
  627. IF ( CONDA ) THEN
  628. OPTWRK = N + MAX( LWRK_CGEQP3, LWCON, LWRK_CGESVD )
  629. ELSE
  630. OPTWRK = N + MAX( LWRK_CGEQP3, LWRK_CGESVD )
  631. END IF
  632. END IF
  633. ELSE
  634. * .. minimal and optimal sizes of the complex workspace if the
  635. * full SVD is requested
  636. IF ( RTRANS ) THEN
  637. MINWRK = MAX( LWQP3, LWSVD, LWUNQ )
  638. IF ( CONDA ) MINWRK = MAX( MINWRK, LWCON )
  639. MINWRK = MINWRK + N
  640. IF ( WNTVA ) THEN
  641. * .. minimal workspace length for N x N/2 CGEQRF
  642. LWQRF = MAX( N/2, 1 )
  643. * .. minimal workspace lengt for N/2 x N/2 CGESVD
  644. LWSVD2 = MAX( 3 * (N/2), 1 )
  645. LWUNQ2 = MAX( N, 1 )
  646. MINWRK2 = MAX( LWQP3, N/2+LWQRF, N/2+LWSVD2,
  647. $ N/2+LWUNQ2, LWUNQ )
  648. IF ( CONDA ) MINWRK2 = MAX( MINWRK2, LWCON )
  649. MINWRK2 = N + MINWRK2
  650. MINWRK = MAX( MINWRK, MINWRK2 )
  651. END IF
  652. ELSE
  653. MINWRK = MAX( LWQP3, LWSVD, LWUNQ )
  654. IF ( CONDA ) MINWRK = MAX( MINWRK, LWCON )
  655. MINWRK = MINWRK + N
  656. IF ( WNTVA ) THEN
  657. * .. minimal workspace length for N/2 x N CGELQF
  658. LWLQF = MAX( N/2, 1 )
  659. LWSVD2 = MAX( 3 * (N/2), 1 )
  660. LWUNLQ = MAX( N , 1 )
  661. MINWRK2 = MAX( LWQP3, N/2+LWLQF, N/2+LWSVD2,
  662. $ N/2+LWUNLQ, LWUNQ )
  663. IF ( CONDA ) MINWRK2 = MAX( MINWRK2, LWCON )
  664. MINWRK2 = N + MINWRK2
  665. MINWRK = MAX( MINWRK, MINWRK2 )
  666. END IF
  667. END IF
  668. IF ( LQUERY ) THEN
  669. IF ( RTRANS ) THEN
  670. CALL CGESVD( 'O', 'A', N, N, A, LDA, S, U, LDU,
  671. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  672. LWRK_CGESVD = INT( CDUMMY(1) )
  673. OPTWRK = MAX(LWRK_CGEQP3,LWRK_CGESVD,LWRK_CUNMQR)
  674. IF ( CONDA ) OPTWRK = MAX( OPTWRK, LWCON )
  675. OPTWRK = N + OPTWRK
  676. IF ( WNTVA ) THEN
  677. CALL CGEQRF(N,N/2,U,LDU,CDUMMY,CDUMMY,-1,IERR)
  678. LWRK_CGEQRF = INT( CDUMMY(1) )
  679. CALL CGESVD( 'S', 'O', N/2,N/2, V,LDV, S, U,LDU,
  680. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  681. LWRK_CGESVD2 = INT( CDUMMY(1) )
  682. CALL CUNMQR( 'R', 'C', N, N, N/2, U, LDU, CDUMMY,
  683. $ V, LDV, CDUMMY, -1, IERR )
  684. LWRK_CUNMQR2 = INT( CDUMMY(1) )
  685. OPTWRK2 = MAX( LWRK_CGEQP3, N/2+LWRK_CGEQRF,
  686. $ N/2+LWRK_CGESVD2, N/2+LWRK_CUNMQR2 )
  687. IF ( CONDA ) OPTWRK2 = MAX( OPTWRK2, LWCON )
  688. OPTWRK2 = N + OPTWRK2
  689. OPTWRK = MAX( OPTWRK, OPTWRK2 )
  690. END IF
  691. ELSE
  692. CALL CGESVD( 'S', 'O', N, N, A, LDA, S, U, LDU,
  693. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  694. LWRK_CGESVD = INT( CDUMMY(1) )
  695. OPTWRK = MAX(LWRK_CGEQP3,LWRK_CGESVD,LWRK_CUNMQR)
  696. IF ( CONDA ) OPTWRK = MAX( OPTWRK, LWCON )
  697. OPTWRK = N + OPTWRK
  698. IF ( WNTVA ) THEN
  699. CALL CGELQF(N/2,N,U,LDU,CDUMMY,CDUMMY,-1,IERR)
  700. LWRK_CGELQF = INT( CDUMMY(1) )
  701. CALL CGESVD( 'S','O', N/2,N/2, V, LDV, S, U, LDU,
  702. $ V, LDV, CDUMMY, -1, RDUMMY, IERR )
  703. LWRK_CGESVD2 = INT( CDUMMY(1) )
  704. CALL CUNMLQ( 'R', 'N', N, N, N/2, U, LDU, CDUMMY,
  705. $ V, LDV, CDUMMY,-1,IERR )
  706. LWRK_CUNMLQ = INT( CDUMMY(1) )
  707. OPTWRK2 = MAX( LWRK_CGEQP3, N/2+LWRK_CGELQF,
  708. $ N/2+LWRK_CGESVD2, N/2+LWRK_CUNMLQ )
  709. IF ( CONDA ) OPTWRK2 = MAX( OPTWRK2, LWCON )
  710. OPTWRK2 = N + OPTWRK2
  711. OPTWRK = MAX( OPTWRK, OPTWRK2 )
  712. END IF
  713. END IF
  714. END IF
  715. END IF
  716. *
  717. MINWRK = MAX( 2, MINWRK )
  718. OPTWRK = MAX( 2, OPTWRK )
  719. IF ( LCWORK .LT. MINWRK .AND. (.NOT.LQUERY) ) INFO = -19
  720. *
  721. END IF
  722. *
  723. IF (INFO .EQ. 0 .AND. LRWORK .LT. RMINWRK .AND. .NOT. LQUERY) THEN
  724. INFO = -21
  725. END IF
  726. IF( INFO.NE.0 ) THEN
  727. CALL XERBLA( 'CGESVDQ', -INFO )
  728. RETURN
  729. ELSE IF ( LQUERY ) THEN
  730. *
  731. * Return optimal workspace
  732. *
  733. IWORK(1) = IMINWRK
  734. CWORK(1) = OPTWRK
  735. CWORK(2) = MINWRK
  736. RWORK(1) = RMINWRK
  737. RETURN
  738. END IF
  739. *
  740. * Quick return if the matrix is void.
  741. *
  742. IF( ( M.EQ.0 ) .OR. ( N.EQ.0 ) ) THEN
  743. * .. all output is void.
  744. RETURN
  745. END IF
  746. *
  747. BIG = SLAMCH('O')
  748. ASCALED = .FALSE.
  749. IF ( ROWPRM ) THEN
  750. * .. reordering the rows in decreasing sequence in the
  751. * ell-infinity norm - this enhances numerical robustness in
  752. * the case of differently scaled rows.
  753. DO 1904 p = 1, M
  754. * RWORK(p) = ABS( A(p,ICAMAX(N,A(p,1),LDA)) )
  755. * [[CLANGE will return NaN if an entry of the p-th row is Nan]]
  756. RWORK(p) = CLANGE( 'M', 1, N, A(p,1), LDA, RDUMMY )
  757. * .. check for NaN's and Inf's
  758. IF ( ( RWORK(p) .NE. RWORK(p) ) .OR.
  759. $ ( (RWORK(p)*ZERO) .NE. ZERO ) ) THEN
  760. INFO = - 8
  761. CALL XERBLA( 'CGESVDQ', -INFO )
  762. RETURN
  763. END IF
  764. 1904 CONTINUE
  765. DO 1952 p = 1, M - 1
  766. q = ISAMAX( M-p+1, RWORK(p), 1 ) + p - 1
  767. IWORK(N+p) = q
  768. IF ( p .NE. q ) THEN
  769. RTMP = RWORK(p)
  770. RWORK(p) = RWORK(q)
  771. RWORK(q) = RTMP
  772. END IF
  773. 1952 CONTINUE
  774. *
  775. IF ( RWORK(1) .EQ. ZERO ) THEN
  776. * Quick return: A is the M x N zero matrix.
  777. NUMRANK = 0
  778. CALL SLASET( 'G', N, 1, ZERO, ZERO, S, N )
  779. IF ( WNTUS ) CALL CLASET('G', M, N, CZERO, CONE, U, LDU)
  780. IF ( WNTUA ) CALL CLASET('G', M, M, CZERO, CONE, U, LDU)
  781. IF ( WNTVA ) CALL CLASET('G', N, N, CZERO, CONE, V, LDV)
  782. IF ( WNTUF ) THEN
  783. CALL CLASET( 'G', N, 1, CZERO, CZERO, CWORK, N )
  784. CALL CLASET( 'G', M, N, CZERO, CONE, U, LDU )
  785. END IF
  786. DO 5001 p = 1, N
  787. IWORK(p) = p
  788. 5001 CONTINUE
  789. IF ( ROWPRM ) THEN
  790. DO 5002 p = N + 1, N + M - 1
  791. IWORK(p) = p - N
  792. 5002 CONTINUE
  793. END IF
  794. IF ( CONDA ) RWORK(1) = -1
  795. RWORK(2) = -1
  796. RETURN
  797. END IF
  798. *
  799. IF ( RWORK(1) .GT. BIG / SQRT(REAL(M)) ) THEN
  800. * .. to prevent overflow in the QR factorization, scale the
  801. * matrix by 1/sqrt(M) if too large entry detected
  802. CALL CLASCL('G',0,0,SQRT(REAL(M)),ONE, M,N, A,LDA, IERR)
  803. ASCALED = .TRUE.
  804. END IF
  805. CALL CLASWP( N, A, LDA, 1, M-1, IWORK(N+1), 1 )
  806. END IF
  807. *
  808. * .. At this stage, preemptive scaling is done only to avoid column
  809. * norms overflows during the QR factorization. The SVD procedure should
  810. * have its own scaling to save the singular values from overflows and
  811. * underflows. That depends on the SVD procedure.
  812. *
  813. IF ( .NOT.ROWPRM ) THEN
  814. RTMP = CLANGE( 'M', M, N, A, LDA, RWORK )
  815. IF ( ( RTMP .NE. RTMP ) .OR.
  816. $ ( (RTMP*ZERO) .NE. ZERO ) ) THEN
  817. INFO = - 8
  818. CALL XERBLA( 'CGESVDQ', -INFO )
  819. RETURN
  820. END IF
  821. IF ( RTMP .GT. BIG / SQRT(REAL(M)) ) THEN
  822. * .. to prevent overflow in the QR factorization, scale the
  823. * matrix by 1/sqrt(M) if too large entry detected
  824. CALL CLASCL('G',0,0, SQRT(REAL(M)),ONE, M,N, A,LDA, IERR)
  825. ASCALED = .TRUE.
  826. END IF
  827. END IF
  828. *
  829. * .. QR factorization with column pivoting
  830. *
  831. * A * P = Q * [ R ]
  832. * [ 0 ]
  833. *
  834. DO 1963 p = 1, N
  835. * .. all columns are free columns
  836. IWORK(p) = 0
  837. 1963 CONTINUE
  838. CALL CGEQP3( M, N, A, LDA, IWORK, CWORK, CWORK(N+1), LCWORK-N,
  839. $ RWORK, IERR )
  840. *
  841. * If the user requested accuracy level allows truncation in the
  842. * computed upper triangular factor, the matrix R is examined and,
  843. * if possible, replaced with its leading upper trapezoidal part.
  844. *
  845. EPSLN = SLAMCH('E')
  846. SFMIN = SLAMCH('S')
  847. * SMALL = SFMIN / EPSLN
  848. NR = N
  849. *
  850. IF ( ACCLA ) THEN
  851. *
  852. * Standard absolute error bound suffices. All sigma_i with
  853. * sigma_i < N*EPS*||A||_F are flushed to zero. This is an
  854. * aggressive enforcement of lower numerical rank by introducing a
  855. * backward error of the order of N*EPS*||A||_F.
  856. NR = 1
  857. RTMP = SQRT(REAL(N))*EPSLN
  858. DO 3001 p = 2, N
  859. IF ( ABS(A(p,p)) .LT. (RTMP*ABS(A(1,1))) ) GO TO 3002
  860. NR = NR + 1
  861. 3001 CONTINUE
  862. 3002 CONTINUE
  863. *
  864. ELSEIF ( ACCLM ) THEN
  865. * .. similarly as above, only slightly more gentle (less aggressive).
  866. * Sudden drop on the diagonal of R is used as the criterion for being
  867. * close-to-rank-deficient. The threshold is set to EPSLN=SLAMCH('E').
  868. * [[This can be made more flexible by replacing this hard-coded value
  869. * with a user specified threshold.]] Also, the values that underflow
  870. * will be truncated.
  871. NR = 1
  872. DO 3401 p = 2, N
  873. IF ( ( ABS(A(p,p)) .LT. (EPSLN*ABS(A(p-1,p-1))) ) .OR.
  874. $ ( ABS(A(p,p)) .LT. SFMIN ) ) GO TO 3402
  875. NR = NR + 1
  876. 3401 CONTINUE
  877. 3402 CONTINUE
  878. *
  879. ELSE
  880. * .. RRQR not authorized to determine numerical rank except in the
  881. * obvious case of zero pivots.
  882. * .. inspect R for exact zeros on the diagonal;
  883. * R(i,i)=0 => R(i:N,i:N)=0.
  884. NR = 1
  885. DO 3501 p = 2, N
  886. IF ( ABS(A(p,p)) .EQ. ZERO ) GO TO 3502
  887. NR = NR + 1
  888. 3501 CONTINUE
  889. 3502 CONTINUE
  890. *
  891. IF ( CONDA ) THEN
  892. * Estimate the scaled condition number of A. Use the fact that it is
  893. * the same as the scaled condition number of R.
  894. * .. V is used as workspace
  895. CALL CLACPY( 'U', N, N, A, LDA, V, LDV )
  896. * Only the leading NR x NR submatrix of the triangular factor
  897. * is considered. Only if NR=N will this give a reliable error
  898. * bound. However, even for NR < N, this can be used on an
  899. * expert level and obtain useful information in the sense of
  900. * perturbation theory.
  901. DO 3053 p = 1, NR
  902. RTMP = SCNRM2( p, V(1,p), 1 )
  903. CALL CSSCAL( p, ONE/RTMP, V(1,p), 1 )
  904. 3053 CONTINUE
  905. IF ( .NOT. ( LSVEC .OR. RSVEC ) ) THEN
  906. CALL CPOCON( 'U', NR, V, LDV, ONE, RTMP,
  907. $ CWORK, RWORK, IERR )
  908. ELSE
  909. CALL CPOCON( 'U', NR, V, LDV, ONE, RTMP,
  910. $ CWORK(N+1), RWORK, IERR )
  911. END IF
  912. SCONDA = ONE / SQRT(RTMP)
  913. * For NR=N, SCONDA is an estimate of SQRT(||(R^* * R)^(-1)||_1),
  914. * N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
  915. * See the reference [1] for more details.
  916. END IF
  917. *
  918. ENDIF
  919. *
  920. IF ( WNTUR ) THEN
  921. N1 = NR
  922. ELSE IF ( WNTUS .OR. WNTUF) THEN
  923. N1 = N
  924. ELSE IF ( WNTUA ) THEN
  925. N1 = M
  926. END IF
  927. *
  928. IF ( .NOT. ( RSVEC .OR. LSVEC ) ) THEN
  929. *.......................................................................
  930. * .. only the singular values are requested
  931. *.......................................................................
  932. IF ( RTRANS ) THEN
  933. *
  934. * .. compute the singular values of R**H = [A](1:NR,1:N)**H
  935. * .. set the lower triangle of [A] to [A](1:NR,1:N)**H and
  936. * the upper triangle of [A] to zero.
  937. DO 1146 p = 1, MIN( N, NR )
  938. A(p,p) = CONJG(A(p,p))
  939. DO 1147 q = p + 1, N
  940. A(q,p) = CONJG(A(p,q))
  941. IF ( q .LE. NR ) A(p,q) = CZERO
  942. 1147 CONTINUE
  943. 1146 CONTINUE
  944. *
  945. CALL CGESVD( 'N', 'N', N, NR, A, LDA, S, U, LDU,
  946. $ V, LDV, CWORK, LCWORK, RWORK, INFO )
  947. *
  948. ELSE
  949. *
  950. * .. compute the singular values of R = [A](1:NR,1:N)
  951. *
  952. IF ( NR .GT. 1 )
  953. $ CALL CLASET( 'L', NR-1,NR-1, CZERO,CZERO, A(2,1), LDA )
  954. CALL CGESVD( 'N', 'N', NR, N, A, LDA, S, U, LDU,
  955. $ V, LDV, CWORK, LCWORK, RWORK, INFO )
  956. *
  957. END IF
  958. *
  959. ELSE IF ( LSVEC .AND. ( .NOT. RSVEC) ) THEN
  960. *.......................................................................
  961. * .. the singular values and the left singular vectors requested
  962. *.......................................................................""""""""
  963. IF ( RTRANS ) THEN
  964. * .. apply CGESVD to R**H
  965. * .. copy R**H into [U] and overwrite [U] with the right singular
  966. * vectors of R
  967. DO 1192 p = 1, NR
  968. DO 1193 q = p, N
  969. U(q,p) = CONJG(A(p,q))
  970. 1193 CONTINUE
  971. 1192 CONTINUE
  972. IF ( NR .GT. 1 )
  973. $ CALL CLASET( 'U', NR-1,NR-1, CZERO,CZERO, U(1,2), LDU )
  974. * .. the left singular vectors not computed, the NR right singular
  975. * vectors overwrite [U](1:NR,1:NR) as conjugate transposed. These
  976. * will be pre-multiplied by Q to build the left singular vectors of A.
  977. CALL CGESVD( 'N', 'O', N, NR, U, LDU, S, U, LDU,
  978. $ U, LDU, CWORK(N+1), LCWORK-N, RWORK, INFO )
  979. *
  980. DO 1119 p = 1, NR
  981. U(p,p) = CONJG(U(p,p))
  982. DO 1120 q = p + 1, NR
  983. CTMP = CONJG(U(q,p))
  984. U(q,p) = CONJG(U(p,q))
  985. U(p,q) = CTMP
  986. 1120 CONTINUE
  987. 1119 CONTINUE
  988. *
  989. ELSE
  990. * .. apply CGESVD to R
  991. * .. copy R into [U] and overwrite [U] with the left singular vectors
  992. CALL CLACPY( 'U', NR, N, A, LDA, U, LDU )
  993. IF ( NR .GT. 1 )
  994. $ CALL CLASET( 'L', NR-1, NR-1, CZERO, CZERO, U(2,1), LDU )
  995. * .. the right singular vectors not computed, the NR left singular
  996. * vectors overwrite [U](1:NR,1:NR)
  997. CALL CGESVD( 'O', 'N', NR, N, U, LDU, S, U, LDU,
  998. $ V, LDV, CWORK(N+1), LCWORK-N, RWORK, INFO )
  999. * .. now [U](1:NR,1:NR) contains the NR left singular vectors of
  1000. * R. These will be pre-multiplied by Q to build the left singular
  1001. * vectors of A.
  1002. END IF
  1003. *
  1004. * .. assemble the left singular vector matrix U of dimensions
  1005. * (M x NR) or (M x N) or (M x M).
  1006. IF ( ( NR .LT. M ) .AND. ( .NOT.WNTUF ) ) THEN
  1007. CALL CLASET('A', M-NR, NR, CZERO, CZERO, U(NR+1,1), LDU)
  1008. IF ( NR .LT. N1 ) THEN
  1009. CALL CLASET( 'A',NR,N1-NR,CZERO,CZERO,U(1,NR+1), LDU )
  1010. CALL CLASET( 'A',M-NR,N1-NR,CZERO,CONE,
  1011. $ U(NR+1,NR+1), LDU )
  1012. END IF
  1013. END IF
  1014. *
  1015. * The Q matrix from the first QRF is built into the left singular
  1016. * vectors matrix U.
  1017. *
  1018. IF ( .NOT.WNTUF )
  1019. $ CALL CUNMQR( 'L', 'N', M, N1, N, A, LDA, CWORK, U,
  1020. $ LDU, CWORK(N+1), LCWORK-N, IERR )
  1021. IF ( ROWPRM .AND. .NOT.WNTUF )
  1022. $ CALL CLASWP( N1, U, LDU, 1, M-1, IWORK(N+1), -1 )
  1023. *
  1024. ELSE IF ( RSVEC .AND. ( .NOT. LSVEC ) ) THEN
  1025. *.......................................................................
  1026. * .. the singular values and the right singular vectors requested
  1027. *.......................................................................
  1028. IF ( RTRANS ) THEN
  1029. * .. apply CGESVD to R**H
  1030. * .. copy R**H into V and overwrite V with the left singular vectors
  1031. DO 1165 p = 1, NR
  1032. DO 1166 q = p, N
  1033. V(q,p) = CONJG(A(p,q))
  1034. 1166 CONTINUE
  1035. 1165 CONTINUE
  1036. IF ( NR .GT. 1 )
  1037. $ CALL CLASET( 'U', NR-1,NR-1, CZERO,CZERO, V(1,2), LDV )
  1038. * .. the left singular vectors of R**H overwrite V, the right singular
  1039. * vectors not computed
  1040. IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
  1041. CALL CGESVD( 'O', 'N', N, NR, V, LDV, S, U, LDU,
  1042. $ U, LDU, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1043. *
  1044. DO 1121 p = 1, NR
  1045. V(p,p) = CONJG(V(p,p))
  1046. DO 1122 q = p + 1, NR
  1047. CTMP = CONJG(V(q,p))
  1048. V(q,p) = CONJG(V(p,q))
  1049. V(p,q) = CTMP
  1050. 1122 CONTINUE
  1051. 1121 CONTINUE
  1052. *
  1053. IF ( NR .LT. N ) THEN
  1054. DO 1103 p = 1, NR
  1055. DO 1104 q = NR + 1, N
  1056. V(p,q) = CONJG(V(q,p))
  1057. 1104 CONTINUE
  1058. 1103 CONTINUE
  1059. END IF
  1060. CALL CLAPMT( .FALSE., NR, N, V, LDV, IWORK )
  1061. ELSE
  1062. * .. need all N right singular vectors and NR < N
  1063. * [!] This is simple implementation that augments [V](1:N,1:NR)
  1064. * by padding a zero block. In the case NR << N, a more efficient
  1065. * way is to first use the QR factorization. For more details
  1066. * how to implement this, see the " FULL SVD " branch.
  1067. CALL CLASET('G', N, N-NR, CZERO, CZERO, V(1,NR+1), LDV)
  1068. CALL CGESVD( 'O', 'N', N, N, V, LDV, S, U, LDU,
  1069. $ U, LDU, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1070. *
  1071. DO 1123 p = 1, N
  1072. V(p,p) = CONJG(V(p,p))
  1073. DO 1124 q = p + 1, N
  1074. CTMP = CONJG(V(q,p))
  1075. V(q,p) = CONJG(V(p,q))
  1076. V(p,q) = CTMP
  1077. 1124 CONTINUE
  1078. 1123 CONTINUE
  1079. CALL CLAPMT( .FALSE., N, N, V, LDV, IWORK )
  1080. END IF
  1081. *
  1082. ELSE
  1083. * .. aply CGESVD to R
  1084. * .. copy R into V and overwrite V with the right singular vectors
  1085. CALL CLACPY( 'U', NR, N, A, LDA, V, LDV )
  1086. IF ( NR .GT. 1 )
  1087. $ CALL CLASET( 'L', NR-1, NR-1, CZERO, CZERO, V(2,1), LDV )
  1088. * .. the right singular vectors overwrite V, the NR left singular
  1089. * vectors stored in U(1:NR,1:NR)
  1090. IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
  1091. CALL CGESVD( 'N', 'O', NR, N, V, LDV, S, U, LDU,
  1092. $ V, LDV, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1093. CALL CLAPMT( .FALSE., NR, N, V, LDV, IWORK )
  1094. * .. now [V](1:NR,1:N) contains V(1:N,1:NR)**H
  1095. ELSE
  1096. * .. need all N right singular vectors and NR < N
  1097. * [!] This is simple implementation that augments [V](1:NR,1:N)
  1098. * by padding a zero block. In the case NR << N, a more efficient
  1099. * way is to first use the LQ factorization. For more details
  1100. * how to implement this, see the " FULL SVD " branch.
  1101. CALL CLASET('G', N-NR, N, CZERO,CZERO, V(NR+1,1), LDV)
  1102. CALL CGESVD( 'N', 'O', N, N, V, LDV, S, U, LDU,
  1103. $ V, LDV, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1104. CALL CLAPMT( .FALSE., N, N, V, LDV, IWORK )
  1105. END IF
  1106. * .. now [V] contains the adjoint of the matrix of the right singular
  1107. * vectors of A.
  1108. END IF
  1109. *
  1110. ELSE
  1111. *.......................................................................
  1112. * .. FULL SVD requested
  1113. *.......................................................................
  1114. IF ( RTRANS ) THEN
  1115. *
  1116. * .. apply CGESVD to R**H [[this option is left for R&D&T]]
  1117. *
  1118. IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
  1119. * .. copy R**H into [V] and overwrite [V] with the left singular
  1120. * vectors of R**H
  1121. DO 1168 p = 1, NR
  1122. DO 1169 q = p, N
  1123. V(q,p) = CONJG(A(p,q))
  1124. 1169 CONTINUE
  1125. 1168 CONTINUE
  1126. IF ( NR .GT. 1 )
  1127. $ CALL CLASET( 'U', NR-1,NR-1, CZERO,CZERO, V(1,2), LDV )
  1128. *
  1129. * .. the left singular vectors of R**H overwrite [V], the NR right
  1130. * singular vectors of R**H stored in [U](1:NR,1:NR) as conjugate
  1131. * transposed
  1132. CALL CGESVD( 'O', 'A', N, NR, V, LDV, S, V, LDV,
  1133. $ U, LDU, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1134. * .. assemble V
  1135. DO 1115 p = 1, NR
  1136. V(p,p) = CONJG(V(p,p))
  1137. DO 1116 q = p + 1, NR
  1138. CTMP = CONJG(V(q,p))
  1139. V(q,p) = CONJG(V(p,q))
  1140. V(p,q) = CTMP
  1141. 1116 CONTINUE
  1142. 1115 CONTINUE
  1143. IF ( NR .LT. N ) THEN
  1144. DO 1101 p = 1, NR
  1145. DO 1102 q = NR+1, N
  1146. V(p,q) = CONJG(V(q,p))
  1147. 1102 CONTINUE
  1148. 1101 CONTINUE
  1149. END IF
  1150. CALL CLAPMT( .FALSE., NR, N, V, LDV, IWORK )
  1151. *
  1152. DO 1117 p = 1, NR
  1153. U(p,p) = CONJG(U(p,p))
  1154. DO 1118 q = p + 1, NR
  1155. CTMP = CONJG(U(q,p))
  1156. U(q,p) = CONJG(U(p,q))
  1157. U(p,q) = CTMP
  1158. 1118 CONTINUE
  1159. 1117 CONTINUE
  1160. *
  1161. IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
  1162. CALL CLASET('A', M-NR,NR, CZERO,CZERO, U(NR+1,1), LDU)
  1163. IF ( NR .LT. N1 ) THEN
  1164. CALL CLASET('A',NR,N1-NR,CZERO,CZERO,U(1,NR+1),LDU)
  1165. CALL CLASET( 'A',M-NR,N1-NR,CZERO,CONE,
  1166. $ U(NR+1,NR+1), LDU )
  1167. END IF
  1168. END IF
  1169. *
  1170. ELSE
  1171. * .. need all N right singular vectors and NR < N
  1172. * .. copy R**H into [V] and overwrite [V] with the left singular
  1173. * vectors of R**H
  1174. * [[The optimal ratio N/NR for using QRF instead of padding
  1175. * with zeros. Here hard coded to 2; it must be at least
  1176. * two due to work space constraints.]]
  1177. * OPTRATIO = ILAENV(6, 'CGESVD', 'S' // 'O', NR,N,0,0)
  1178. * OPTRATIO = MAX( OPTRATIO, 2 )
  1179. OPTRATIO = 2
  1180. IF ( OPTRATIO*NR .GT. N ) THEN
  1181. DO 1198 p = 1, NR
  1182. DO 1199 q = p, N
  1183. V(q,p) = CONJG(A(p,q))
  1184. 1199 CONTINUE
  1185. 1198 CONTINUE
  1186. IF ( NR .GT. 1 )
  1187. $ CALL CLASET('U',NR-1,NR-1, CZERO,CZERO, V(1,2),LDV)
  1188. *
  1189. CALL CLASET('A',N,N-NR,CZERO,CZERO,V(1,NR+1),LDV)
  1190. CALL CGESVD( 'O', 'A', N, N, V, LDV, S, V, LDV,
  1191. $ U, LDU, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1192. *
  1193. DO 1113 p = 1, N
  1194. V(p,p) = CONJG(V(p,p))
  1195. DO 1114 q = p + 1, N
  1196. CTMP = CONJG(V(q,p))
  1197. V(q,p) = CONJG(V(p,q))
  1198. V(p,q) = CTMP
  1199. 1114 CONTINUE
  1200. 1113 CONTINUE
  1201. CALL CLAPMT( .FALSE., N, N, V, LDV, IWORK )
  1202. * .. assemble the left singular vector matrix U of dimensions
  1203. * (M x N1), i.e. (M x N) or (M x M).
  1204. *
  1205. DO 1111 p = 1, N
  1206. U(p,p) = CONJG(U(p,p))
  1207. DO 1112 q = p + 1, N
  1208. CTMP = CONJG(U(q,p))
  1209. U(q,p) = CONJG(U(p,q))
  1210. U(p,q) = CTMP
  1211. 1112 CONTINUE
  1212. 1111 CONTINUE
  1213. *
  1214. IF ( ( N .LT. M ) .AND. .NOT.(WNTUF)) THEN
  1215. CALL CLASET('A',M-N,N,CZERO,CZERO,U(N+1,1),LDU)
  1216. IF ( N .LT. N1 ) THEN
  1217. CALL CLASET('A',N,N1-N,CZERO,CZERO,U(1,N+1),LDU)
  1218. CALL CLASET('A',M-N,N1-N,CZERO,CONE,
  1219. $ U(N+1,N+1), LDU )
  1220. END IF
  1221. END IF
  1222. ELSE
  1223. * .. copy R**H into [U] and overwrite [U] with the right
  1224. * singular vectors of R
  1225. DO 1196 p = 1, NR
  1226. DO 1197 q = p, N
  1227. U(q,NR+p) = CONJG(A(p,q))
  1228. 1197 CONTINUE
  1229. 1196 CONTINUE
  1230. IF ( NR .GT. 1 )
  1231. $ CALL CLASET('U',NR-1,NR-1,CZERO,CZERO,U(1,NR+2),LDU)
  1232. CALL CGEQRF( N, NR, U(1,NR+1), LDU, CWORK(N+1),
  1233. $ CWORK(N+NR+1), LCWORK-N-NR, IERR )
  1234. DO 1143 p = 1, NR
  1235. DO 1144 q = 1, N
  1236. V(q,p) = CONJG(U(p,NR+q))
  1237. 1144 CONTINUE
  1238. 1143 CONTINUE
  1239. CALL CLASET('U',NR-1,NR-1,CZERO,CZERO,V(1,2),LDV)
  1240. CALL CGESVD( 'S', 'O', NR, NR, V, LDV, S, U, LDU,
  1241. $ V,LDV, CWORK(N+NR+1),LCWORK-N-NR,RWORK, INFO )
  1242. CALL CLASET('A',N-NR,NR,CZERO,CZERO,V(NR+1,1),LDV)
  1243. CALL CLASET('A',NR,N-NR,CZERO,CZERO,V(1,NR+1),LDV)
  1244. CALL CLASET('A',N-NR,N-NR,CZERO,CONE,V(NR+1,NR+1),LDV)
  1245. CALL CUNMQR('R','C', N, N, NR, U(1,NR+1), LDU,
  1246. $ CWORK(N+1),V,LDV,CWORK(N+NR+1),LCWORK-N-NR,IERR)
  1247. CALL CLAPMT( .FALSE., N, N, V, LDV, IWORK )
  1248. * .. assemble the left singular vector matrix U of dimensions
  1249. * (M x NR) or (M x N) or (M x M).
  1250. IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
  1251. CALL CLASET('A',M-NR,NR,CZERO,CZERO,U(NR+1,1),LDU)
  1252. IF ( NR .LT. N1 ) THEN
  1253. CALL CLASET('A',NR,N1-NR,CZERO,CZERO,U(1,NR+1),LDU)
  1254. CALL CLASET( 'A',M-NR,N1-NR,CZERO,CONE,
  1255. $ U(NR+1,NR+1),LDU)
  1256. END IF
  1257. END IF
  1258. END IF
  1259. END IF
  1260. *
  1261. ELSE
  1262. *
  1263. * .. apply CGESVD to R [[this is the recommended option]]
  1264. *
  1265. IF ( WNTVR .OR. ( NR .EQ. N ) ) THEN
  1266. * .. copy R into [V] and overwrite V with the right singular vectors
  1267. CALL CLACPY( 'U', NR, N, A, LDA, V, LDV )
  1268. IF ( NR .GT. 1 )
  1269. $ CALL CLASET( 'L', NR-1,NR-1, CZERO,CZERO, V(2,1), LDV )
  1270. * .. the right singular vectors of R overwrite [V], the NR left
  1271. * singular vectors of R stored in [U](1:NR,1:NR)
  1272. CALL CGESVD( 'S', 'O', NR, N, V, LDV, S, U, LDU,
  1273. $ V, LDV, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1274. CALL CLAPMT( .FALSE., NR, N, V, LDV, IWORK )
  1275. * .. now [V](1:NR,1:N) contains V(1:N,1:NR)**H
  1276. * .. assemble the left singular vector matrix U of dimensions
  1277. * (M x NR) or (M x N) or (M x M).
  1278. IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
  1279. CALL CLASET('A', M-NR,NR, CZERO,CZERO, U(NR+1,1), LDU)
  1280. IF ( NR .LT. N1 ) THEN
  1281. CALL CLASET('A',NR,N1-NR,CZERO,CZERO,U(1,NR+1),LDU)
  1282. CALL CLASET( 'A',M-NR,N1-NR,CZERO,CONE,
  1283. $ U(NR+1,NR+1), LDU )
  1284. END IF
  1285. END IF
  1286. *
  1287. ELSE
  1288. * .. need all N right singular vectors and NR < N
  1289. * .. the requested number of the left singular vectors
  1290. * is then N1 (N or M)
  1291. * [[The optimal ratio N/NR for using LQ instead of padding
  1292. * with zeros. Here hard coded to 2; it must be at least
  1293. * two due to work space constraints.]]
  1294. * OPTRATIO = ILAENV(6, 'CGESVD', 'S' // 'O', NR,N,0,0)
  1295. * OPTRATIO = MAX( OPTRATIO, 2 )
  1296. OPTRATIO = 2
  1297. IF ( OPTRATIO * NR .GT. N ) THEN
  1298. CALL CLACPY( 'U', NR, N, A, LDA, V, LDV )
  1299. IF ( NR .GT. 1 )
  1300. $ CALL CLASET('L', NR-1,NR-1, CZERO,CZERO, V(2,1),LDV)
  1301. * .. the right singular vectors of R overwrite [V], the NR left
  1302. * singular vectors of R stored in [U](1:NR,1:NR)
  1303. CALL CLASET('A', N-NR,N, CZERO,CZERO, V(NR+1,1),LDV)
  1304. CALL CGESVD( 'S', 'O', N, N, V, LDV, S, U, LDU,
  1305. $ V, LDV, CWORK(N+1), LCWORK-N, RWORK, INFO )
  1306. CALL CLAPMT( .FALSE., N, N, V, LDV, IWORK )
  1307. * .. now [V] contains the adjoint of the matrix of the right
  1308. * singular vectors of A. The leading N left singular vectors
  1309. * are in [U](1:N,1:N)
  1310. * .. assemble the left singular vector matrix U of dimensions
  1311. * (M x N1), i.e. (M x N) or (M x M).
  1312. IF ( ( N .LT. M ) .AND. .NOT.(WNTUF)) THEN
  1313. CALL CLASET('A',M-N,N,CZERO,CZERO,U(N+1,1),LDU)
  1314. IF ( N .LT. N1 ) THEN
  1315. CALL CLASET('A',N,N1-N,CZERO,CZERO,U(1,N+1),LDU)
  1316. CALL CLASET( 'A',M-N,N1-N,CZERO,CONE,
  1317. $ U(N+1,N+1), LDU )
  1318. END IF
  1319. END IF
  1320. ELSE
  1321. CALL CLACPY( 'U', NR, N, A, LDA, U(NR+1,1), LDU )
  1322. IF ( NR .GT. 1 )
  1323. $ CALL CLASET('L',NR-1,NR-1,CZERO,CZERO,U(NR+2,1),LDU)
  1324. CALL CGELQF( NR, N, U(NR+1,1), LDU, CWORK(N+1),
  1325. $ CWORK(N+NR+1), LCWORK-N-NR, IERR )
  1326. CALL CLACPY('L',NR,NR,U(NR+1,1),LDU,V,LDV)
  1327. IF ( NR .GT. 1 )
  1328. $ CALL CLASET('U',NR-1,NR-1,CZERO,CZERO,V(1,2),LDV)
  1329. CALL CGESVD( 'S', 'O', NR, NR, V, LDV, S, U, LDU,
  1330. $ V, LDV, CWORK(N+NR+1), LCWORK-N-NR, RWORK, INFO )
  1331. CALL CLASET('A',N-NR,NR,CZERO,CZERO,V(NR+1,1),LDV)
  1332. CALL CLASET('A',NR,N-NR,CZERO,CZERO,V(1,NR+1),LDV)
  1333. CALL CLASET('A',N-NR,N-NR,CZERO,CONE,V(NR+1,NR+1),LDV)
  1334. CALL CUNMLQ('R','N',N,N,NR,U(NR+1,1),LDU,CWORK(N+1),
  1335. $ V, LDV, CWORK(N+NR+1),LCWORK-N-NR,IERR)
  1336. CALL CLAPMT( .FALSE., N, N, V, LDV, IWORK )
  1337. * .. assemble the left singular vector matrix U of dimensions
  1338. * (M x NR) or (M x N) or (M x M).
  1339. IF ( ( NR .LT. M ) .AND. .NOT.(WNTUF)) THEN
  1340. CALL CLASET('A',M-NR,NR,CZERO,CZERO,U(NR+1,1),LDU)
  1341. IF ( NR .LT. N1 ) THEN
  1342. CALL CLASET('A',NR,N1-NR,CZERO,CZERO,U(1,NR+1),LDU)
  1343. CALL CLASET( 'A',M-NR,N1-NR,CZERO,CONE,
  1344. $ U(NR+1,NR+1), LDU )
  1345. END IF
  1346. END IF
  1347. END IF
  1348. END IF
  1349. * .. end of the "R**H or R" branch
  1350. END IF
  1351. *
  1352. * The Q matrix from the first QRF is built into the left singular
  1353. * vectors matrix U.
  1354. *
  1355. IF ( .NOT. WNTUF )
  1356. $ CALL CUNMQR( 'L', 'N', M, N1, N, A, LDA, CWORK, U,
  1357. $ LDU, CWORK(N+1), LCWORK-N, IERR )
  1358. IF ( ROWPRM .AND. .NOT.WNTUF )
  1359. $ CALL CLASWP( N1, U, LDU, 1, M-1, IWORK(N+1), -1 )
  1360. *
  1361. * ... end of the "full SVD" branch
  1362. END IF
  1363. *
  1364. * Check whether some singular values are returned as zeros, e.g.
  1365. * due to underflow, and update the numerical rank.
  1366. p = NR
  1367. DO 4001 q = p, 1, -1
  1368. IF ( S(q) .GT. ZERO ) GO TO 4002
  1369. NR = NR - 1
  1370. 4001 CONTINUE
  1371. 4002 CONTINUE
  1372. *
  1373. * .. if numerical rank deficiency is detected, the truncated
  1374. * singular values are set to zero.
  1375. IF ( NR .LT. N ) CALL SLASET( 'G', N-NR,1, ZERO,ZERO, S(NR+1), N )
  1376. * .. undo scaling; this may cause overflow in the largest singular
  1377. * values.
  1378. IF ( ASCALED )
  1379. $ CALL SLASCL( 'G',0,0, ONE,SQRT(REAL(M)), NR,1, S, N, IERR )
  1380. IF ( CONDA ) RWORK(1) = SCONDA
  1381. RWORK(2) = p - NR
  1382. * .. p-NR is the number of singular values that are computed as
  1383. * exact zeros in CGESVD() applied to the (possibly truncated)
  1384. * full row rank triangular (trapezoidal) factor of A.
  1385. NUMRANK = NR
  1386. *
  1387. RETURN
  1388. *
  1389. * End of CGESVDQ
  1390. *
  1391. END