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

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