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.

sgesvdq.f 58 kB

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