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.

dgesvj.f 68 kB

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