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

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