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

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