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.

sgesvj.f 68 kB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598
  1. *> \brief \b SGESVJ
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download SGESVJ + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgesvj.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgesvj.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgesvj.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE SGESVJ( JOBA, JOBU, JOBV, M, N, A, LDA, SVA, MV, V,
  22. * LDV, WORK, LWORK, INFO )
  23. *
  24. * .. Scalar Arguments ..
  25. * INTEGER INFO, LDA, LDV, LWORK, M, MV, N
  26. * CHARACTER*1 JOBA, JOBU, JOBV
  27. * ..
  28. * .. Array Arguments ..
  29. * REAL A( LDA, * ), SVA( N ), V( LDV, * ),
  30. * $ WORK( LWORK )
  31. * ..
  32. *
  33. *
  34. *> \par Purpose:
  35. * =============
  36. *>
  37. *> \verbatim
  38. *>
  39. *> SGESVJ computes the singular value decomposition (SVD) of a real
  40. *> M-by-N matrix A, where M >= N. The SVD of A is written as
  41. *> [++] [xx] [x0] [xx]
  42. *> A = U * SIGMA * V^t, [++] = [xx] * [ox] * [xx]
  43. *> [++] [xx]
  44. *> where SIGMA is an N-by-N diagonal matrix, U is an M-by-N orthonormal
  45. *> matrix, and V is an N-by-N orthogonal matrix. The diagonal elements
  46. *> of SIGMA are the singular values of A. The columns of U and V are the
  47. *> left and the right singular vectors of A, respectively.
  48. *> SGESVJ can sometimes compute tiny singular values and their singular vectors much
  49. *> more accurately than other SVD routines, see below under Further Details.
  50. *> \endverbatim
  51. *
  52. * Arguments:
  53. * ==========
  54. *
  55. *> \param[in] JOBA
  56. *> \verbatim
  57. *> JOBA is CHARACTER*1
  58. *> Specifies the structure of A.
  59. *> = 'L': The input matrix A is lower triangular;
  60. *> = 'U': The input matrix A is upper triangular;
  61. *> = 'G': The input matrix A is general M-by-N matrix, M >= N.
  62. *> \endverbatim
  63. *>
  64. *> \param[in] JOBU
  65. *> \verbatim
  66. *> JOBU is CHARACTER*1
  67. *> Specifies whether to compute the left singular vectors
  68. *> (columns of U):
  69. *> = 'U': The left singular vectors corresponding to the nonzero
  70. *> singular values are computed and returned in the leading
  71. *> columns of A. See more details in the description of A.
  72. *> The default numerical orthogonality threshold is set to
  73. *> approximately TOL=CTOL*EPS, CTOL=SQRT(M), EPS=SLAMCH('E').
  74. *> = 'C': Analogous to JOBU='U', except that user can control the
  75. *> level of numerical orthogonality of the computed left
  76. *> singular vectors. TOL can be set to TOL = CTOL*EPS, where
  77. *> CTOL is given on input in the array WORK.
  78. *> No CTOL smaller than ONE is allowed. CTOL greater
  79. *> than 1 / EPS is meaningless. The option 'C'
  80. *> can be used if M*EPS is satisfactory orthogonality
  81. *> of the computed left singular vectors, so CTOL=M could
  82. *> save few sweeps of Jacobi rotations.
  83. *> See the descriptions of A and WORK(1).
  84. *> = 'N': The matrix U is not computed. However, see the
  85. *> description of A.
  86. *> \endverbatim
  87. *>
  88. *> \param[in] JOBV
  89. *> \verbatim
  90. *> JOBV is CHARACTER*1
  91. *> Specifies whether to compute the right singular vectors, that
  92. *> is, the matrix V:
  93. *> = 'V': the matrix V is computed and returned in the array V
  94. *> = 'A': the Jacobi rotations are applied to the MV-by-N
  95. *> array V. In other words, the right singular vector
  96. *> matrix V is not computed explicitly; instead it is
  97. *> applied to an MV-by-N matrix initially stored in the
  98. *> first MV rows of V.
  99. *> = 'N': the matrix V is not computed and the array V is not
  100. *> referenced
  101. *> \endverbatim
  102. *>
  103. *> \param[in] M
  104. *> \verbatim
  105. *> M is INTEGER
  106. *> The number of rows of the input matrix A. 1/SLAMCH('E') > M >= 0.
  107. *> \endverbatim
  108. *>
  109. *> \param[in] N
  110. *> \verbatim
  111. *> N is INTEGER
  112. *> The number of columns of the input matrix A.
  113. *> M >= N >= 0.
  114. *> \endverbatim
  115. *>
  116. *> \param[in,out] A
  117. *> \verbatim
  118. *> A is REAL array, dimension (LDA,N)
  119. *> On entry, the M-by-N matrix A.
  120. *> On exit,
  121. *> If JOBU = 'U' .OR. JOBU = 'C':
  122. *> If INFO = 0:
  123. *> RANKA orthonormal columns of U are returned in the
  124. *> leading RANKA columns of the array A. Here RANKA <= N
  125. *> is the number of computed singular values of A that are
  126. *> above the underflow threshold SLAMCH('S'). The singular
  127. *> vectors corresponding to underflowed or zero singular
  128. *> values are not computed. The value of RANKA is returned
  129. *> in the array WORK as RANKA=NINT(WORK(2)). Also see the
  130. *> descriptions of SVA and WORK. The computed columns of U
  131. *> are mutually numerically orthogonal up to approximately
  132. *> TOL=SQRT(M)*EPS (default); or TOL=CTOL*EPS (JOBU = 'C'),
  133. *> see the description of JOBU.
  134. *> If INFO > 0,
  135. *> the procedure SGESVJ did not converge in the given number
  136. *> of iterations (sweeps). In that case, the computed
  137. *> columns of U may not be orthogonal up to TOL. The output
  138. *> U (stored in A), SIGMA (given by the computed singular
  139. *> values in SVA(1:N)) and V is still a decomposition of the
  140. *> input matrix A in the sense that the residual
  141. *> ||A-SCALE*U*SIGMA*V^T||_2 / ||A||_2 is small.
  142. *> If JOBU = 'N':
  143. *> If INFO = 0:
  144. *> Note that the left singular vectors are 'for free' in the
  145. *> one-sided Jacobi SVD algorithm. However, if only the
  146. *> singular values are needed, the level of numerical
  147. *> orthogonality of U is not an issue and iterations are
  148. *> stopped when the columns of the iterated matrix are
  149. *> numerically orthogonal up to approximately M*EPS. Thus,
  150. *> on exit, A contains the columns of U scaled with the
  151. *> corresponding singular values.
  152. *> If INFO > 0:
  153. *> the procedure SGESVJ did not converge in the given number
  154. *> of iterations (sweeps).
  155. *> \endverbatim
  156. *>
  157. *> \param[in] LDA
  158. *> \verbatim
  159. *> LDA is INTEGER
  160. *> The leading dimension of the array A. LDA >= max(1,M).
  161. *> \endverbatim
  162. *>
  163. *> \param[out] SVA
  164. *> \verbatim
  165. *> SVA is REAL array, dimension (N)
  166. *> On exit,
  167. *> If INFO = 0 :
  168. *> depending on the value SCALE = WORK(1), we have:
  169. *> If SCALE = ONE:
  170. *> SVA(1:N) contains the computed singular values of A.
  171. *> During the computation SVA contains the Euclidean column
  172. *> norms of the iterated matrices in the array A.
  173. *> If SCALE .NE. ONE:
  174. *> The singular values of A are SCALE*SVA(1:N), and this
  175. *> factored representation is due to the fact that some of the
  176. *> singular values of A might underflow or overflow.
  177. *>
  178. *> If INFO > 0 :
  179. *> the procedure SGESVJ did not converge in the given number of
  180. *> iterations (sweeps) and SCALE*SVA(1:N) may not be accurate.
  181. *> \endverbatim
  182. *>
  183. *> \param[in] MV
  184. *> \verbatim
  185. *> MV is INTEGER
  186. *> If JOBV = 'A', then the product of Jacobi rotations in SGESVJ
  187. *> is applied to the first MV rows of V. See the description of JOBV.
  188. *> \endverbatim
  189. *>
  190. *> \param[in,out] V
  191. *> \verbatim
  192. *> V is REAL array, dimension (LDV,N)
  193. *> If JOBV = 'V', then V contains on exit the N-by-N matrix of
  194. *> the right singular vectors;
  195. *> If JOBV = 'A', then V contains the product of the computed right
  196. *> singular vector matrix and the initial matrix in
  197. *> the array V.
  198. *> If JOBV = 'N', then V is not referenced.
  199. *> \endverbatim
  200. *>
  201. *> \param[in] LDV
  202. *> \verbatim
  203. *> LDV is INTEGER
  204. *> The leading dimension of the array V, LDV >= 1.
  205. *> If JOBV = 'V', then LDV >= max(1,N).
  206. *> If JOBV = 'A', then LDV >= max(1,MV) .
  207. *> \endverbatim
  208. *>
  209. *> \param[in,out] WORK
  210. *> \verbatim
  211. *> WORK is REAL array, dimension (LWORK)
  212. *> On entry,
  213. *> If JOBU = 'C' :
  214. *> WORK(1) = CTOL, where CTOL defines the threshold for convergence.
  215. *> The process stops if all columns of A are mutually
  216. *> orthogonal up to CTOL*EPS, EPS=SLAMCH('E').
  217. *> It is required that CTOL >= ONE, i.e. it is not
  218. *> allowed to force the routine to obtain orthogonality
  219. *> below EPSILON.
  220. *> On exit,
  221. *> WORK(1) = SCALE is the scaling factor such that SCALE*SVA(1:N)
  222. *> are the computed singular vcalues of A.
  223. *> (See description of SVA().)
  224. *> WORK(2) = NINT(WORK(2)) is the number of the computed nonzero
  225. *> singular values.
  226. *> WORK(3) = NINT(WORK(3)) is the number of the computed singular
  227. *> values that are larger than the underflow threshold.
  228. *> WORK(4) = NINT(WORK(4)) is the number of sweeps of Jacobi
  229. *> rotations needed for numerical convergence.
  230. *> WORK(5) = max_{i.NE.j} |COS(A(:,i),A(:,j))| in the last sweep.
  231. *> This is useful information in cases when SGESVJ did
  232. *> not converge, as it can be used to estimate whether
  233. *> the output is still useful and for post festum analysis.
  234. *> WORK(6) = the largest absolute value over all sines of the
  235. *> Jacobi rotation angles in the last sweep. It can be
  236. *> useful for a post festum analysis.
  237. *> \endverbatim
  238. *>
  239. *> \param[in] LWORK
  240. *> \verbatim
  241. *> LWORK is INTEGER
  242. *> length of WORK, WORK >= MAX(6,M+N)
  243. *> \endverbatim
  244. *>
  245. *> \param[out] INFO
  246. *> \verbatim
  247. *> INFO is INTEGER
  248. *> = 0: successful exit.
  249. *> < 0: if INFO = -i, then the i-th argument had an illegal value
  250. *> > 0: SGESVJ did not converge in the maximal allowed number (30)
  251. *> of sweeps. The output may still be useful. See the
  252. *> description of WORK.
  253. *> \endverbatim
  254. *
  255. * Authors:
  256. * ========
  257. *
  258. *> \author Univ. of Tennessee
  259. *> \author Univ. of California Berkeley
  260. *> \author Univ. of Colorado Denver
  261. *> \author NAG Ltd.
  262. *
  263. *> \ingroup realGEcomputational
  264. *
  265. *> \par Further Details:
  266. * =====================
  267. *>
  268. *> The orthogonal N-by-N matrix V is obtained as a product of Jacobi plane
  269. *> rotations. The rotations are implemented as fast scaled rotations of
  270. *> Anda and Park [1]. In the case of underflow of the Jacobi angle, a
  271. *> modified Jacobi transformation of Drmac [4] is used. Pivot strategy uses
  272. *> column interchanges of de Rijk [2]. The relative accuracy of the computed
  273. *> singular values and the accuracy of the computed singular vectors (in
  274. *> angle metric) is as guaranteed by the theory of Demmel and Veselic [3].
  275. *> The condition number that determines the accuracy in the full rank case
  276. *> is essentially min_{D=diag} kappa(A*D), where kappa(.) is the
  277. *> spectral condition number. The best performance of this Jacobi SVD
  278. *> procedure is achieved if used in an accelerated version of Drmac and
  279. *> Veselic [5,6], and it is the kernel routine in the SIGMA library [7].
  280. *> Some tuning parameters (marked with [TP]) are available for the
  281. *> implementer. \n
  282. *> The computational range for the nonzero singular values is the machine
  283. *> number interval ( UNDERFLOW , OVERFLOW ). In extreme cases, even
  284. *> denormalized singular values can be computed with the corresponding
  285. *> gradual loss of accurate digits.
  286. *>
  287. *> \par Contributors:
  288. * ==================
  289. *>
  290. *> Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
  291. *>
  292. *> \par References:
  293. * ================
  294. *>
  295. *> [1] A. A. Anda and H. Park: Fast plane rotations with dynamic scaling. \n
  296. *> SIAM J. matrix Anal. Appl., Vol. 15 (1994), pp. 162-174. \n\n
  297. *> [2] P. P. M. De Rijk: A one-sided Jacobi algorithm for computing the
  298. *> singular value decomposition on a vector computer. \n
  299. *> SIAM J. Sci. Stat. Comp., Vol. 10 (1998), pp. 359-371. \n\n
  300. *> [3] J. Demmel and K. Veselic: Jacobi method is more accurate than QR. \n
  301. *> [4] Z. Drmac: Implementation of Jacobi rotations for accurate singular
  302. *> value computation in floating point arithmetic. \n
  303. *> SIAM J. Sci. Comp., Vol. 18 (1997), pp. 1200-1222. \n\n
  304. *> [5] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm I. \n
  305. *> SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1322-1342. \n
  306. *> LAPACK Working note 169. \n\n
  307. *> [6] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm II. \n
  308. *> SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1343-1362. \n
  309. *> LAPACK Working note 170. \n\n
  310. *> [7] Z. Drmac: SIGMA - mathematical software library for accurate SVD, PSV,
  311. *> QSVD, (H,K)-SVD computations.\n
  312. *> Department of Mathematics, University of Zagreb, 2008.
  313. *>
  314. *> \par Bugs, Examples and Comments:
  315. * =================================
  316. *>
  317. *> Please report all bugs and send interesting test examples and comments to
  318. *> drmac@math.hr. Thank you.
  319. *
  320. * =====================================================================
  321. SUBROUTINE SGESVJ( JOBA, JOBU, JOBV, M, N, A, LDA, SVA, MV, V,
  322. $ LDV, WORK, LWORK, INFO )
  323. *
  324. * -- LAPACK computational routine --
  325. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  326. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  327. *
  328. * .. Scalar Arguments ..
  329. INTEGER INFO, LDA, LDV, LWORK, M, MV, N
  330. CHARACTER*1 JOBA, JOBU, JOBV
  331. * ..
  332. * .. Array Arguments ..
  333. REAL A( LDA, * ), SVA( N ), V( LDV, * ),
  334. $ WORK( LWORK )
  335. * ..
  336. *
  337. * =====================================================================
  338. *
  339. * .. Local Parameters ..
  340. REAL ZERO, HALF, ONE
  341. PARAMETER ( ZERO = 0.0E0, HALF = 0.5E0, ONE = 1.0E0)
  342. INTEGER NSWEEP
  343. PARAMETER ( NSWEEP = 30 )
  344. * ..
  345. * .. Local Scalars ..
  346. REAL AAPP, AAPP0, AAPQ, AAQQ, APOAQ, AQOAP, BIG,
  347. $ BIGTHETA, CS, CTOL, EPSLN, LARGE, MXAAPQ,
  348. $ MXSINJ, ROOTBIG, ROOTEPS, ROOTSFMIN, ROOTTOL,
  349. $ SKL, SFMIN, SMALL, SN, T, TEMP1, THETA,
  350. $ THSIGN, TOL
  351. INTEGER BLSKIP, EMPTSW, i, ibr, IERR, igl, IJBLSK, ir1,
  352. $ ISWROT, jbc, jgl, KBL, LKAHEAD, MVL, N2, N34,
  353. $ N4, NBL, NOTROT, p, PSKIPPED, q, ROWSKIP,
  354. $ SWBAND
  355. LOGICAL APPLV, GOSCALE, LOWER, LSVEC, NOSCALE, ROTOK,
  356. $ RSVEC, UCTOL, UPPER
  357. * ..
  358. * .. Local Arrays ..
  359. REAL FASTR( 5 )
  360. * ..
  361. * .. Intrinsic Functions ..
  362. INTRINSIC ABS, MAX, MIN, FLOAT, SIGN, SQRT
  363. * ..
  364. * .. External Functions ..
  365. * ..
  366. * from BLAS
  367. REAL SDOT, SNRM2
  368. EXTERNAL SDOT, SNRM2
  369. INTEGER ISAMAX
  370. EXTERNAL ISAMAX
  371. * from LAPACK
  372. REAL SLAMCH
  373. EXTERNAL SLAMCH
  374. LOGICAL LSAME
  375. EXTERNAL LSAME
  376. * ..
  377. * .. External Subroutines ..
  378. * ..
  379. * from BLAS
  380. EXTERNAL SAXPY, SCOPY, SROTM, SSCAL, SSWAP
  381. * from LAPACK
  382. EXTERNAL SLASCL, SLASET, SLASSQ, XERBLA
  383. *
  384. EXTERNAL SGSVJ0, SGSVJ1
  385. * ..
  386. * .. Executable Statements ..
  387. *
  388. * Test the input arguments
  389. *
  390. LSVEC = LSAME( JOBU, 'U' )
  391. UCTOL = LSAME( JOBU, 'C' )
  392. RSVEC = LSAME( JOBV, 'V' )
  393. APPLV = LSAME( JOBV, 'A' )
  394. UPPER = LSAME( JOBA, 'U' )
  395. LOWER = LSAME( JOBA, 'L' )
  396. *
  397. IF( .NOT.( UPPER .OR. LOWER .OR. LSAME( JOBA, 'G' ) ) ) THEN
  398. INFO = -1
  399. ELSE IF( .NOT.( LSVEC .OR. UCTOL .OR. LSAME( JOBU, 'N' ) ) ) THEN
  400. INFO = -2
  401. ELSE IF( .NOT.( RSVEC .OR. APPLV .OR. LSAME( JOBV, 'N' ) ) ) THEN
  402. INFO = -3
  403. ELSE IF( M.LT.0 ) THEN
  404. INFO = -4
  405. ELSE IF( ( N.LT.0 ) .OR. ( N.GT.M ) ) THEN
  406. INFO = -5
  407. ELSE IF( LDA.LT.M ) THEN
  408. INFO = -7
  409. ELSE IF( MV.LT.0 ) THEN
  410. INFO = -9
  411. ELSE IF( ( RSVEC .AND. ( LDV.LT.N ) ) .OR.
  412. $ ( APPLV .AND. ( LDV.LT.MV ) ) ) THEN
  413. INFO = -11
  414. ELSE IF( UCTOL .AND. ( WORK( 1 ).LE.ONE ) ) THEN
  415. INFO = -12
  416. ELSE IF( LWORK.LT.MAX( M+N, 6 ) ) THEN
  417. INFO = -13
  418. ELSE
  419. INFO = 0
  420. END IF
  421. *
  422. * #:(
  423. IF( INFO.NE.0 ) THEN
  424. CALL XERBLA( 'SGESVJ', -INFO )
  425. RETURN
  426. END IF
  427. *
  428. * #:) Quick return for void matrix
  429. *
  430. IF( ( M.EQ.0 ) .OR. ( N.EQ.0 ) )RETURN
  431. *
  432. * Set numerical parameters
  433. * The stopping criterion for Jacobi rotations is
  434. *
  435. * max_{i<>j}|A(:,i)^T * A(:,j)|/(||A(:,i)||*||A(:,j)||) < CTOL*EPS
  436. *
  437. * where EPS is the round-off and CTOL is defined as follows:
  438. *
  439. IF( UCTOL ) THEN
  440. * ... user controlled
  441. CTOL = WORK( 1 )
  442. ELSE
  443. * ... default
  444. IF( LSVEC .OR. RSVEC .OR. APPLV ) THEN
  445. CTOL = SQRT( FLOAT( M ) )
  446. ELSE
  447. CTOL = FLOAT( M )
  448. END IF
  449. END IF
  450. * ... and the machine dependent parameters are
  451. *[!] (Make sure that SLAMCH() works properly on the target machine.)
  452. *
  453. EPSLN = SLAMCH( 'Epsilon' )
  454. ROOTEPS = SQRT( EPSLN )
  455. SFMIN = SLAMCH( 'SafeMinimum' )
  456. ROOTSFMIN = SQRT( SFMIN )
  457. SMALL = SFMIN / EPSLN
  458. BIG = SLAMCH( 'Overflow' )
  459. * BIG = ONE / SFMIN
  460. ROOTBIG = ONE / ROOTSFMIN
  461. LARGE = BIG / SQRT( FLOAT( M*N ) )
  462. BIGTHETA = ONE / ROOTEPS
  463. *
  464. TOL = CTOL*EPSLN
  465. ROOTTOL = SQRT( TOL )
  466. *
  467. IF( FLOAT( M )*EPSLN.GE.ONE ) THEN
  468. INFO = -4
  469. CALL XERBLA( 'SGESVJ', -INFO )
  470. RETURN
  471. END IF
  472. *
  473. * Initialize the right singular vector matrix.
  474. *
  475. IF( RSVEC ) THEN
  476. MVL = N
  477. CALL SLASET( 'A', MVL, N, ZERO, ONE, V, LDV )
  478. ELSE IF( APPLV ) THEN
  479. MVL = MV
  480. END IF
  481. RSVEC = RSVEC .OR. APPLV
  482. *
  483. * Initialize SVA( 1:N ) = ( ||A e_i||_2, i = 1:N )
  484. *(!) If necessary, scale A to protect the largest singular value
  485. * from overflow. It is possible that saving the largest singular
  486. * value destroys the information about the small ones.
  487. * This initial scaling is almost minimal in the sense that the
  488. * goal is to make sure that no column norm overflows, and that
  489. * SQRT(N)*max_i SVA(i) does not overflow. If INFinite entries
  490. * in A are detected, the procedure returns with INFO=-6.
  491. *
  492. SKL = ONE / SQRT( FLOAT( M )*FLOAT( N ) )
  493. NOSCALE = .TRUE.
  494. GOSCALE = .TRUE.
  495. *
  496. IF( LOWER ) THEN
  497. * the input matrix is M-by-N lower triangular (trapezoidal)
  498. DO 1874 p = 1, N
  499. AAPP = ZERO
  500. AAQQ = ONE
  501. CALL SLASSQ( M-p+1, A( p, p ), 1, AAPP, AAQQ )
  502. IF( AAPP.GT.BIG ) THEN
  503. INFO = -6
  504. CALL XERBLA( 'SGESVJ', -INFO )
  505. RETURN
  506. END IF
  507. AAQQ = SQRT( AAQQ )
  508. IF( ( AAPP.LT.( BIG / AAQQ ) ) .AND. NOSCALE ) THEN
  509. SVA( p ) = AAPP*AAQQ
  510. ELSE
  511. NOSCALE = .FALSE.
  512. SVA( p ) = AAPP*( AAQQ*SKL )
  513. IF( GOSCALE ) THEN
  514. GOSCALE = .FALSE.
  515. DO 1873 q = 1, p - 1
  516. SVA( q ) = SVA( q )*SKL
  517. 1873 CONTINUE
  518. END IF
  519. END IF
  520. 1874 CONTINUE
  521. ELSE IF( UPPER ) THEN
  522. * the input matrix is M-by-N upper triangular (trapezoidal)
  523. DO 2874 p = 1, N
  524. AAPP = ZERO
  525. AAQQ = ONE
  526. CALL SLASSQ( p, A( 1, p ), 1, AAPP, AAQQ )
  527. IF( AAPP.GT.BIG ) THEN
  528. INFO = -6
  529. CALL XERBLA( 'SGESVJ', -INFO )
  530. RETURN
  531. END IF
  532. AAQQ = SQRT( AAQQ )
  533. IF( ( AAPP.LT.( BIG / AAQQ ) ) .AND. NOSCALE ) THEN
  534. SVA( p ) = AAPP*AAQQ
  535. ELSE
  536. NOSCALE = .FALSE.
  537. SVA( p ) = AAPP*( AAQQ*SKL )
  538. IF( GOSCALE ) THEN
  539. GOSCALE = .FALSE.
  540. DO 2873 q = 1, p - 1
  541. SVA( q ) = SVA( q )*SKL
  542. 2873 CONTINUE
  543. END IF
  544. END IF
  545. 2874 CONTINUE
  546. ELSE
  547. * the input matrix is M-by-N general dense
  548. DO 3874 p = 1, N
  549. AAPP = ZERO
  550. AAQQ = ONE
  551. CALL SLASSQ( M, A( 1, p ), 1, AAPP, AAQQ )
  552. IF( AAPP.GT.BIG ) THEN
  553. INFO = -6
  554. CALL XERBLA( 'SGESVJ', -INFO )
  555. RETURN
  556. END IF
  557. AAQQ = SQRT( AAQQ )
  558. IF( ( AAPP.LT.( BIG / AAQQ ) ) .AND. NOSCALE ) THEN
  559. SVA( p ) = AAPP*AAQQ
  560. ELSE
  561. NOSCALE = .FALSE.
  562. SVA( p ) = AAPP*( AAQQ*SKL )
  563. IF( GOSCALE ) THEN
  564. GOSCALE = .FALSE.
  565. DO 3873 q = 1, p - 1
  566. SVA( q ) = SVA( q )*SKL
  567. 3873 CONTINUE
  568. END IF
  569. END IF
  570. 3874 CONTINUE
  571. END IF
  572. *
  573. IF( NOSCALE )SKL = ONE
  574. *
  575. * Move the smaller part of the spectrum from the underflow threshold
  576. *(!) Start by determining the position of the nonzero entries of the
  577. * array SVA() relative to ( SFMIN, BIG ).
  578. *
  579. AAPP = ZERO
  580. AAQQ = BIG
  581. DO 4781 p = 1, N
  582. IF( SVA( p ).NE.ZERO )AAQQ = MIN( AAQQ, SVA( p ) )
  583. AAPP = MAX( AAPP, SVA( p ) )
  584. 4781 CONTINUE
  585. *
  586. * #:) Quick return for zero matrix
  587. *
  588. IF( AAPP.EQ.ZERO ) THEN
  589. IF( LSVEC )CALL SLASET( 'G', M, N, ZERO, ONE, A, LDA )
  590. WORK( 1 ) = ONE
  591. WORK( 2 ) = ZERO
  592. WORK( 3 ) = ZERO
  593. WORK( 4 ) = ZERO
  594. WORK( 5 ) = ZERO
  595. WORK( 6 ) = ZERO
  596. RETURN
  597. END IF
  598. *
  599. * #:) Quick return for one-column matrix
  600. *
  601. IF( N.EQ.1 ) THEN
  602. IF( LSVEC )CALL SLASCL( 'G', 0, 0, SVA( 1 ), SKL, M, 1,
  603. $ A( 1, 1 ), LDA, IERR )
  604. WORK( 1 ) = ONE / SKL
  605. IF( SVA( 1 ).GE.SFMIN ) THEN
  606. WORK( 2 ) = ONE
  607. ELSE
  608. WORK( 2 ) = ZERO
  609. END IF
  610. WORK( 3 ) = ZERO
  611. WORK( 4 ) = ZERO
  612. WORK( 5 ) = ZERO
  613. WORK( 6 ) = ZERO
  614. RETURN
  615. END IF
  616. *
  617. * Protect small singular values from underflow, and try to
  618. * avoid underflows/overflows in computing Jacobi rotations.
  619. *
  620. SN = SQRT( SFMIN / EPSLN )
  621. TEMP1 = SQRT( BIG / FLOAT( N ) )
  622. IF( ( AAPP.LE.SN ) .OR. ( AAQQ.GE.TEMP1 ) .OR.
  623. $ ( ( SN.LE.AAQQ ) .AND. ( AAPP.LE.TEMP1 ) ) ) THEN
  624. TEMP1 = MIN( BIG, TEMP1 / AAPP )
  625. * AAQQ = AAQQ*TEMP1
  626. * AAPP = AAPP*TEMP1
  627. ELSE IF( ( AAQQ.LE.SN ) .AND. ( AAPP.LE.TEMP1 ) ) THEN
  628. TEMP1 = MIN( SN / AAQQ, BIG / ( AAPP*SQRT( FLOAT( N ) ) ) )
  629. * AAQQ = AAQQ*TEMP1
  630. * AAPP = AAPP*TEMP1
  631. ELSE IF( ( AAQQ.GE.SN ) .AND. ( AAPP.GE.TEMP1 ) ) THEN
  632. TEMP1 = MAX( SN / AAQQ, TEMP1 / AAPP )
  633. * AAQQ = AAQQ*TEMP1
  634. * AAPP = AAPP*TEMP1
  635. ELSE IF( ( AAQQ.LE.SN ) .AND. ( AAPP.GE.TEMP1 ) ) THEN
  636. TEMP1 = MIN( SN / AAQQ, BIG / ( SQRT( FLOAT( N ) )*AAPP ) )
  637. * AAQQ = AAQQ*TEMP1
  638. * AAPP = AAPP*TEMP1
  639. ELSE
  640. TEMP1 = ONE
  641. END IF
  642. *
  643. * Scale, if necessary
  644. *
  645. IF( TEMP1.NE.ONE ) THEN
  646. CALL SLASCL( 'G', 0, 0, ONE, TEMP1, N, 1, SVA, N, IERR )
  647. END IF
  648. SKL = TEMP1*SKL
  649. IF( SKL.NE.ONE ) THEN
  650. CALL SLASCL( JOBA, 0, 0, ONE, SKL, M, N, A, LDA, IERR )
  651. SKL = ONE / SKL
  652. END IF
  653. *
  654. * Row-cyclic Jacobi SVD algorithm with column pivoting
  655. *
  656. EMPTSW = ( N*( N-1 ) ) / 2
  657. NOTROT = 0
  658. FASTR( 1 ) = ZERO
  659. *
  660. * A is represented in factored form A = A * diag(WORK), where diag(WORK)
  661. * is initialized to identity. WORK is updated during fast scaled
  662. * rotations.
  663. *
  664. DO 1868 q = 1, N
  665. WORK( q ) = ONE
  666. 1868 CONTINUE
  667. *
  668. *
  669. SWBAND = 3
  670. *[TP] SWBAND is a tuning parameter [TP]. It is meaningful and effective
  671. * if SGESVJ is used as a computational routine in the preconditioned
  672. * Jacobi SVD algorithm SGESVJ. For sweeps i=1:SWBAND the procedure
  673. * works on pivots inside a band-like region around the diagonal.
  674. * The boundaries are determined dynamically, based on the number of
  675. * pivots above a threshold.
  676. *
  677. KBL = MIN( 8, N )
  678. *[TP] KBL is a tuning parameter that defines the tile size in the
  679. * tiling of the p-q loops of pivot pairs. In general, an optimal
  680. * value of KBL depends on the matrix dimensions and on the
  681. * parameters of the computer's memory.
  682. *
  683. NBL = N / KBL
  684. IF( ( NBL*KBL ).NE.N )NBL = NBL + 1
  685. *
  686. BLSKIP = KBL**2
  687. *[TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL.
  688. *
  689. ROWSKIP = MIN( 5, KBL )
  690. *[TP] ROWSKIP is a tuning parameter.
  691. *
  692. LKAHEAD = 1
  693. *[TP] LKAHEAD is a tuning parameter.
  694. *
  695. * Quasi block transformations, using the lower (upper) triangular
  696. * structure of the input matrix. The quasi-block-cycling usually
  697. * invokes cubic convergence. Big part of this cycle is done inside
  698. * canonical subspaces of dimensions less than M.
  699. *
  700. IF( ( LOWER .OR. UPPER ) .AND. ( N.GT.MAX( 64, 4*KBL ) ) ) THEN
  701. *[TP] The number of partition levels and the actual partition are
  702. * tuning parameters.
  703. N4 = N / 4
  704. N2 = N / 2
  705. N34 = 3*N4
  706. IF( APPLV ) THEN
  707. q = 0
  708. ELSE
  709. q = 1
  710. END IF
  711. *
  712. IF( LOWER ) THEN
  713. *
  714. * This works very well on lower triangular matrices, in particular
  715. * in the framework of the preconditioned Jacobi SVD (xGEJSV).
  716. * The idea is simple:
  717. * [+ 0 0 0] Note that Jacobi transformations of [0 0]
  718. * [+ + 0 0] [0 0]
  719. * [+ + x 0] actually work on [x 0] [x 0]
  720. * [+ + x x] [x x]. [x x]
  721. *
  722. CALL SGSVJ0( JOBV, M-N34, N-N34, A( N34+1, N34+1 ), LDA,
  723. $ WORK( N34+1 ), SVA( N34+1 ), MVL,
  724. $ V( N34*q+1, N34+1 ), LDV, EPSLN, SFMIN, TOL,
  725. $ 2, WORK( N+1 ), LWORK-N, IERR )
  726. *
  727. CALL SGSVJ0( JOBV, M-N2, N34-N2, A( N2+1, N2+1 ), LDA,
  728. $ WORK( N2+1 ), SVA( N2+1 ), MVL,
  729. $ V( N2*q+1, N2+1 ), LDV, EPSLN, SFMIN, TOL, 2,
  730. $ WORK( N+1 ), LWORK-N, IERR )
  731. *
  732. CALL SGSVJ1( JOBV, M-N2, N-N2, N4, A( N2+1, N2+1 ), LDA,
  733. $ WORK( N2+1 ), SVA( N2+1 ), MVL,
  734. $ V( N2*q+1, N2+1 ), LDV, EPSLN, SFMIN, TOL, 1,
  735. $ WORK( N+1 ), LWORK-N, IERR )
  736. *
  737. CALL SGSVJ0( JOBV, M-N4, N2-N4, A( N4+1, N4+1 ), LDA,
  738. $ WORK( N4+1 ), SVA( N4+1 ), MVL,
  739. $ V( N4*q+1, N4+1 ), LDV, EPSLN, SFMIN, TOL, 1,
  740. $ WORK( N+1 ), LWORK-N, IERR )
  741. *
  742. CALL SGSVJ0( JOBV, M, N4, A, LDA, WORK, SVA, MVL, V, LDV,
  743. $ EPSLN, SFMIN, TOL, 1, WORK( N+1 ), LWORK-N,
  744. $ IERR )
  745. *
  746. CALL SGSVJ1( JOBV, M, N2, N4, A, LDA, WORK, SVA, MVL, V,
  747. $ LDV, EPSLN, SFMIN, TOL, 1, WORK( N+1 ),
  748. $ LWORK-N, IERR )
  749. *
  750. *
  751. ELSE IF( UPPER ) THEN
  752. *
  753. *
  754. CALL SGSVJ0( JOBV, N4, N4, A, LDA, WORK, SVA, MVL, V, LDV,
  755. $ EPSLN, SFMIN, TOL, 2, WORK( N+1 ), LWORK-N,
  756. $ IERR )
  757. *
  758. CALL SGSVJ0( JOBV, N2, N4, A( 1, N4+1 ), LDA, WORK( N4+1 ),
  759. $ SVA( N4+1 ), MVL, V( N4*q+1, N4+1 ), LDV,
  760. $ EPSLN, SFMIN, TOL, 1, WORK( N+1 ), LWORK-N,
  761. $ IERR )
  762. *
  763. CALL SGSVJ1( JOBV, N2, N2, N4, A, LDA, WORK, SVA, MVL, V,
  764. $ LDV, EPSLN, SFMIN, TOL, 1, WORK( N+1 ),
  765. $ LWORK-N, IERR )
  766. *
  767. CALL SGSVJ0( JOBV, N2+N4, N4, A( 1, N2+1 ), LDA,
  768. $ WORK( N2+1 ), SVA( N2+1 ), MVL,
  769. $ V( N2*q+1, N2+1 ), LDV, EPSLN, SFMIN, TOL, 1,
  770. $ WORK( N+1 ), LWORK-N, IERR )
  771. END IF
  772. *
  773. END IF
  774. *
  775. * .. Row-cyclic pivot strategy with de Rijk's pivoting ..
  776. *
  777. DO 1993 i = 1, NSWEEP
  778. *
  779. * .. go go go ...
  780. *
  781. MXAAPQ = ZERO
  782. MXSINJ = ZERO
  783. ISWROT = 0
  784. *
  785. NOTROT = 0
  786. PSKIPPED = 0
  787. *
  788. * Each sweep is unrolled using KBL-by-KBL tiles over the pivot pairs
  789. * 1 <= p < q <= N. This is the first step toward a blocked implementation
  790. * of the rotations. New implementation, based on block transformations,
  791. * is under development.
  792. *
  793. DO 2000 ibr = 1, NBL
  794. *
  795. igl = ( ibr-1 )*KBL + 1
  796. *
  797. DO 1002 ir1 = 0, MIN( LKAHEAD, NBL-ibr )
  798. *
  799. igl = igl + ir1*KBL
  800. *
  801. DO 2001 p = igl, MIN( igl+KBL-1, N-1 )
  802. *
  803. * .. de Rijk's pivoting
  804. *
  805. q = ISAMAX( N-p+1, SVA( p ), 1 ) + p - 1
  806. IF( p.NE.q ) THEN
  807. CALL SSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
  808. IF( RSVEC )CALL SSWAP( MVL, V( 1, p ), 1,
  809. $ V( 1, q ), 1 )
  810. TEMP1 = SVA( p )
  811. SVA( p ) = SVA( q )
  812. SVA( q ) = TEMP1
  813. TEMP1 = WORK( p )
  814. WORK( p ) = WORK( q )
  815. WORK( q ) = TEMP1
  816. END IF
  817. *
  818. IF( ir1.EQ.0 ) THEN
  819. *
  820. * Column norms are periodically updated by explicit
  821. * norm computation.
  822. * Caveat:
  823. * Unfortunately, some BLAS implementations compute SNRM2(M,A(1,p),1)
  824. * as SQRT(SDOT(M,A(1,p),1,A(1,p),1)), which may cause the result to
  825. * overflow for ||A(:,p)||_2 > SQRT(overflow_threshold), and to
  826. * underflow for ||A(:,p)||_2 < SQRT(underflow_threshold).
  827. * Hence, SNRM2 cannot be trusted, not even in the case when
  828. * the true norm is far from the under(over)flow boundaries.
  829. * If properly implemented SNRM2 is available, the IF-THEN-ELSE
  830. * below should read "AAPP = SNRM2( M, A(1,p), 1 ) * WORK(p)".
  831. *
  832. IF( ( SVA( p ).LT.ROOTBIG ) .AND.
  833. $ ( SVA( p ).GT.ROOTSFMIN ) ) THEN
  834. SVA( p ) = SNRM2( M, A( 1, p ), 1 )*WORK( p )
  835. ELSE
  836. TEMP1 = ZERO
  837. AAPP = ONE
  838. CALL SLASSQ( M, A( 1, p ), 1, TEMP1, AAPP )
  839. SVA( p ) = TEMP1*SQRT( AAPP )*WORK( p )
  840. END IF
  841. AAPP = SVA( p )
  842. ELSE
  843. AAPP = SVA( p )
  844. END IF
  845. *
  846. IF( AAPP.GT.ZERO ) THEN
  847. *
  848. PSKIPPED = 0
  849. *
  850. DO 2002 q = p + 1, MIN( igl+KBL-1, N )
  851. *
  852. AAQQ = SVA( q )
  853. *
  854. IF( AAQQ.GT.ZERO ) THEN
  855. *
  856. AAPP0 = AAPP
  857. IF( AAQQ.GE.ONE ) THEN
  858. ROTOK = ( SMALL*AAPP ).LE.AAQQ
  859. IF( AAPP.LT.( BIG / AAQQ ) ) THEN
  860. AAPQ = ( SDOT( M, A( 1, p ), 1, A( 1,
  861. $ q ), 1 )*WORK( p )*WORK( q ) /
  862. $ AAQQ ) / AAPP
  863. ELSE
  864. CALL SCOPY( M, A( 1, p ), 1,
  865. $ WORK( N+1 ), 1 )
  866. CALL SLASCL( 'G', 0, 0, AAPP,
  867. $ WORK( p ), M, 1,
  868. $ WORK( N+1 ), LDA, IERR )
  869. AAPQ = SDOT( M, WORK( N+1 ), 1,
  870. $ A( 1, q ), 1 )*WORK( q ) / AAQQ
  871. END IF
  872. ELSE
  873. ROTOK = AAPP.LE.( AAQQ / SMALL )
  874. IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
  875. AAPQ = ( SDOT( M, A( 1, p ), 1, A( 1,
  876. $ q ), 1 )*WORK( p )*WORK( q ) /
  877. $ AAQQ ) / AAPP
  878. ELSE
  879. CALL SCOPY( M, A( 1, q ), 1,
  880. $ WORK( N+1 ), 1 )
  881. CALL SLASCL( 'G', 0, 0, AAQQ,
  882. $ WORK( q ), M, 1,
  883. $ WORK( N+1 ), LDA, IERR )
  884. AAPQ = SDOT( M, WORK( N+1 ), 1,
  885. $ A( 1, p ), 1 )*WORK( p ) / AAPP
  886. END IF
  887. END IF
  888. *
  889. MXAAPQ = MAX( MXAAPQ, ABS( AAPQ ) )
  890. *
  891. * TO rotate or NOT to rotate, THAT is the question ...
  892. *
  893. IF( ABS( AAPQ ).GT.TOL ) THEN
  894. *
  895. * .. rotate
  896. *[RTD] ROTATED = ROTATED + ONE
  897. *
  898. IF( ir1.EQ.0 ) THEN
  899. NOTROT = 0
  900. PSKIPPED = 0
  901. ISWROT = ISWROT + 1
  902. END IF
  903. *
  904. IF( ROTOK ) THEN
  905. *
  906. AQOAP = AAQQ / AAPP
  907. APOAQ = AAPP / AAQQ
  908. THETA = -HALF*ABS( AQOAP-APOAQ ) / AAPQ
  909. *
  910. IF( ABS( THETA ).GT.BIGTHETA ) THEN
  911. *
  912. T = HALF / THETA
  913. FASTR( 3 ) = T*WORK( p ) / WORK( q )
  914. FASTR( 4 ) = -T*WORK( q ) /
  915. $ WORK( p )
  916. CALL SROTM( M, A( 1, p ), 1,
  917. $ A( 1, q ), 1, FASTR )
  918. IF( RSVEC )CALL SROTM( MVL,
  919. $ V( 1, p ), 1,
  920. $ V( 1, q ), 1,
  921. $ FASTR )
  922. SVA( q ) = AAQQ*SQRT( MAX( ZERO,
  923. $ ONE+T*APOAQ*AAPQ ) )
  924. AAPP = AAPP*SQRT( MAX( ZERO,
  925. $ ONE-T*AQOAP*AAPQ ) )
  926. MXSINJ = MAX( MXSINJ, ABS( T ) )
  927. *
  928. ELSE
  929. *
  930. * .. choose correct signum for THETA and rotate
  931. *
  932. THSIGN = -SIGN( ONE, AAPQ )
  933. T = ONE / ( THETA+THSIGN*
  934. $ SQRT( ONE+THETA*THETA ) )
  935. CS = SQRT( ONE / ( ONE+T*T ) )
  936. SN = T*CS
  937. *
  938. MXSINJ = MAX( MXSINJ, ABS( SN ) )
  939. SVA( q ) = AAQQ*SQRT( MAX( ZERO,
  940. $ ONE+T*APOAQ*AAPQ ) )
  941. AAPP = AAPP*SQRT( MAX( ZERO,
  942. $ ONE-T*AQOAP*AAPQ ) )
  943. *
  944. APOAQ = WORK( p ) / WORK( q )
  945. AQOAP = WORK( q ) / WORK( p )
  946. IF( WORK( p ).GE.ONE ) THEN
  947. IF( WORK( q ).GE.ONE ) THEN
  948. FASTR( 3 ) = T*APOAQ
  949. FASTR( 4 ) = -T*AQOAP
  950. WORK( p ) = WORK( p )*CS
  951. WORK( q ) = WORK( q )*CS
  952. CALL SROTM( M, A( 1, p ), 1,
  953. $ A( 1, q ), 1,
  954. $ FASTR )
  955. IF( RSVEC )CALL SROTM( MVL,
  956. $ V( 1, p ), 1, V( 1, q ),
  957. $ 1, FASTR )
  958. ELSE
  959. CALL SAXPY( M, -T*AQOAP,
  960. $ A( 1, q ), 1,
  961. $ A( 1, p ), 1 )
  962. CALL SAXPY( M, CS*SN*APOAQ,
  963. $ A( 1, p ), 1,
  964. $ A( 1, q ), 1 )
  965. WORK( p ) = WORK( p )*CS
  966. WORK( q ) = WORK( q ) / CS
  967. IF( RSVEC ) THEN
  968. CALL SAXPY( MVL, -T*AQOAP,
  969. $ V( 1, q ), 1,
  970. $ V( 1, p ), 1 )
  971. CALL SAXPY( MVL,
  972. $ CS*SN*APOAQ,
  973. $ V( 1, p ), 1,
  974. $ V( 1, q ), 1 )
  975. END IF
  976. END IF
  977. ELSE
  978. IF( WORK( q ).GE.ONE ) THEN
  979. CALL SAXPY( M, T*APOAQ,
  980. $ A( 1, p ), 1,
  981. $ A( 1, q ), 1 )
  982. CALL SAXPY( M, -CS*SN*AQOAP,
  983. $ A( 1, q ), 1,
  984. $ A( 1, p ), 1 )
  985. WORK( p ) = WORK( p ) / CS
  986. WORK( q ) = WORK( q )*CS
  987. IF( RSVEC ) THEN
  988. CALL SAXPY( MVL, T*APOAQ,
  989. $ V( 1, p ), 1,
  990. $ V( 1, q ), 1 )
  991. CALL SAXPY( MVL,
  992. $ -CS*SN*AQOAP,
  993. $ V( 1, q ), 1,
  994. $ V( 1, p ), 1 )
  995. END IF
  996. ELSE
  997. IF( WORK( p ).GE.WORK( q ) )
  998. $ THEN
  999. CALL SAXPY( M, -T*AQOAP,
  1000. $ A( 1, q ), 1,
  1001. $ A( 1, p ), 1 )
  1002. CALL SAXPY( M, CS*SN*APOAQ,
  1003. $ A( 1, p ), 1,
  1004. $ A( 1, q ), 1 )
  1005. WORK( p ) = WORK( p )*CS
  1006. WORK( q ) = WORK( q ) / CS
  1007. IF( RSVEC ) THEN
  1008. CALL SAXPY( MVL,
  1009. $ -T*AQOAP,
  1010. $ V( 1, q ), 1,
  1011. $ V( 1, p ), 1 )
  1012. CALL SAXPY( MVL,
  1013. $ CS*SN*APOAQ,
  1014. $ V( 1, p ), 1,
  1015. $ V( 1, q ), 1 )
  1016. END IF
  1017. ELSE
  1018. CALL SAXPY( M, T*APOAQ,
  1019. $ A( 1, p ), 1,
  1020. $ A( 1, q ), 1 )
  1021. CALL SAXPY( M,
  1022. $ -CS*SN*AQOAP,
  1023. $ A( 1, q ), 1,
  1024. $ A( 1, p ), 1 )
  1025. WORK( p ) = WORK( p ) / CS
  1026. WORK( q ) = WORK( q )*CS
  1027. IF( RSVEC ) THEN
  1028. CALL SAXPY( MVL,
  1029. $ T*APOAQ, V( 1, p ),
  1030. $ 1, V( 1, q ), 1 )
  1031. CALL SAXPY( MVL,
  1032. $ -CS*SN*AQOAP,
  1033. $ V( 1, q ), 1,
  1034. $ V( 1, p ), 1 )
  1035. END IF
  1036. END IF
  1037. END IF
  1038. END IF
  1039. END IF
  1040. *
  1041. ELSE
  1042. * .. have to use modified Gram-Schmidt like transformation
  1043. CALL SCOPY( M, A( 1, p ), 1,
  1044. $ WORK( N+1 ), 1 )
  1045. CALL SLASCL( 'G', 0, 0, AAPP, ONE, M,
  1046. $ 1, WORK( N+1 ), LDA,
  1047. $ IERR )
  1048. CALL SLASCL( 'G', 0, 0, AAQQ, ONE, M,
  1049. $ 1, A( 1, q ), LDA, IERR )
  1050. TEMP1 = -AAPQ*WORK( p ) / WORK( q )
  1051. CALL SAXPY( M, TEMP1, WORK( N+1 ), 1,
  1052. $ A( 1, q ), 1 )
  1053. CALL SLASCL( 'G', 0, 0, ONE, AAQQ, M,
  1054. $ 1, A( 1, q ), LDA, IERR )
  1055. SVA( q ) = AAQQ*SQRT( MAX( ZERO,
  1056. $ ONE-AAPQ*AAPQ ) )
  1057. MXSINJ = MAX( MXSINJ, SFMIN )
  1058. END IF
  1059. * END IF ROTOK THEN ... ELSE
  1060. *
  1061. * In the case of cancellation in updating SVA(q), SVA(p)
  1062. * recompute SVA(q), SVA(p).
  1063. *
  1064. IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
  1065. $ THEN
  1066. IF( ( AAQQ.LT.ROOTBIG ) .AND.
  1067. $ ( AAQQ.GT.ROOTSFMIN ) ) THEN
  1068. SVA( q ) = SNRM2( M, A( 1, q ), 1 )*
  1069. $ WORK( q )
  1070. ELSE
  1071. T = ZERO
  1072. AAQQ = ONE
  1073. CALL SLASSQ( M, A( 1, q ), 1, T,
  1074. $ AAQQ )
  1075. SVA( q ) = T*SQRT( AAQQ )*WORK( q )
  1076. END IF
  1077. END IF
  1078. IF( ( AAPP / AAPP0 ).LE.ROOTEPS ) THEN
  1079. IF( ( AAPP.LT.ROOTBIG ) .AND.
  1080. $ ( AAPP.GT.ROOTSFMIN ) ) THEN
  1081. AAPP = SNRM2( M, A( 1, p ), 1 )*
  1082. $ WORK( p )
  1083. ELSE
  1084. T = ZERO
  1085. AAPP = ONE
  1086. CALL SLASSQ( M, A( 1, p ), 1, T,
  1087. $ AAPP )
  1088. AAPP = T*SQRT( AAPP )*WORK( p )
  1089. END IF
  1090. SVA( p ) = AAPP
  1091. END IF
  1092. *
  1093. ELSE
  1094. * A(:,p) and A(:,q) already numerically orthogonal
  1095. IF( ir1.EQ.0 )NOTROT = NOTROT + 1
  1096. *[RTD] SKIPPED = SKIPPED + 1
  1097. PSKIPPED = PSKIPPED + 1
  1098. END IF
  1099. ELSE
  1100. * A(:,q) is zero column
  1101. IF( ir1.EQ.0 )NOTROT = NOTROT + 1
  1102. PSKIPPED = PSKIPPED + 1
  1103. END IF
  1104. *
  1105. IF( ( i.LE.SWBAND ) .AND.
  1106. $ ( PSKIPPED.GT.ROWSKIP ) ) THEN
  1107. IF( ir1.EQ.0 )AAPP = -AAPP
  1108. NOTROT = 0
  1109. GO TO 2103
  1110. END IF
  1111. *
  1112. 2002 CONTINUE
  1113. * END q-LOOP
  1114. *
  1115. 2103 CONTINUE
  1116. * bailed out of q-loop
  1117. *
  1118. SVA( p ) = AAPP
  1119. *
  1120. ELSE
  1121. SVA( p ) = AAPP
  1122. IF( ( ir1.EQ.0 ) .AND. ( AAPP.EQ.ZERO ) )
  1123. $ NOTROT = NOTROT + MIN( igl+KBL-1, N ) - p
  1124. END IF
  1125. *
  1126. 2001 CONTINUE
  1127. * end of the p-loop
  1128. * end of doing the block ( ibr, ibr )
  1129. 1002 CONTINUE
  1130. * end of ir1-loop
  1131. *
  1132. * ... go to the off diagonal blocks
  1133. *
  1134. igl = ( ibr-1 )*KBL + 1
  1135. *
  1136. DO 2010 jbc = ibr + 1, NBL
  1137. *
  1138. jgl = ( jbc-1 )*KBL + 1
  1139. *
  1140. * doing the block at ( ibr, jbc )
  1141. *
  1142. IJBLSK = 0
  1143. DO 2100 p = igl, MIN( igl+KBL-1, N )
  1144. *
  1145. AAPP = SVA( p )
  1146. IF( AAPP.GT.ZERO ) THEN
  1147. *
  1148. PSKIPPED = 0
  1149. *
  1150. DO 2200 q = jgl, MIN( jgl+KBL-1, N )
  1151. *
  1152. AAQQ = SVA( q )
  1153. IF( AAQQ.GT.ZERO ) THEN
  1154. AAPP0 = AAPP
  1155. *
  1156. * .. M x 2 Jacobi SVD ..
  1157. *
  1158. * Safe Gram matrix computation
  1159. *
  1160. IF( AAQQ.GE.ONE ) THEN
  1161. IF( AAPP.GE.AAQQ ) THEN
  1162. ROTOK = ( SMALL*AAPP ).LE.AAQQ
  1163. ELSE
  1164. ROTOK = ( SMALL*AAQQ ).LE.AAPP
  1165. END IF
  1166. IF( AAPP.LT.( BIG / AAQQ ) ) THEN
  1167. AAPQ = ( SDOT( M, A( 1, p ), 1, A( 1,
  1168. $ q ), 1 )*WORK( p )*WORK( q ) /
  1169. $ AAQQ ) / AAPP
  1170. ELSE
  1171. CALL SCOPY( M, A( 1, p ), 1,
  1172. $ WORK( N+1 ), 1 )
  1173. CALL SLASCL( 'G', 0, 0, AAPP,
  1174. $ WORK( p ), M, 1,
  1175. $ WORK( N+1 ), LDA, IERR )
  1176. AAPQ = SDOT( M, WORK( N+1 ), 1,
  1177. $ A( 1, q ), 1 )*WORK( q ) / AAQQ
  1178. END IF
  1179. ELSE
  1180. IF( AAPP.GE.AAQQ ) THEN
  1181. ROTOK = AAPP.LE.( AAQQ / SMALL )
  1182. ELSE
  1183. ROTOK = AAQQ.LE.( AAPP / SMALL )
  1184. END IF
  1185. IF( AAPP.GT.( SMALL / AAQQ ) ) THEN
  1186. AAPQ = ( SDOT( M, A( 1, p ), 1, A( 1,
  1187. $ q ), 1 )*WORK( p )*WORK( q ) /
  1188. $ AAQQ ) / AAPP
  1189. ELSE
  1190. CALL SCOPY( M, A( 1, q ), 1,
  1191. $ WORK( N+1 ), 1 )
  1192. CALL SLASCL( 'G', 0, 0, AAQQ,
  1193. $ WORK( q ), M, 1,
  1194. $ WORK( N+1 ), LDA, IERR )
  1195. AAPQ = SDOT( M, WORK( N+1 ), 1,
  1196. $ A( 1, p ), 1 )*WORK( p ) / AAPP
  1197. END IF
  1198. END IF
  1199. *
  1200. MXAAPQ = MAX( MXAAPQ, ABS( AAPQ ) )
  1201. *
  1202. * TO rotate or NOT to rotate, THAT is the question ...
  1203. *
  1204. IF( ABS( AAPQ ).GT.TOL ) THEN
  1205. NOTROT = 0
  1206. *[RTD] ROTATED = ROTATED + 1
  1207. PSKIPPED = 0
  1208. ISWROT = ISWROT + 1
  1209. *
  1210. IF( ROTOK ) THEN
  1211. *
  1212. AQOAP = AAQQ / AAPP
  1213. APOAQ = AAPP / AAQQ
  1214. THETA = -HALF*ABS( AQOAP-APOAQ ) / AAPQ
  1215. IF( AAQQ.GT.AAPP0 )THETA = -THETA
  1216. *
  1217. IF( ABS( THETA ).GT.BIGTHETA ) THEN
  1218. T = HALF / THETA
  1219. FASTR( 3 ) = T*WORK( p ) / WORK( q )
  1220. FASTR( 4 ) = -T*WORK( q ) /
  1221. $ WORK( p )
  1222. CALL SROTM( M, A( 1, p ), 1,
  1223. $ A( 1, q ), 1, FASTR )
  1224. IF( RSVEC )CALL SROTM( MVL,
  1225. $ V( 1, p ), 1,
  1226. $ V( 1, q ), 1,
  1227. $ FASTR )
  1228. SVA( q ) = AAQQ*SQRT( MAX( ZERO,
  1229. $ ONE+T*APOAQ*AAPQ ) )
  1230. AAPP = AAPP*SQRT( MAX( ZERO,
  1231. $ ONE-T*AQOAP*AAPQ ) )
  1232. MXSINJ = MAX( MXSINJ, ABS( T ) )
  1233. ELSE
  1234. *
  1235. * .. choose correct signum for THETA and rotate
  1236. *
  1237. THSIGN = -SIGN( ONE, AAPQ )
  1238. IF( AAQQ.GT.AAPP0 )THSIGN = -THSIGN
  1239. T = ONE / ( THETA+THSIGN*
  1240. $ SQRT( ONE+THETA*THETA ) )
  1241. CS = SQRT( ONE / ( ONE+T*T ) )
  1242. SN = T*CS
  1243. MXSINJ = MAX( MXSINJ, ABS( SN ) )
  1244. SVA( q ) = AAQQ*SQRT( MAX( ZERO,
  1245. $ ONE+T*APOAQ*AAPQ ) )
  1246. AAPP = AAPP*SQRT( MAX( ZERO,
  1247. $ ONE-T*AQOAP*AAPQ ) )
  1248. *
  1249. APOAQ = WORK( p ) / WORK( q )
  1250. AQOAP = WORK( q ) / WORK( p )
  1251. IF( WORK( p ).GE.ONE ) THEN
  1252. *
  1253. IF( WORK( q ).GE.ONE ) THEN
  1254. FASTR( 3 ) = T*APOAQ
  1255. FASTR( 4 ) = -T*AQOAP
  1256. WORK( p ) = WORK( p )*CS
  1257. WORK( q ) = WORK( q )*CS
  1258. CALL SROTM( M, A( 1, p ), 1,
  1259. $ A( 1, q ), 1,
  1260. $ FASTR )
  1261. IF( RSVEC )CALL SROTM( MVL,
  1262. $ V( 1, p ), 1, V( 1, q ),
  1263. $ 1, FASTR )
  1264. ELSE
  1265. CALL SAXPY( M, -T*AQOAP,
  1266. $ A( 1, q ), 1,
  1267. $ A( 1, p ), 1 )
  1268. CALL SAXPY( M, CS*SN*APOAQ,
  1269. $ A( 1, p ), 1,
  1270. $ A( 1, q ), 1 )
  1271. IF( RSVEC ) THEN
  1272. CALL SAXPY( MVL, -T*AQOAP,
  1273. $ V( 1, q ), 1,
  1274. $ V( 1, p ), 1 )
  1275. CALL SAXPY( MVL,
  1276. $ CS*SN*APOAQ,
  1277. $ V( 1, p ), 1,
  1278. $ V( 1, q ), 1 )
  1279. END IF
  1280. WORK( p ) = WORK( p )*CS
  1281. WORK( q ) = WORK( q ) / CS
  1282. END IF
  1283. ELSE
  1284. IF( WORK( q ).GE.ONE ) THEN
  1285. CALL SAXPY( M, T*APOAQ,
  1286. $ A( 1, p ), 1,
  1287. $ A( 1, q ), 1 )
  1288. CALL SAXPY( M, -CS*SN*AQOAP,
  1289. $ A( 1, q ), 1,
  1290. $ A( 1, p ), 1 )
  1291. IF( RSVEC ) THEN
  1292. CALL SAXPY( MVL, T*APOAQ,
  1293. $ V( 1, p ), 1,
  1294. $ V( 1, q ), 1 )
  1295. CALL SAXPY( MVL,
  1296. $ -CS*SN*AQOAP,
  1297. $ V( 1, q ), 1,
  1298. $ V( 1, p ), 1 )
  1299. END IF
  1300. WORK( p ) = WORK( p ) / CS
  1301. WORK( q ) = WORK( q )*CS
  1302. ELSE
  1303. IF( WORK( p ).GE.WORK( q ) )
  1304. $ THEN
  1305. CALL SAXPY( M, -T*AQOAP,
  1306. $ A( 1, q ), 1,
  1307. $ A( 1, p ), 1 )
  1308. CALL SAXPY( M, CS*SN*APOAQ,
  1309. $ A( 1, p ), 1,
  1310. $ A( 1, q ), 1 )
  1311. WORK( p ) = WORK( p )*CS
  1312. WORK( q ) = WORK( q ) / CS
  1313. IF( RSVEC ) THEN
  1314. CALL SAXPY( MVL,
  1315. $ -T*AQOAP,
  1316. $ V( 1, q ), 1,
  1317. $ V( 1, p ), 1 )
  1318. CALL SAXPY( MVL,
  1319. $ CS*SN*APOAQ,
  1320. $ V( 1, p ), 1,
  1321. $ V( 1, q ), 1 )
  1322. END IF
  1323. ELSE
  1324. CALL SAXPY( M, T*APOAQ,
  1325. $ A( 1, p ), 1,
  1326. $ A( 1, q ), 1 )
  1327. CALL SAXPY( M,
  1328. $ -CS*SN*AQOAP,
  1329. $ A( 1, q ), 1,
  1330. $ A( 1, p ), 1 )
  1331. WORK( p ) = WORK( p ) / CS
  1332. WORK( q ) = WORK( q )*CS
  1333. IF( RSVEC ) THEN
  1334. CALL SAXPY( MVL,
  1335. $ T*APOAQ, V( 1, p ),
  1336. $ 1, V( 1, q ), 1 )
  1337. CALL SAXPY( MVL,
  1338. $ -CS*SN*AQOAP,
  1339. $ V( 1, q ), 1,
  1340. $ V( 1, p ), 1 )
  1341. END IF
  1342. END IF
  1343. END IF
  1344. END IF
  1345. END IF
  1346. *
  1347. ELSE
  1348. IF( AAPP.GT.AAQQ ) THEN
  1349. CALL SCOPY( M, A( 1, p ), 1,
  1350. $ WORK( N+1 ), 1 )
  1351. CALL SLASCL( 'G', 0, 0, AAPP, ONE,
  1352. $ M, 1, WORK( N+1 ), LDA,
  1353. $ IERR )
  1354. CALL SLASCL( 'G', 0, 0, AAQQ, ONE,
  1355. $ M, 1, A( 1, q ), LDA,
  1356. $ IERR )
  1357. TEMP1 = -AAPQ*WORK( p ) / WORK( q )
  1358. CALL SAXPY( M, TEMP1, WORK( N+1 ),
  1359. $ 1, A( 1, q ), 1 )
  1360. CALL SLASCL( 'G', 0, 0, ONE, AAQQ,
  1361. $ M, 1, A( 1, q ), LDA,
  1362. $ IERR )
  1363. SVA( q ) = AAQQ*SQRT( MAX( ZERO,
  1364. $ ONE-AAPQ*AAPQ ) )
  1365. MXSINJ = MAX( MXSINJ, SFMIN )
  1366. ELSE
  1367. CALL SCOPY( M, A( 1, q ), 1,
  1368. $ WORK( N+1 ), 1 )
  1369. CALL SLASCL( 'G', 0, 0, AAQQ, ONE,
  1370. $ M, 1, WORK( N+1 ), LDA,
  1371. $ IERR )
  1372. CALL SLASCL( 'G', 0, 0, AAPP, ONE,
  1373. $ M, 1, A( 1, p ), LDA,
  1374. $ IERR )
  1375. TEMP1 = -AAPQ*WORK( q ) / WORK( p )
  1376. CALL SAXPY( M, TEMP1, WORK( N+1 ),
  1377. $ 1, A( 1, p ), 1 )
  1378. CALL SLASCL( 'G', 0, 0, ONE, AAPP,
  1379. $ M, 1, A( 1, p ), LDA,
  1380. $ IERR )
  1381. SVA( p ) = AAPP*SQRT( MAX( ZERO,
  1382. $ ONE-AAPQ*AAPQ ) )
  1383. MXSINJ = MAX( MXSINJ, SFMIN )
  1384. END IF
  1385. END IF
  1386. * END IF ROTOK THEN ... ELSE
  1387. *
  1388. * In the case of cancellation in updating SVA(q)
  1389. * .. recompute SVA(q)
  1390. IF( ( SVA( q ) / AAQQ )**2.LE.ROOTEPS )
  1391. $ THEN
  1392. IF( ( AAQQ.LT.ROOTBIG ) .AND.
  1393. $ ( AAQQ.GT.ROOTSFMIN ) ) THEN
  1394. SVA( q ) = SNRM2( M, A( 1, q ), 1 )*
  1395. $ WORK( q )
  1396. ELSE
  1397. T = ZERO
  1398. AAQQ = ONE
  1399. CALL SLASSQ( M, A( 1, q ), 1, T,
  1400. $ AAQQ )
  1401. SVA( q ) = T*SQRT( AAQQ )*WORK( q )
  1402. END IF
  1403. END IF
  1404. IF( ( AAPP / AAPP0 )**2.LE.ROOTEPS ) THEN
  1405. IF( ( AAPP.LT.ROOTBIG ) .AND.
  1406. $ ( AAPP.GT.ROOTSFMIN ) ) THEN
  1407. AAPP = SNRM2( M, A( 1, p ), 1 )*
  1408. $ WORK( p )
  1409. ELSE
  1410. T = ZERO
  1411. AAPP = ONE
  1412. CALL SLASSQ( M, A( 1, p ), 1, T,
  1413. $ AAPP )
  1414. AAPP = T*SQRT( AAPP )*WORK( p )
  1415. END IF
  1416. SVA( p ) = AAPP
  1417. END IF
  1418. * end of OK rotation
  1419. ELSE
  1420. NOTROT = NOTROT + 1
  1421. *[RTD] SKIPPED = SKIPPED + 1
  1422. PSKIPPED = PSKIPPED + 1
  1423. IJBLSK = IJBLSK + 1
  1424. END IF
  1425. ELSE
  1426. NOTROT = NOTROT + 1
  1427. PSKIPPED = PSKIPPED + 1
  1428. IJBLSK = IJBLSK + 1
  1429. END IF
  1430. *
  1431. IF( ( i.LE.SWBAND ) .AND. ( IJBLSK.GE.BLSKIP ) )
  1432. $ THEN
  1433. SVA( p ) = AAPP
  1434. NOTROT = 0
  1435. GO TO 2011
  1436. END IF
  1437. IF( ( i.LE.SWBAND ) .AND.
  1438. $ ( PSKIPPED.GT.ROWSKIP ) ) THEN
  1439. AAPP = -AAPP
  1440. NOTROT = 0
  1441. GO TO 2203
  1442. END IF
  1443. *
  1444. 2200 CONTINUE
  1445. * end of the q-loop
  1446. 2203 CONTINUE
  1447. *
  1448. SVA( p ) = AAPP
  1449. *
  1450. ELSE
  1451. *
  1452. IF( AAPP.EQ.ZERO )NOTROT = NOTROT +
  1453. $ MIN( jgl+KBL-1, N ) - jgl + 1
  1454. IF( AAPP.LT.ZERO )NOTROT = 0
  1455. *
  1456. END IF
  1457. *
  1458. 2100 CONTINUE
  1459. * end of the p-loop
  1460. 2010 CONTINUE
  1461. * end of the jbc-loop
  1462. 2011 CONTINUE
  1463. *2011 bailed out of the jbc-loop
  1464. DO 2012 p = igl, MIN( igl+KBL-1, N )
  1465. SVA( p ) = ABS( SVA( p ) )
  1466. 2012 CONTINUE
  1467. ***
  1468. 2000 CONTINUE
  1469. *2000 :: end of the ibr-loop
  1470. *
  1471. * .. update SVA(N)
  1472. IF( ( SVA( N ).LT.ROOTBIG ) .AND. ( SVA( N ).GT.ROOTSFMIN ) )
  1473. $ THEN
  1474. SVA( N ) = SNRM2( M, A( 1, N ), 1 )*WORK( N )
  1475. ELSE
  1476. T = ZERO
  1477. AAPP = ONE
  1478. CALL SLASSQ( M, A( 1, N ), 1, T, AAPP )
  1479. SVA( N ) = T*SQRT( AAPP )*WORK( N )
  1480. END IF
  1481. *
  1482. * Additional steering devices
  1483. *
  1484. IF( ( i.LT.SWBAND ) .AND. ( ( MXAAPQ.LE.ROOTTOL ) .OR.
  1485. $ ( ISWROT.LE.N ) ) )SWBAND = i
  1486. *
  1487. IF( ( i.GT.SWBAND+1 ) .AND. ( MXAAPQ.LT.SQRT( FLOAT( N ) )*
  1488. $ TOL ) .AND. ( FLOAT( N )*MXAAPQ*MXSINJ.LT.TOL ) ) THEN
  1489. GO TO 1994
  1490. END IF
  1491. *
  1492. IF( NOTROT.GE.EMPTSW )GO TO 1994
  1493. *
  1494. 1993 CONTINUE
  1495. * end i=1:NSWEEP loop
  1496. *
  1497. * #:( Reaching this point means that the procedure has not converged.
  1498. INFO = NSWEEP - 1
  1499. GO TO 1995
  1500. *
  1501. 1994 CONTINUE
  1502. * #:) Reaching this point means numerical convergence after the i-th
  1503. * sweep.
  1504. *
  1505. INFO = 0
  1506. * #:) INFO = 0 confirms successful iterations.
  1507. 1995 CONTINUE
  1508. *
  1509. * Sort the singular values and find how many are above
  1510. * the underflow threshold.
  1511. *
  1512. N2 = 0
  1513. N4 = 0
  1514. DO 5991 p = 1, N - 1
  1515. q = ISAMAX( N-p+1, SVA( p ), 1 ) + p - 1
  1516. IF( p.NE.q ) THEN
  1517. TEMP1 = SVA( p )
  1518. SVA( p ) = SVA( q )
  1519. SVA( q ) = TEMP1
  1520. TEMP1 = WORK( p )
  1521. WORK( p ) = WORK( q )
  1522. WORK( q ) = TEMP1
  1523. CALL SSWAP( M, A( 1, p ), 1, A( 1, q ), 1 )
  1524. IF( RSVEC )CALL SSWAP( MVL, V( 1, p ), 1, V( 1, q ), 1 )
  1525. END IF
  1526. IF( SVA( p ).NE.ZERO ) THEN
  1527. N4 = N4 + 1
  1528. IF( SVA( p )*SKL.GT.SFMIN )N2 = N2 + 1
  1529. END IF
  1530. 5991 CONTINUE
  1531. IF( SVA( N ).NE.ZERO ) THEN
  1532. N4 = N4 + 1
  1533. IF( SVA( N )*SKL.GT.SFMIN )N2 = N2 + 1
  1534. END IF
  1535. *
  1536. * Normalize the left singular vectors.
  1537. *
  1538. IF( LSVEC .OR. UCTOL ) THEN
  1539. DO 1998 p = 1, N2
  1540. CALL SSCAL( M, WORK( p ) / SVA( p ), A( 1, p ), 1 )
  1541. 1998 CONTINUE
  1542. END IF
  1543. *
  1544. * Scale the product of Jacobi rotations (assemble the fast rotations).
  1545. *
  1546. IF( RSVEC ) THEN
  1547. IF( APPLV ) THEN
  1548. DO 2398 p = 1, N
  1549. CALL SSCAL( MVL, WORK( p ), V( 1, p ), 1 )
  1550. 2398 CONTINUE
  1551. ELSE
  1552. DO 2399 p = 1, N
  1553. TEMP1 = ONE / SNRM2( MVL, V( 1, p ), 1 )
  1554. CALL SSCAL( MVL, TEMP1, V( 1, p ), 1 )
  1555. 2399 CONTINUE
  1556. END IF
  1557. END IF
  1558. *
  1559. * Undo scaling, if necessary (and possible).
  1560. IF( ( ( SKL.GT.ONE ) .AND. ( SVA( 1 ).LT.( BIG / SKL ) ) )
  1561. $ .OR. ( ( SKL.LT.ONE ) .AND. ( SVA( MAX( N2, 1 ) ) .GT.
  1562. $ ( SFMIN / SKL ) ) ) ) THEN
  1563. DO 2400 p = 1, N
  1564. SVA( P ) = SKL*SVA( P )
  1565. 2400 CONTINUE
  1566. SKL = ONE
  1567. END IF
  1568. *
  1569. WORK( 1 ) = SKL
  1570. * The singular values of A are SKL*SVA(1:N). If SKL.NE.ONE
  1571. * then some of the singular values may overflow or underflow and
  1572. * the spectrum is given in this factored representation.
  1573. *
  1574. WORK( 2 ) = FLOAT( N4 )
  1575. * N4 is the number of computed nonzero singular values of A.
  1576. *
  1577. WORK( 3 ) = FLOAT( N2 )
  1578. * N2 is the number of singular values of A greater than SFMIN.
  1579. * If N2<N, SVA(N2:N) contains ZEROS and/or denormalized numbers
  1580. * that may carry some information.
  1581. *
  1582. WORK( 4 ) = FLOAT( i )
  1583. * i is the index of the last sweep before declaring convergence.
  1584. *
  1585. WORK( 5 ) = MXAAPQ
  1586. * MXAAPQ is the largest absolute value of scaled pivots in the
  1587. * last sweep
  1588. *
  1589. WORK( 6 ) = MXSINJ
  1590. * MXSINJ is the largest absolute value of the sines of Jacobi angles
  1591. * in the last sweep
  1592. *
  1593. RETURN
  1594. * ..
  1595. * .. END OF SGESVJ
  1596. * ..
  1597. END