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.

sgejsv.f 72 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779
  1. *> \brief \b SGEJSV
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download SGEJSV + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgejsv.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgejsv.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgejsv.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE SGEJSV( JOBA, JOBU, JOBV, JOBR, JOBT, JOBP,
  22. * M, N, A, LDA, SVA, U, LDU, V, LDV,
  23. * WORK, LWORK, IWORK, INFO )
  24. *
  25. * .. Scalar Arguments ..
  26. * IMPLICIT NONE
  27. * INTEGER INFO, LDA, LDU, LDV, LWORK, M, N
  28. * ..
  29. * .. Array Arguments ..
  30. * REAL A( LDA, * ), SVA( N ), U( LDU, * ), V( LDV, * ),
  31. * $ WORK( LWORK )
  32. * INTEGER IWORK( * )
  33. * CHARACTER*1 JOBA, JOBP, JOBR, JOBT, JOBU, JOBV
  34. * ..
  35. *
  36. *
  37. *> \par Purpose:
  38. * =============
  39. *>
  40. *> \verbatim
  41. *>
  42. *> SGEJSV computes the singular value decomposition (SVD) of a real M-by-N
  43. *> matrix [A], where M >= N. The SVD of [A] is written as
  44. *>
  45. *> [A] = [U] * [SIGMA] * [V]^t,
  46. *>
  47. *> where [SIGMA] is an N-by-N (M-by-N) matrix which is zero except for its N
  48. *> diagonal elements, [U] is an M-by-N (or M-by-M) orthonormal matrix, and
  49. *> [V] is an N-by-N orthogonal matrix. The diagonal elements of [SIGMA] are
  50. *> the singular values of [A]. The columns of [U] and [V] are the left and
  51. *> the right singular vectors of [A], respectively. The matrices [U] and [V]
  52. *> are computed and stored in the arrays U and V, respectively. The diagonal
  53. *> of [SIGMA] is computed and stored in the array SVA.
  54. *> SGEJSV can sometimes compute tiny singular values and their singular vectors much
  55. *> more accurately than other SVD routines, see below under Further Details.
  56. *> \endverbatim
  57. *
  58. * Arguments:
  59. * ==========
  60. *
  61. *> \param[in] JOBA
  62. *> \verbatim
  63. *> JOBA is CHARACTER*1
  64. *> Specifies the level of accuracy:
  65. *> = 'C': This option works well (high relative accuracy) if A = B * D,
  66. *> with well-conditioned B and arbitrary diagonal matrix D.
  67. *> The accuracy cannot be spoiled by COLUMN scaling. The
  68. *> accuracy of the computed output depends on the condition of
  69. *> B, and the procedure aims at the best theoretical accuracy.
  70. *> The relative error max_{i=1:N}|d sigma_i| / sigma_i is
  71. *> bounded by f(M,N)*epsilon* cond(B), independent of D.
  72. *> The input matrix is preprocessed with the QRF with column
  73. *> pivoting. This initial preprocessing and preconditioning by
  74. *> a rank revealing QR factorization is common for all values of
  75. *> JOBA. Additional actions are specified as follows:
  76. *> = 'E': Computation as with 'C' with an additional estimate of the
  77. *> condition number of B. It provides a realistic error bound.
  78. *> = 'F': If A = D1 * C * D2 with ill-conditioned diagonal scalings
  79. *> D1, D2, and well-conditioned matrix C, this option gives
  80. *> higher accuracy than the 'C' option. If the structure of the
  81. *> input matrix is not known, and relative accuracy is
  82. *> desirable, then this option is advisable. The input matrix A
  83. *> is preprocessed with QR factorization with FULL (row and
  84. *> column) pivoting.
  85. *> = 'G' Computation as with 'F' with an additional estimate of the
  86. *> condition number of B, where A=D*B. If A has heavily weighted
  87. *> rows, then using this condition number gives too pessimistic
  88. *> error bound.
  89. *> = 'A': Small singular values are the noise and the matrix is treated
  90. *> as numerically rank defficient. The error in the computed
  91. *> singular values is bounded by f(m,n)*epsilon*||A||.
  92. *> The computed SVD A = U * S * V^t restores A up to
  93. *> f(m,n)*epsilon*||A||.
  94. *> This gives the procedure the licence to discard (set to zero)
  95. *> all singular values below N*epsilon*||A||.
  96. *> = 'R': Similar as in 'A'. Rank revealing property of the initial
  97. *> QR factorization is used do reveal (using triangular factor)
  98. *> a gap sigma_{r+1} < epsilon * sigma_r in which case the
  99. *> numerical RANK is declared to be r. The SVD is computed with
  100. *> absolute error bounds, but more accurately than with 'A'.
  101. *> \endverbatim
  102. *>
  103. *> \param[in] JOBU
  104. *> \verbatim
  105. *> JOBU is CHARACTER*1
  106. *> Specifies whether to compute the columns of U:
  107. *> = 'U': N columns of U are returned in the array U.
  108. *> = 'F': full set of M left sing. vectors is returned in the array U.
  109. *> = 'W': U may be used as workspace of length M*N. See the description
  110. *> of U.
  111. *> = 'N': U is not computed.
  112. *> \endverbatim
  113. *>
  114. *> \param[in] JOBV
  115. *> \verbatim
  116. *> JOBV is CHARACTER*1
  117. *> Specifies whether to compute the matrix V:
  118. *> = 'V': N columns of V are returned in the array V; Jacobi rotations
  119. *> are not explicitly accumulated.
  120. *> = 'J': N columns of V are returned in the array V, but they are
  121. *> computed as the product of Jacobi rotations. This option is
  122. *> allowed only if JOBU .NE. 'N', i.e. in computing the full SVD.
  123. *> = 'W': V may be used as workspace of length N*N. See the description
  124. *> of V.
  125. *> = 'N': V is not computed.
  126. *> \endverbatim
  127. *>
  128. *> \param[in] JOBR
  129. *> \verbatim
  130. *> JOBR is CHARACTER*1
  131. *> Specifies the RANGE for the singular values. Issues the licence to
  132. *> set to zero small positive singular values if they are outside
  133. *> specified range. If A .NE. 0 is scaled so that the largest singular
  134. *> value of c*A is around SQRT(BIG), BIG=SLAMCH('O'), then JOBR issues
  135. *> the licence to kill columns of A whose norm in c*A is less than
  136. *> SQRT(SFMIN) (for JOBR.EQ.'R'), or less than SMALL=SFMIN/EPSLN,
  137. *> where SFMIN=SLAMCH('S'), EPSLN=SLAMCH('E').
  138. *> = 'N': Do not kill small columns of c*A. This option assumes that
  139. *> BLAS and QR factorizations and triangular solvers are
  140. *> implemented to work in that range. If the condition of A
  141. *> is greater than BIG, use SGESVJ.
  142. *> = 'R': RESTRICTED range for sigma(c*A) is [SQRT(SFMIN), SQRT(BIG)]
  143. *> (roughly, as described above). This option is recommended.
  144. *> ===========================
  145. *> For computing the singular values in the FULL range [SFMIN,BIG]
  146. *> use SGESVJ.
  147. *> \endverbatim
  148. *>
  149. *> \param[in] JOBT
  150. *> \verbatim
  151. *> JOBT is CHARACTER*1
  152. *> If the matrix is square then the procedure may determine to use
  153. *> transposed A if A^t seems to be better with respect to convergence.
  154. *> If the matrix is not square, JOBT is ignored. This is subject to
  155. *> changes in the future.
  156. *> The decision is based on two values of entropy over the adjoint
  157. *> orbit of A^t * A. See the descriptions of WORK(6) and WORK(7).
  158. *> = 'T': transpose if entropy test indicates possibly faster
  159. *> convergence of Jacobi process if A^t is taken as input. If A is
  160. *> replaced with A^t, then the row pivoting is included automatically.
  161. *> = 'N': do not speculate.
  162. *> This option can be used to compute only the singular values, or the
  163. *> full SVD (U, SIGMA and V). For only one set of singular vectors
  164. *> (U or V), the caller should provide both U and V, as one of the
  165. *> matrices is used as workspace if the matrix A is transposed.
  166. *> The implementer can easily remove this constraint and make the
  167. *> code more complicated. See the descriptions of U and V.
  168. *> \endverbatim
  169. *>
  170. *> \param[in] JOBP
  171. *> \verbatim
  172. *> JOBP is CHARACTER*1
  173. *> Issues the licence to introduce structured perturbations to drown
  174. *> denormalized numbers. This licence should be active if the
  175. *> denormals are poorly implemented, causing slow computation,
  176. *> especially in cases of fast convergence (!). For details see [1,2].
  177. *> For the sake of simplicity, this perturbations are included only
  178. *> when the full SVD or only the singular values are requested. The
  179. *> implementer/user can easily add the perturbation for the cases of
  180. *> computing one set of singular vectors.
  181. *> = 'P': introduce perturbation
  182. *> = 'N': do not perturb
  183. *> \endverbatim
  184. *>
  185. *> \param[in] M
  186. *> \verbatim
  187. *> M is INTEGER
  188. *> The number of rows of the input matrix A. M >= 0.
  189. *> \endverbatim
  190. *>
  191. *> \param[in] N
  192. *> \verbatim
  193. *> N is INTEGER
  194. *> The number of columns of the input matrix A. M >= N >= 0.
  195. *> \endverbatim
  196. *>
  197. *> \param[in,out] A
  198. *> \verbatim
  199. *> A is REAL array, dimension (LDA,N)
  200. *> On entry, the M-by-N matrix A.
  201. *> \endverbatim
  202. *>
  203. *> \param[in] LDA
  204. *> \verbatim
  205. *> LDA is INTEGER
  206. *> The leading dimension of the array A. LDA >= max(1,M).
  207. *> \endverbatim
  208. *>
  209. *> \param[out] SVA
  210. *> \verbatim
  211. *> SVA is REAL array, dimension (N)
  212. *> On exit,
  213. *> - For WORK(1)/WORK(2) = ONE: The singular values of A. During the
  214. *> computation SVA contains Euclidean column norms of the
  215. *> iterated matrices in the array A.
  216. *> - For WORK(1) .NE. WORK(2): The singular values of A are
  217. *> (WORK(1)/WORK(2)) * SVA(1:N). This factored form is used if
  218. *> sigma_max(A) overflows or if small singular values have been
  219. *> saved from underflow by scaling the input matrix A.
  220. *> - If JOBR='R' then some of the singular values may be returned
  221. *> as exact zeros obtained by "set to zero" because they are
  222. *> below the numerical rank threshold or are denormalized numbers.
  223. *> \endverbatim
  224. *>
  225. *> \param[out] U
  226. *> \verbatim
  227. *> U is REAL array, dimension ( LDU, N )
  228. *> If JOBU = 'U', then U contains on exit the M-by-N matrix of
  229. *> the left singular vectors.
  230. *> If JOBU = 'F', then U contains on exit the M-by-M matrix of
  231. *> the left singular vectors, including an ONB
  232. *> of the orthogonal complement of the Range(A).
  233. *> If JOBU = 'W' .AND. (JOBV.EQ.'V' .AND. JOBT.EQ.'T' .AND. M.EQ.N),
  234. *> then U is used as workspace if the procedure
  235. *> replaces A with A^t. In that case, [V] is computed
  236. *> in U as left singular vectors of A^t and then
  237. *> copied back to the V array. This 'W' option is just
  238. *> a reminder to the caller that in this case U is
  239. *> reserved as workspace of length N*N.
  240. *> If JOBU = 'N' U is not referenced.
  241. *> \endverbatim
  242. *>
  243. *> \param[in] LDU
  244. *> \verbatim
  245. *> LDU is INTEGER
  246. *> The leading dimension of the array U, LDU >= 1.
  247. *> IF JOBU = 'U' or 'F' or 'W', then LDU >= M.
  248. *> \endverbatim
  249. *>
  250. *> \param[out] V
  251. *> \verbatim
  252. *> V is REAL array, dimension ( LDV, N )
  253. *> If JOBV = 'V', 'J' then V contains on exit the N-by-N matrix of
  254. *> the right singular vectors;
  255. *> If JOBV = 'W', AND (JOBU.EQ.'U' AND JOBT.EQ.'T' AND M.EQ.N),
  256. *> then V is used as workspace if the pprocedure
  257. *> replaces A with A^t. In that case, [U] is computed
  258. *> in V as right singular vectors of A^t and then
  259. *> copied back to the U array. This 'W' option is just
  260. *> a reminder to the caller that in this case V is
  261. *> reserved as workspace of length N*N.
  262. *> If JOBV = 'N' V is not referenced.
  263. *> \endverbatim
  264. *>
  265. *> \param[in] LDV
  266. *> \verbatim
  267. *> LDV is INTEGER
  268. *> The leading dimension of the array V, LDV >= 1.
  269. *> If JOBV = 'V' or 'J' or 'W', then LDV >= N.
  270. *> \endverbatim
  271. *>
  272. *> \param[out] WORK
  273. *> \verbatim
  274. *> WORK is REAL array, dimension at least LWORK.
  275. *> On exit,
  276. *> WORK(1) = SCALE = WORK(2) / WORK(1) is the scaling factor such
  277. *> that SCALE*SVA(1:N) are the computed singular values
  278. *> of A. (See the description of SVA().)
  279. *> WORK(2) = See the description of WORK(1).
  280. *> WORK(3) = SCONDA is an estimate for the condition number of
  281. *> column equilibrated A. (If JOBA .EQ. 'E' or 'G')
  282. *> SCONDA is an estimate of SQRT(||(R^t * R)^(-1)||_1).
  283. *> It is computed using SPOCON. It holds
  284. *> N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
  285. *> where R is the triangular factor from the QRF of A.
  286. *> However, if R is truncated and the numerical rank is
  287. *> determined to be strictly smaller than N, SCONDA is
  288. *> returned as -1, thus indicating that the smallest
  289. *> singular values might be lost.
  290. *>
  291. *> If full SVD is needed, the following two condition numbers are
  292. *> useful for the analysis of the algorithm. They are provied for
  293. *> a developer/implementer who is familiar with the details of
  294. *> the method.
  295. *>
  296. *> WORK(4) = an estimate of the scaled condition number of the
  297. *> triangular factor in the first QR factorization.
  298. *> WORK(5) = an estimate of the scaled condition number of the
  299. *> triangular factor in the second QR factorization.
  300. *> The following two parameters are computed if JOBT .EQ. 'T'.
  301. *> They are provided for a developer/implementer who is familiar
  302. *> with the details of the method.
  303. *>
  304. *> WORK(6) = the entropy of A^t*A :: this is the Shannon entropy
  305. *> of diag(A^t*A) / Trace(A^t*A) taken as point in the
  306. *> probability simplex.
  307. *> WORK(7) = the entropy of A*A^t.
  308. *> \endverbatim
  309. *>
  310. *> \param[in] LWORK
  311. *> \verbatim
  312. *> LWORK is INTEGER
  313. *> Length of WORK to confirm proper allocation of work space.
  314. *> LWORK depends on the job:
  315. *>
  316. *> If only SIGMA is needed ( JOBU.EQ.'N', JOBV.EQ.'N' ) and
  317. *> -> .. no scaled condition estimate required (JOBE.EQ.'N'):
  318. *> LWORK >= max(2*M+N,4*N+1,7). This is the minimal requirement.
  319. *> ->> For optimal performance (blocked code) the optimal value
  320. *> is LWORK >= max(2*M+N,3*N+(N+1)*NB,7). Here NB is the optimal
  321. *> block size for DGEQP3 and DGEQRF.
  322. *> In general, optimal LWORK is computed as
  323. *> LWORK >= max(2*M+N,N+LWORK(DGEQP3),N+LWORK(DGEQRF), 7).
  324. *> -> .. an estimate of the scaled condition number of A is
  325. *> required (JOBA='E', 'G'). In this case, LWORK is the maximum
  326. *> of the above and N*N+4*N, i.e. LWORK >= max(2*M+N,N*N+4*N,7).
  327. *> ->> For optimal performance (blocked code) the optimal value
  328. *> is LWORK >= max(2*M+N,3*N+(N+1)*NB, N*N+4*N, 7).
  329. *> In general, the optimal length LWORK is computed as
  330. *> LWORK >= max(2*M+N,N+LWORK(DGEQP3),N+LWORK(DGEQRF),
  331. *> N+N*N+LWORK(DPOCON),7).
  332. *>
  333. *> If SIGMA and the right singular vectors are needed (JOBV.EQ.'V'),
  334. *> -> the minimal requirement is LWORK >= max(2*M+N,4*N+1,7).
  335. *> -> For optimal performance, LWORK >= max(2*M+N,3*N+(N+1)*NB,7),
  336. *> where NB is the optimal block size for DGEQP3, DGEQRF, DGELQ,
  337. *> DORMLQ. In general, the optimal length LWORK is computed as
  338. *> LWORK >= max(2*M+N,N+LWORK(DGEQP3), N+LWORK(DPOCON),
  339. *> N+LWORK(DGELQ), 2*N+LWORK(DGEQRF), N+LWORK(DORMLQ)).
  340. *>
  341. *> If SIGMA and the left singular vectors are needed
  342. *> -> the minimal requirement is LWORK >= max(2*M+N,4*N+1,7).
  343. *> -> For optimal performance:
  344. *> if JOBU.EQ.'U' :: LWORK >= max(2*M+N,3*N+(N+1)*NB,7),
  345. *> if JOBU.EQ.'F' :: LWORK >= max(2*M+N,3*N+(N+1)*NB,N+M*NB,7),
  346. *> where NB is the optimal block size for DGEQP3, DGEQRF, DORMQR.
  347. *> In general, the optimal length LWORK is computed as
  348. *> LWORK >= max(2*M+N,N+LWORK(DGEQP3),N+LWORK(DPOCON),
  349. *> 2*N+LWORK(DGEQRF), N+LWORK(DORMQR)).
  350. *> Here LWORK(DORMQR) equals N*NB (for JOBU.EQ.'U') or
  351. *> M*NB (for JOBU.EQ.'F').
  352. *>
  353. *> If the full SVD is needed: (JOBU.EQ.'U' or JOBU.EQ.'F') and
  354. *> -> if JOBV.EQ.'V'
  355. *> the minimal requirement is LWORK >= max(2*M+N,6*N+2*N*N).
  356. *> -> if JOBV.EQ.'J' the minimal requirement is
  357. *> LWORK >= max(2*M+N, 4*N+N*N,2*N+N*N+6).
  358. *> -> For optimal performance, LWORK should be additionally
  359. *> larger than N+M*NB, where NB is the optimal block size
  360. *> for DORMQR.
  361. *> \endverbatim
  362. *>
  363. *> \param[out] IWORK
  364. *> \verbatim
  365. *> IWORK is INTEGER array, dimension M+3*N.
  366. *> On exit,
  367. *> IWORK(1) = the numerical rank determined after the initial
  368. *> QR factorization with pivoting. See the descriptions
  369. *> of JOBA and JOBR.
  370. *> IWORK(2) = the number of the computed nonzero singular values
  371. *> IWORK(3) = if nonzero, a warning message:
  372. *> If IWORK(3).EQ.1 then some of the column norms of A
  373. *> were denormalized floats. The requested high accuracy
  374. *> is not warranted by the data.
  375. *> \endverbatim
  376. *>
  377. *> \param[out] INFO
  378. *> \verbatim
  379. *> INFO is INTEGER
  380. *> < 0 : if INFO = -i, then the i-th argument had an illegal value.
  381. *> = 0 : successfull exit;
  382. *> > 0 : SGEJSV did not converge in the maximal allowed number
  383. *> of sweeps. The computed values may be inaccurate.
  384. *> \endverbatim
  385. *
  386. * Authors:
  387. * ========
  388. *
  389. *> \author Univ. of Tennessee
  390. *> \author Univ. of California Berkeley
  391. *> \author Univ. of Colorado Denver
  392. *> \author NAG Ltd.
  393. *
  394. *> \date November 2015
  395. *
  396. *> \ingroup realGEsing
  397. *
  398. *> \par Further Details:
  399. * =====================
  400. *>
  401. *> \verbatim
  402. *>
  403. *> SGEJSV implements a preconditioned Jacobi SVD algorithm. It uses SGEQP3,
  404. *> SGEQRF, and SGELQF as preprocessors and preconditioners. Optionally, an
  405. *> additional row pivoting can be used as a preprocessor, which in some
  406. *> cases results in much higher accuracy. An example is matrix A with the
  407. *> structure A = D1 * C * D2, where D1, D2 are arbitrarily ill-conditioned
  408. *> diagonal matrices and C is well-conditioned matrix. In that case, complete
  409. *> pivoting in the first QR factorizations provides accuracy dependent on the
  410. *> condition number of C, and independent of D1, D2. Such higher accuracy is
  411. *> not completely understood theoretically, but it works well in practice.
  412. *> Further, if A can be written as A = B*D, with well-conditioned B and some
  413. *> diagonal D, then the high accuracy is guaranteed, both theoretically and
  414. *> in software, independent of D. For more details see [1], [2].
  415. *> The computational range for the singular values can be the full range
  416. *> ( UNDERFLOW,OVERFLOW ), provided that the machine arithmetic and the BLAS
  417. *> & LAPACK routines called by SGEJSV are implemented to work in that range.
  418. *> If that is not the case, then the restriction for safe computation with
  419. *> the singular values in the range of normalized IEEE numbers is that the
  420. *> spectral condition number kappa(A)=sigma_max(A)/sigma_min(A) does not
  421. *> overflow. This code (SGEJSV) is best used in this restricted range,
  422. *> meaning that singular values of magnitude below ||A||_2 / SLAMCH('O') are
  423. *> returned as zeros. See JOBR for details on this.
  424. *> Further, this implementation is somewhat slower than the one described
  425. *> in [1,2] due to replacement of some non-LAPACK components, and because
  426. *> the choice of some tuning parameters in the iterative part (SGESVJ) is
  427. *> left to the implementer on a particular machine.
  428. *> The rank revealing QR factorization (in this code: SGEQP3) should be
  429. *> implemented as in [3]. We have a new version of SGEQP3 under development
  430. *> that is more robust than the current one in LAPACK, with a cleaner cut in
  431. *> rank defficient cases. It will be available in the SIGMA library [4].
  432. *> If M is much larger than N, it is obvious that the inital QRF with
  433. *> column pivoting can be preprocessed by the QRF without pivoting. That
  434. *> well known trick is not used in SGEJSV because in some cases heavy row
  435. *> weighting can be treated with complete pivoting. The overhead in cases
  436. *> M much larger than N is then only due to pivoting, but the benefits in
  437. *> terms of accuracy have prevailed. The implementer/user can incorporate
  438. *> this extra QRF step easily. The implementer can also improve data movement
  439. *> (matrix transpose, matrix copy, matrix transposed copy) - this
  440. *> implementation of SGEJSV uses only the simplest, naive data movement.
  441. *> \endverbatim
  442. *
  443. *> \par Contributors:
  444. * ==================
  445. *>
  446. *> Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany)
  447. *
  448. *> \par References:
  449. * ================
  450. *>
  451. *> \verbatim
  452. *>
  453. *> [1] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm I.
  454. *> SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1322-1342.
  455. *> LAPACK Working note 169.
  456. *> [2] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm II.
  457. *> SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1343-1362.
  458. *> LAPACK Working note 170.
  459. *> [3] Z. Drmac and Z. Bujanovic: On the failure of rank-revealing QR
  460. *> factorization software - a case study.
  461. *> ACM Trans. math. Softw. Vol. 35, No 2 (2008), pp. 1-28.
  462. *> LAPACK Working note 176.
  463. *> [4] Z. Drmac: SIGMA - mathematical software library for accurate SVD, PSV,
  464. *> QSVD, (H,K)-SVD computations.
  465. *> Department of Mathematics, University of Zagreb, 2008.
  466. *> \endverbatim
  467. *
  468. *> \par Bugs, examples and comments:
  469. * =================================
  470. *>
  471. *> Please report all bugs and send interesting examples and/or comments to
  472. *> drmac@math.hr. Thank you.
  473. *>
  474. * =====================================================================
  475. SUBROUTINE SGEJSV( JOBA, JOBU, JOBV, JOBR, JOBT, JOBP,
  476. $ M, N, A, LDA, SVA, U, LDU, V, LDV,
  477. $ WORK, LWORK, IWORK, INFO )
  478. *
  479. * -- LAPACK computational routine (version 3.6.0) --
  480. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  481. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  482. * November 2015
  483. *
  484. * .. Scalar Arguments ..
  485. IMPLICIT NONE
  486. INTEGER INFO, LDA, LDU, LDV, LWORK, M, N
  487. * ..
  488. * .. Array Arguments ..
  489. REAL A( LDA, * ), SVA( N ), U( LDU, * ), V( LDV, * ),
  490. $ WORK( LWORK )
  491. INTEGER IWORK( * )
  492. CHARACTER*1 JOBA, JOBP, JOBR, JOBT, JOBU, JOBV
  493. * ..
  494. *
  495. * ===========================================================================
  496. *
  497. * .. Local Parameters ..
  498. REAL ZERO, ONE
  499. PARAMETER ( ZERO = 0.0E0, ONE = 1.0E0 )
  500. * ..
  501. * .. Local Scalars ..
  502. REAL AAPP, AAQQ, AATMAX, AATMIN, BIG, BIG1, COND_OK,
  503. $ CONDR1, CONDR2, ENTRA, ENTRAT, EPSLN, MAXPRJ, SCALEM,
  504. $ SCONDA, SFMIN, SMALL, TEMP1, USCAL1, USCAL2, XSC
  505. INTEGER IERR, N1, NR, NUMRANK, p, q, WARNING
  506. LOGICAL ALMORT, DEFR, ERREST, GOSCAL, JRACC, KILL, LSVEC,
  507. $ L2ABER, L2KILL, L2PERT, L2RANK, L2TRAN,
  508. $ NOSCAL, ROWPIV, RSVEC, TRANSP
  509. * ..
  510. * .. Intrinsic Functions ..
  511. INTRINSIC ABS, ALOG, MAX, MIN, FLOAT, NINT, SIGN, SQRT
  512. * ..
  513. * .. External Functions ..
  514. REAL SLAMCH, SNRM2
  515. INTEGER ISAMAX
  516. LOGICAL LSAME
  517. EXTERNAL ISAMAX, LSAME, SLAMCH, SNRM2
  518. * ..
  519. * .. External Subroutines ..
  520. EXTERNAL SCOPY, SGELQF, SGEQP3, SGEQRF, SLACPY, SLASCL,
  521. $ SLASET, SLASSQ, SLASWP, SORGQR, SORMLQ,
  522. $ SORMQR, SPOCON, SSCAL, SSWAP, STRSM, XERBLA
  523. *
  524. EXTERNAL SGESVJ
  525. * ..
  526. *
  527. * Test the input arguments
  528. *
  529. LSVEC = LSAME( JOBU, 'U' ) .OR. LSAME( JOBU, 'F' )
  530. JRACC = LSAME( JOBV, 'J' )
  531. RSVEC = LSAME( JOBV, 'V' ) .OR. JRACC
  532. ROWPIV = LSAME( JOBA, 'F' ) .OR. LSAME( JOBA, 'G' )
  533. L2RANK = LSAME( JOBA, 'R' )
  534. L2ABER = LSAME( JOBA, 'A' )
  535. ERREST = LSAME( JOBA, 'E' ) .OR. LSAME( JOBA, 'G' )
  536. L2TRAN = LSAME( JOBT, 'T' )
  537. L2KILL = LSAME( JOBR, 'R' )
  538. DEFR = LSAME( JOBR, 'N' )
  539. L2PERT = LSAME( JOBP, 'P' )
  540. *
  541. IF ( .NOT.(ROWPIV .OR. L2RANK .OR. L2ABER .OR.
  542. $ ERREST .OR. LSAME( JOBA, 'C' ) )) THEN
  543. INFO = - 1
  544. ELSE IF ( .NOT.( LSVEC .OR. LSAME( JOBU, 'N' ) .OR.
  545. $ LSAME( JOBU, 'W' )) ) THEN
  546. INFO = - 2
  547. ELSE IF ( .NOT.( RSVEC .OR. LSAME( JOBV, 'N' ) .OR.
  548. $ LSAME( JOBV, 'W' )) .OR. ( JRACC .AND. (.NOT.LSVEC) ) ) THEN
  549. INFO = - 3
  550. ELSE IF ( .NOT. ( L2KILL .OR. DEFR ) ) THEN
  551. INFO = - 4
  552. ELSE IF ( .NOT. ( L2TRAN .OR. LSAME( JOBT, 'N' ) ) ) THEN
  553. INFO = - 5
  554. ELSE IF ( .NOT. ( L2PERT .OR. LSAME( JOBP, 'N' ) ) ) THEN
  555. INFO = - 6
  556. ELSE IF ( M .LT. 0 ) THEN
  557. INFO = - 7
  558. ELSE IF ( ( N .LT. 0 ) .OR. ( N .GT. M ) ) THEN
  559. INFO = - 8
  560. ELSE IF ( LDA .LT. M ) THEN
  561. INFO = - 10
  562. ELSE IF ( LSVEC .AND. ( LDU .LT. M ) ) THEN
  563. INFO = - 13
  564. ELSE IF ( RSVEC .AND. ( LDV .LT. N ) ) THEN
  565. INFO = - 14
  566. ELSE IF ( (.NOT.(LSVEC .OR. RSVEC .OR. ERREST).AND.
  567. $ (LWORK .LT. MAX(7,4*N+1,2*M+N))) .OR.
  568. $ (.NOT.(LSVEC .OR. RSVEC) .AND. ERREST .AND.
  569. $ (LWORK .LT. MAX(7,4*N+N*N,2*M+N))) .OR.
  570. $ (LSVEC .AND. (.NOT.RSVEC) .AND. (LWORK .LT. MAX(7,2*M+N,4*N+1)))
  571. $ .OR.
  572. $ (RSVEC .AND. (.NOT.LSVEC) .AND. (LWORK .LT. MAX(7,2*M+N,4*N+1)))
  573. $ .OR.
  574. $ (LSVEC .AND. RSVEC .AND. (.NOT.JRACC) .AND.
  575. $ (LWORK.LT.MAX(2*M+N,6*N+2*N*N)))
  576. $ .OR. (LSVEC .AND. RSVEC .AND. JRACC .AND.
  577. $ LWORK.LT.MAX(2*M+N,4*N+N*N,2*N+N*N+6)))
  578. $ THEN
  579. INFO = - 17
  580. ELSE
  581. * #:)
  582. INFO = 0
  583. END IF
  584. *
  585. IF ( INFO .NE. 0 ) THEN
  586. * #:(
  587. CALL XERBLA( 'SGEJSV', - INFO )
  588. RETURN
  589. END IF
  590. *
  591. * Quick return for void matrix (Y3K safe)
  592. * #:)
  593. IF ( ( M .EQ. 0 ) .OR. ( N .EQ. 0 ) ) RETURN
  594. *
  595. * Determine whether the matrix U should be M x N or M x M
  596. *
  597. IF ( LSVEC ) THEN
  598. N1 = N
  599. IF ( LSAME( JOBU, 'F' ) ) N1 = M
  600. END IF
  601. *
  602. * Set numerical parameters
  603. *
  604. *! NOTE: Make sure SLAMCH() does not fail on the target architecture.
  605. *
  606. EPSLN = SLAMCH('Epsilon')
  607. SFMIN = SLAMCH('SafeMinimum')
  608. SMALL = SFMIN / EPSLN
  609. BIG = SLAMCH('O')
  610. * BIG = ONE / SFMIN
  611. *
  612. * Initialize SVA(1:N) = diag( ||A e_i||_2 )_1^N
  613. *
  614. *(!) If necessary, scale SVA() to protect the largest norm from
  615. * overflow. It is possible that this scaling pushes the smallest
  616. * column norm left from the underflow threshold (extreme case).
  617. *
  618. SCALEM = ONE / SQRT(FLOAT(M)*FLOAT(N))
  619. NOSCAL = .TRUE.
  620. GOSCAL = .TRUE.
  621. DO 1874 p = 1, N
  622. AAPP = ZERO
  623. AAQQ = ONE
  624. CALL SLASSQ( M, A(1,p), 1, AAPP, AAQQ )
  625. IF ( AAPP .GT. BIG ) THEN
  626. INFO = - 9
  627. CALL XERBLA( 'SGEJSV', -INFO )
  628. RETURN
  629. END IF
  630. AAQQ = SQRT(AAQQ)
  631. IF ( ( AAPP .LT. (BIG / AAQQ) ) .AND. NOSCAL ) THEN
  632. SVA(p) = AAPP * AAQQ
  633. ELSE
  634. NOSCAL = .FALSE.
  635. SVA(p) = AAPP * ( AAQQ * SCALEM )
  636. IF ( GOSCAL ) THEN
  637. GOSCAL = .FALSE.
  638. CALL SSCAL( p-1, SCALEM, SVA, 1 )
  639. END IF
  640. END IF
  641. 1874 CONTINUE
  642. *
  643. IF ( NOSCAL ) SCALEM = ONE
  644. *
  645. AAPP = ZERO
  646. AAQQ = BIG
  647. DO 4781 p = 1, N
  648. AAPP = MAX( AAPP, SVA(p) )
  649. IF ( SVA(p) .NE. ZERO ) AAQQ = MIN( AAQQ, SVA(p) )
  650. 4781 CONTINUE
  651. *
  652. * Quick return for zero M x N matrix
  653. * #:)
  654. IF ( AAPP .EQ. ZERO ) THEN
  655. IF ( LSVEC ) CALL SLASET( 'G', M, N1, ZERO, ONE, U, LDU )
  656. IF ( RSVEC ) CALL SLASET( 'G', N, N, ZERO, ONE, V, LDV )
  657. WORK(1) = ONE
  658. WORK(2) = ONE
  659. IF ( ERREST ) WORK(3) = ONE
  660. IF ( LSVEC .AND. RSVEC ) THEN
  661. WORK(4) = ONE
  662. WORK(5) = ONE
  663. END IF
  664. IF ( L2TRAN ) THEN
  665. WORK(6) = ZERO
  666. WORK(7) = ZERO
  667. END IF
  668. IWORK(1) = 0
  669. IWORK(2) = 0
  670. IWORK(3) = 0
  671. RETURN
  672. END IF
  673. *
  674. * Issue warning if denormalized column norms detected. Override the
  675. * high relative accuracy request. Issue licence to kill columns
  676. * (set them to zero) whose norm is less than sigma_max / BIG (roughly).
  677. * #:(
  678. WARNING = 0
  679. IF ( AAQQ .LE. SFMIN ) THEN
  680. L2RANK = .TRUE.
  681. L2KILL = .TRUE.
  682. WARNING = 1
  683. END IF
  684. *
  685. * Quick return for one-column matrix
  686. * #:)
  687. IF ( N .EQ. 1 ) THEN
  688. *
  689. IF ( LSVEC ) THEN
  690. CALL SLASCL( 'G',0,0,SVA(1),SCALEM, M,1,A(1,1),LDA,IERR )
  691. CALL SLACPY( 'A', M, 1, A, LDA, U, LDU )
  692. * computing all M left singular vectors of the M x 1 matrix
  693. IF ( N1 .NE. N ) THEN
  694. CALL SGEQRF( M, N, U,LDU, WORK, WORK(N+1),LWORK-N,IERR )
  695. CALL SORGQR( M,N1,1, U,LDU,WORK,WORK(N+1),LWORK-N,IERR )
  696. CALL SCOPY( M, A(1,1), 1, U(1,1), 1 )
  697. END IF
  698. END IF
  699. IF ( RSVEC ) THEN
  700. V(1,1) = ONE
  701. END IF
  702. IF ( SVA(1) .LT. (BIG*SCALEM) ) THEN
  703. SVA(1) = SVA(1) / SCALEM
  704. SCALEM = ONE
  705. END IF
  706. WORK(1) = ONE / SCALEM
  707. WORK(2) = ONE
  708. IF ( SVA(1) .NE. ZERO ) THEN
  709. IWORK(1) = 1
  710. IF ( ( SVA(1) / SCALEM) .GE. SFMIN ) THEN
  711. IWORK(2) = 1
  712. ELSE
  713. IWORK(2) = 0
  714. END IF
  715. ELSE
  716. IWORK(1) = 0
  717. IWORK(2) = 0
  718. END IF
  719. IF ( ERREST ) WORK(3) = ONE
  720. IF ( LSVEC .AND. RSVEC ) THEN
  721. WORK(4) = ONE
  722. WORK(5) = ONE
  723. END IF
  724. IF ( L2TRAN ) THEN
  725. WORK(6) = ZERO
  726. WORK(7) = ZERO
  727. END IF
  728. RETURN
  729. *
  730. END IF
  731. *
  732. TRANSP = .FALSE.
  733. L2TRAN = L2TRAN .AND. ( M .EQ. N )
  734. *
  735. AATMAX = -ONE
  736. AATMIN = BIG
  737. IF ( ROWPIV .OR. L2TRAN ) THEN
  738. *
  739. * Compute the row norms, needed to determine row pivoting sequence
  740. * (in the case of heavily row weighted A, row pivoting is strongly
  741. * advised) and to collect information needed to compare the
  742. * structures of A * A^t and A^t * A (in the case L2TRAN.EQ..TRUE.).
  743. *
  744. IF ( L2TRAN ) THEN
  745. DO 1950 p = 1, M
  746. XSC = ZERO
  747. TEMP1 = ONE
  748. CALL SLASSQ( N, A(p,1), LDA, XSC, TEMP1 )
  749. * SLASSQ gets both the ell_2 and the ell_infinity norm
  750. * in one pass through the vector
  751. WORK(M+N+p) = XSC * SCALEM
  752. WORK(N+p) = XSC * (SCALEM*SQRT(TEMP1))
  753. AATMAX = MAX( AATMAX, WORK(N+p) )
  754. IF (WORK(N+p) .NE. ZERO) AATMIN = MIN(AATMIN,WORK(N+p))
  755. 1950 CONTINUE
  756. ELSE
  757. DO 1904 p = 1, M
  758. WORK(M+N+p) = SCALEM*ABS( A(p,ISAMAX(N,A(p,1),LDA)) )
  759. AATMAX = MAX( AATMAX, WORK(M+N+p) )
  760. AATMIN = MIN( AATMIN, WORK(M+N+p) )
  761. 1904 CONTINUE
  762. END IF
  763. *
  764. END IF
  765. *
  766. * For square matrix A try to determine whether A^t would be better
  767. * input for the preconditioned Jacobi SVD, with faster convergence.
  768. * The decision is based on an O(N) function of the vector of column
  769. * and row norms of A, based on the Shannon entropy. This should give
  770. * the right choice in most cases when the difference actually matters.
  771. * It may fail and pick the slower converging side.
  772. *
  773. ENTRA = ZERO
  774. ENTRAT = ZERO
  775. IF ( L2TRAN ) THEN
  776. *
  777. XSC = ZERO
  778. TEMP1 = ONE
  779. CALL SLASSQ( N, SVA, 1, XSC, TEMP1 )
  780. TEMP1 = ONE / TEMP1
  781. *
  782. ENTRA = ZERO
  783. DO 1113 p = 1, N
  784. BIG1 = ( ( SVA(p) / XSC )**2 ) * TEMP1
  785. IF ( BIG1 .NE. ZERO ) ENTRA = ENTRA + BIG1 * ALOG(BIG1)
  786. 1113 CONTINUE
  787. ENTRA = - ENTRA / ALOG(FLOAT(N))
  788. *
  789. * Now, SVA().^2/Trace(A^t * A) is a point in the probability simplex.
  790. * It is derived from the diagonal of A^t * A. Do the same with the
  791. * diagonal of A * A^t, compute the entropy of the corresponding
  792. * probability distribution. Note that A * A^t and A^t * A have the
  793. * same trace.
  794. *
  795. ENTRAT = ZERO
  796. DO 1114 p = N+1, N+M
  797. BIG1 = ( ( WORK(p) / XSC )**2 ) * TEMP1
  798. IF ( BIG1 .NE. ZERO ) ENTRAT = ENTRAT + BIG1 * ALOG(BIG1)
  799. 1114 CONTINUE
  800. ENTRAT = - ENTRAT / ALOG(FLOAT(M))
  801. *
  802. * Analyze the entropies and decide A or A^t. Smaller entropy
  803. * usually means better input for the algorithm.
  804. *
  805. TRANSP = ( ENTRAT .LT. ENTRA )
  806. *
  807. * If A^t is better than A, transpose A.
  808. *
  809. IF ( TRANSP ) THEN
  810. * In an optimal implementation, this trivial transpose
  811. * should be replaced with faster transpose.
  812. DO 1115 p = 1, N - 1
  813. DO 1116 q = p + 1, N
  814. TEMP1 = A(q,p)
  815. A(q,p) = A(p,q)
  816. A(p,q) = TEMP1
  817. 1116 CONTINUE
  818. 1115 CONTINUE
  819. DO 1117 p = 1, N
  820. WORK(M+N+p) = SVA(p)
  821. SVA(p) = WORK(N+p)
  822. 1117 CONTINUE
  823. TEMP1 = AAPP
  824. AAPP = AATMAX
  825. AATMAX = TEMP1
  826. TEMP1 = AAQQ
  827. AAQQ = AATMIN
  828. AATMIN = TEMP1
  829. KILL = LSVEC
  830. LSVEC = RSVEC
  831. RSVEC = KILL
  832. IF ( LSVEC ) N1 = N
  833. *
  834. ROWPIV = .TRUE.
  835. END IF
  836. *
  837. END IF
  838. * END IF L2TRAN
  839. *
  840. * Scale the matrix so that its maximal singular value remains less
  841. * than SQRT(BIG) -- the matrix is scaled so that its maximal column
  842. * has Euclidean norm equal to SQRT(BIG/N). The only reason to keep
  843. * SQRT(BIG) instead of BIG is the fact that SGEJSV uses LAPACK and
  844. * BLAS routines that, in some implementations, are not capable of
  845. * working in the full interval [SFMIN,BIG] and that they may provoke
  846. * overflows in the intermediate results. If the singular values spread
  847. * from SFMIN to BIG, then SGESVJ will compute them. So, in that case,
  848. * one should use SGESVJ instead of SGEJSV.
  849. *
  850. BIG1 = SQRT( BIG )
  851. TEMP1 = SQRT( BIG / FLOAT(N) )
  852. *
  853. CALL SLASCL( 'G', 0, 0, AAPP, TEMP1, N, 1, SVA, N, IERR )
  854. IF ( AAQQ .GT. (AAPP * SFMIN) ) THEN
  855. AAQQ = ( AAQQ / AAPP ) * TEMP1
  856. ELSE
  857. AAQQ = ( AAQQ * TEMP1 ) / AAPP
  858. END IF
  859. TEMP1 = TEMP1 * SCALEM
  860. CALL SLASCL( 'G', 0, 0, AAPP, TEMP1, M, N, A, LDA, IERR )
  861. *
  862. * To undo scaling at the end of this procedure, multiply the
  863. * computed singular values with USCAL2 / USCAL1.
  864. *
  865. USCAL1 = TEMP1
  866. USCAL2 = AAPP
  867. *
  868. IF ( L2KILL ) THEN
  869. * L2KILL enforces computation of nonzero singular values in
  870. * the restricted range of condition number of the initial A,
  871. * sigma_max(A) / sigma_min(A) approx. SQRT(BIG)/SQRT(SFMIN).
  872. XSC = SQRT( SFMIN )
  873. ELSE
  874. XSC = SMALL
  875. *
  876. * Now, if the condition number of A is too big,
  877. * sigma_max(A) / sigma_min(A) .GT. SQRT(BIG/N) * EPSLN / SFMIN,
  878. * as a precaution measure, the full SVD is computed using SGESVJ
  879. * with accumulated Jacobi rotations. This provides numerically
  880. * more robust computation, at the cost of slightly increased run
  881. * time. Depending on the concrete implementation of BLAS and LAPACK
  882. * (i.e. how they behave in presence of extreme ill-conditioning) the
  883. * implementor may decide to remove this switch.
  884. IF ( ( AAQQ.LT.SQRT(SFMIN) ) .AND. LSVEC .AND. RSVEC ) THEN
  885. JRACC = .TRUE.
  886. END IF
  887. *
  888. END IF
  889. IF ( AAQQ .LT. XSC ) THEN
  890. DO 700 p = 1, N
  891. IF ( SVA(p) .LT. XSC ) THEN
  892. CALL SLASET( 'A', M, 1, ZERO, ZERO, A(1,p), LDA )
  893. SVA(p) = ZERO
  894. END IF
  895. 700 CONTINUE
  896. END IF
  897. *
  898. * Preconditioning using QR factorization with pivoting
  899. *
  900. IF ( ROWPIV ) THEN
  901. * Optional row permutation (Bjoerck row pivoting):
  902. * A result by Cox and Higham shows that the Bjoerck's
  903. * row pivoting combined with standard column pivoting
  904. * has similar effect as Powell-Reid complete pivoting.
  905. * The ell-infinity norms of A are made nonincreasing.
  906. DO 1952 p = 1, M - 1
  907. q = ISAMAX( M-p+1, WORK(M+N+p), 1 ) + p - 1
  908. IWORK(2*N+p) = q
  909. IF ( p .NE. q ) THEN
  910. TEMP1 = WORK(M+N+p)
  911. WORK(M+N+p) = WORK(M+N+q)
  912. WORK(M+N+q) = TEMP1
  913. END IF
  914. 1952 CONTINUE
  915. CALL SLASWP( N, A, LDA, 1, M-1, IWORK(2*N+1), 1 )
  916. END IF
  917. *
  918. * End of the preparation phase (scaling, optional sorting and
  919. * transposing, optional flushing of small columns).
  920. *
  921. * Preconditioning
  922. *
  923. * If the full SVD is needed, the right singular vectors are computed
  924. * from a matrix equation, and for that we need theoretical analysis
  925. * of the Businger-Golub pivoting. So we use SGEQP3 as the first RR QRF.
  926. * In all other cases the first RR QRF can be chosen by other criteria
  927. * (eg speed by replacing global with restricted window pivoting, such
  928. * as in SGEQPX from TOMS # 782). Good results will be obtained using
  929. * SGEQPX with properly (!) chosen numerical parameters.
  930. * Any improvement of SGEQP3 improves overal performance of SGEJSV.
  931. *
  932. * A * P1 = Q1 * [ R1^t 0]^t:
  933. DO 1963 p = 1, N
  934. * .. all columns are free columns
  935. IWORK(p) = 0
  936. 1963 CONTINUE
  937. CALL SGEQP3( M,N,A,LDA, IWORK,WORK, WORK(N+1),LWORK-N, IERR )
  938. *
  939. * The upper triangular matrix R1 from the first QRF is inspected for
  940. * rank deficiency and possibilities for deflation, or possible
  941. * ill-conditioning. Depending on the user specified flag L2RANK,
  942. * the procedure explores possibilities to reduce the numerical
  943. * rank by inspecting the computed upper triangular factor. If
  944. * L2RANK or L2ABER are up, then SGEJSV will compute the SVD of
  945. * A + dA, where ||dA|| <= f(M,N)*EPSLN.
  946. *
  947. NR = 1
  948. IF ( L2ABER ) THEN
  949. * Standard absolute error bound suffices. All sigma_i with
  950. * sigma_i < N*EPSLN*||A|| are flushed to zero. This is an
  951. * agressive enforcement of lower numerical rank by introducing a
  952. * backward error of the order of N*EPSLN*||A||.
  953. TEMP1 = SQRT(FLOAT(N))*EPSLN
  954. DO 3001 p = 2, N
  955. IF ( ABS(A(p,p)) .GE. (TEMP1*ABS(A(1,1))) ) THEN
  956. NR = NR + 1
  957. ELSE
  958. GO TO 3002
  959. END IF
  960. 3001 CONTINUE
  961. 3002 CONTINUE
  962. ELSE IF ( L2RANK ) THEN
  963. * .. similarly as above, only slightly more gentle (less agressive).
  964. * Sudden drop on the diagonal of R1 is used as the criterion for
  965. * close-to-rank-defficient.
  966. TEMP1 = SQRT(SFMIN)
  967. DO 3401 p = 2, N
  968. IF ( ( ABS(A(p,p)) .LT. (EPSLN*ABS(A(p-1,p-1))) ) .OR.
  969. $ ( ABS(A(p,p)) .LT. SMALL ) .OR.
  970. $ ( L2KILL .AND. (ABS(A(p,p)) .LT. TEMP1) ) ) GO TO 3402
  971. NR = NR + 1
  972. 3401 CONTINUE
  973. 3402 CONTINUE
  974. *
  975. ELSE
  976. * The goal is high relative accuracy. However, if the matrix
  977. * has high scaled condition number the relative accuracy is in
  978. * general not feasible. Later on, a condition number estimator
  979. * will be deployed to estimate the scaled condition number.
  980. * Here we just remove the underflowed part of the triangular
  981. * factor. This prevents the situation in which the code is
  982. * working hard to get the accuracy not warranted by the data.
  983. TEMP1 = SQRT(SFMIN)
  984. DO 3301 p = 2, N
  985. IF ( ( ABS(A(p,p)) .LT. SMALL ) .OR.
  986. $ ( L2KILL .AND. (ABS(A(p,p)) .LT. TEMP1) ) ) GO TO 3302
  987. NR = NR + 1
  988. 3301 CONTINUE
  989. 3302 CONTINUE
  990. *
  991. END IF
  992. *
  993. ALMORT = .FALSE.
  994. IF ( NR .EQ. N ) THEN
  995. MAXPRJ = ONE
  996. DO 3051 p = 2, N
  997. TEMP1 = ABS(A(p,p)) / SVA(IWORK(p))
  998. MAXPRJ = MIN( MAXPRJ, TEMP1 )
  999. 3051 CONTINUE
  1000. IF ( MAXPRJ**2 .GE. ONE - FLOAT(N)*EPSLN ) ALMORT = .TRUE.
  1001. END IF
  1002. *
  1003. *
  1004. SCONDA = - ONE
  1005. CONDR1 = - ONE
  1006. CONDR2 = - ONE
  1007. *
  1008. IF ( ERREST ) THEN
  1009. IF ( N .EQ. NR ) THEN
  1010. IF ( RSVEC ) THEN
  1011. * .. V is available as workspace
  1012. CALL SLACPY( 'U', N, N, A, LDA, V, LDV )
  1013. DO 3053 p = 1, N
  1014. TEMP1 = SVA(IWORK(p))
  1015. CALL SSCAL( p, ONE/TEMP1, V(1,p), 1 )
  1016. 3053 CONTINUE
  1017. CALL SPOCON( 'U', N, V, LDV, ONE, TEMP1,
  1018. $ WORK(N+1), IWORK(2*N+M+1), IERR )
  1019. ELSE IF ( LSVEC ) THEN
  1020. * .. U is available as workspace
  1021. CALL SLACPY( 'U', N, N, A, LDA, U, LDU )
  1022. DO 3054 p = 1, N
  1023. TEMP1 = SVA(IWORK(p))
  1024. CALL SSCAL( p, ONE/TEMP1, U(1,p), 1 )
  1025. 3054 CONTINUE
  1026. CALL SPOCON( 'U', N, U, LDU, ONE, TEMP1,
  1027. $ WORK(N+1), IWORK(2*N+M+1), IERR )
  1028. ELSE
  1029. CALL SLACPY( 'U', N, N, A, LDA, WORK(N+1), N )
  1030. DO 3052 p = 1, N
  1031. TEMP1 = SVA(IWORK(p))
  1032. CALL SSCAL( p, ONE/TEMP1, WORK(N+(p-1)*N+1), 1 )
  1033. 3052 CONTINUE
  1034. * .. the columns of R are scaled to have unit Euclidean lengths.
  1035. CALL SPOCON( 'U', N, WORK(N+1), N, ONE, TEMP1,
  1036. $ WORK(N+N*N+1), IWORK(2*N+M+1), IERR )
  1037. END IF
  1038. SCONDA = ONE / SQRT(TEMP1)
  1039. * SCONDA is an estimate of SQRT(||(R^t * R)^(-1)||_1).
  1040. * N^(-1/4) * SCONDA <= ||R^(-1)||_2 <= N^(1/4) * SCONDA
  1041. ELSE
  1042. SCONDA = - ONE
  1043. END IF
  1044. END IF
  1045. *
  1046. L2PERT = L2PERT .AND. ( ABS( A(1,1)/A(NR,NR) ) .GT. SQRT(BIG1) )
  1047. * If there is no violent scaling, artificial perturbation is not needed.
  1048. *
  1049. * Phase 3:
  1050. *
  1051. IF ( .NOT. ( RSVEC .OR. LSVEC ) ) THEN
  1052. *
  1053. * Singular Values only
  1054. *
  1055. * .. transpose A(1:NR,1:N)
  1056. DO 1946 p = 1, MIN( N-1, NR )
  1057. CALL SCOPY( N-p, A(p,p+1), LDA, A(p+1,p), 1 )
  1058. 1946 CONTINUE
  1059. *
  1060. * The following two DO-loops introduce small relative perturbation
  1061. * into the strict upper triangle of the lower triangular matrix.
  1062. * Small entries below the main diagonal are also changed.
  1063. * This modification is useful if the computing environment does not
  1064. * provide/allow FLUSH TO ZERO underflow, for it prevents many
  1065. * annoying denormalized numbers in case of strongly scaled matrices.
  1066. * The perturbation is structured so that it does not introduce any
  1067. * new perturbation of the singular values, and it does not destroy
  1068. * the job done by the preconditioner.
  1069. * The licence for this perturbation is in the variable L2PERT, which
  1070. * should be .FALSE. if FLUSH TO ZERO underflow is active.
  1071. *
  1072. IF ( .NOT. ALMORT ) THEN
  1073. *
  1074. IF ( L2PERT ) THEN
  1075. * XSC = SQRT(SMALL)
  1076. XSC = EPSLN / FLOAT(N)
  1077. DO 4947 q = 1, NR
  1078. TEMP1 = XSC*ABS(A(q,q))
  1079. DO 4949 p = 1, N
  1080. IF ( ( (p.GT.q) .AND. (ABS(A(p,q)).LE.TEMP1) )
  1081. $ .OR. ( p .LT. q ) )
  1082. $ A(p,q) = SIGN( TEMP1, A(p,q) )
  1083. 4949 CONTINUE
  1084. 4947 CONTINUE
  1085. ELSE
  1086. CALL SLASET( 'U', NR-1,NR-1, ZERO,ZERO, A(1,2),LDA )
  1087. END IF
  1088. *
  1089. * .. second preconditioning using the QR factorization
  1090. *
  1091. CALL SGEQRF( N,NR, A,LDA, WORK, WORK(N+1),LWORK-N, IERR )
  1092. *
  1093. * .. and transpose upper to lower triangular
  1094. DO 1948 p = 1, NR - 1
  1095. CALL SCOPY( NR-p, A(p,p+1), LDA, A(p+1,p), 1 )
  1096. 1948 CONTINUE
  1097. *
  1098. END IF
  1099. *
  1100. * Row-cyclic Jacobi SVD algorithm with column pivoting
  1101. *
  1102. * .. again some perturbation (a "background noise") is added
  1103. * to drown denormals
  1104. IF ( L2PERT ) THEN
  1105. * XSC = SQRT(SMALL)
  1106. XSC = EPSLN / FLOAT(N)
  1107. DO 1947 q = 1, NR
  1108. TEMP1 = XSC*ABS(A(q,q))
  1109. DO 1949 p = 1, NR
  1110. IF ( ( (p.GT.q) .AND. (ABS(A(p,q)).LE.TEMP1) )
  1111. $ .OR. ( p .LT. q ) )
  1112. $ A(p,q) = SIGN( TEMP1, A(p,q) )
  1113. 1949 CONTINUE
  1114. 1947 CONTINUE
  1115. ELSE
  1116. CALL SLASET( 'U', NR-1, NR-1, ZERO, ZERO, A(1,2), LDA )
  1117. END IF
  1118. *
  1119. * .. and one-sided Jacobi rotations are started on a lower
  1120. * triangular matrix (plus perturbation which is ignored in
  1121. * the part which destroys triangular form (confusing?!))
  1122. *
  1123. CALL SGESVJ( 'L', 'NoU', 'NoV', NR, NR, A, LDA, SVA,
  1124. $ N, V, LDV, WORK, LWORK, INFO )
  1125. *
  1126. SCALEM = WORK(1)
  1127. NUMRANK = NINT(WORK(2))
  1128. *
  1129. *
  1130. ELSE IF ( RSVEC .AND. ( .NOT. LSVEC ) ) THEN
  1131. *
  1132. * -> Singular Values and Right Singular Vectors <-
  1133. *
  1134. IF ( ALMORT ) THEN
  1135. *
  1136. * .. in this case NR equals N
  1137. DO 1998 p = 1, NR
  1138. CALL SCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
  1139. 1998 CONTINUE
  1140. CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
  1141. *
  1142. CALL SGESVJ( 'L','U','N', N, NR, V,LDV, SVA, NR, A,LDA,
  1143. $ WORK, LWORK, INFO )
  1144. SCALEM = WORK(1)
  1145. NUMRANK = NINT(WORK(2))
  1146. ELSE
  1147. *
  1148. * .. two more QR factorizations ( one QRF is not enough, two require
  1149. * accumulated product of Jacobi rotations, three are perfect )
  1150. *
  1151. CALL SLASET( 'Lower', NR-1, NR-1, ZERO, ZERO, A(2,1), LDA )
  1152. CALL SGELQF( NR, N, A, LDA, WORK, WORK(N+1), LWORK-N, IERR)
  1153. CALL SLACPY( 'Lower', NR, NR, A, LDA, V, LDV )
  1154. CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
  1155. CALL SGEQRF( NR, NR, V, LDV, WORK(N+1), WORK(2*N+1),
  1156. $ LWORK-2*N, IERR )
  1157. DO 8998 p = 1, NR
  1158. CALL SCOPY( NR-p+1, V(p,p), LDV, V(p,p), 1 )
  1159. 8998 CONTINUE
  1160. CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
  1161. *
  1162. CALL SGESVJ( 'Lower', 'U','N', NR, NR, V,LDV, SVA, NR, U,
  1163. $ LDU, WORK(N+1), LWORK-N, INFO )
  1164. SCALEM = WORK(N+1)
  1165. NUMRANK = NINT(WORK(N+2))
  1166. IF ( NR .LT. N ) THEN
  1167. CALL SLASET( 'A',N-NR, NR, ZERO,ZERO, V(NR+1,1), LDV )
  1168. CALL SLASET( 'A',NR, N-NR, ZERO,ZERO, V(1,NR+1), LDV )
  1169. CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE, V(NR+1,NR+1), LDV )
  1170. END IF
  1171. *
  1172. CALL SORMLQ( 'Left', 'Transpose', N, N, NR, A, LDA, WORK,
  1173. $ V, LDV, WORK(N+1), LWORK-N, IERR )
  1174. *
  1175. END IF
  1176. *
  1177. DO 8991 p = 1, N
  1178. CALL SCOPY( N, V(p,1), LDV, A(IWORK(p),1), LDA )
  1179. 8991 CONTINUE
  1180. CALL SLACPY( 'All', N, N, A, LDA, V, LDV )
  1181. *
  1182. IF ( TRANSP ) THEN
  1183. CALL SLACPY( 'All', N, N, V, LDV, U, LDU )
  1184. END IF
  1185. *
  1186. ELSE IF ( LSVEC .AND. ( .NOT. RSVEC ) ) THEN
  1187. *
  1188. * .. Singular Values and Left Singular Vectors ..
  1189. *
  1190. * .. second preconditioning step to avoid need to accumulate
  1191. * Jacobi rotations in the Jacobi iterations.
  1192. DO 1965 p = 1, NR
  1193. CALL SCOPY( N-p+1, A(p,p), LDA, U(p,p), 1 )
  1194. 1965 CONTINUE
  1195. CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, U(1,2), LDU )
  1196. *
  1197. CALL SGEQRF( N, NR, U, LDU, WORK(N+1), WORK(2*N+1),
  1198. $ LWORK-2*N, IERR )
  1199. *
  1200. DO 1967 p = 1, NR - 1
  1201. CALL SCOPY( NR-p, U(p,p+1), LDU, U(p+1,p), 1 )
  1202. 1967 CONTINUE
  1203. CALL SLASET( 'Upper', NR-1, NR-1, ZERO, ZERO, U(1,2), LDU )
  1204. *
  1205. CALL SGESVJ( 'Lower', 'U', 'N', NR,NR, U, LDU, SVA, NR, A,
  1206. $ LDA, WORK(N+1), LWORK-N, INFO )
  1207. SCALEM = WORK(N+1)
  1208. NUMRANK = NINT(WORK(N+2))
  1209. *
  1210. IF ( NR .LT. M ) THEN
  1211. CALL SLASET( 'A', M-NR, NR,ZERO, ZERO, U(NR+1,1), LDU )
  1212. IF ( NR .LT. N1 ) THEN
  1213. CALL SLASET( 'A',NR, N1-NR, ZERO, ZERO, U(1,NR+1), LDU )
  1214. CALL SLASET( 'A',M-NR,N1-NR,ZERO,ONE,U(NR+1,NR+1), LDU )
  1215. END IF
  1216. END IF
  1217. *
  1218. CALL SORMQR( 'Left', 'No Tr', M, N1, N, A, LDA, WORK, U,
  1219. $ LDU, WORK(N+1), LWORK-N, IERR )
  1220. *
  1221. IF ( ROWPIV )
  1222. $ CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
  1223. *
  1224. DO 1974 p = 1, N1
  1225. XSC = ONE / SNRM2( M, U(1,p), 1 )
  1226. CALL SSCAL( M, XSC, U(1,p), 1 )
  1227. 1974 CONTINUE
  1228. *
  1229. IF ( TRANSP ) THEN
  1230. CALL SLACPY( 'All', N, N, U, LDU, V, LDV )
  1231. END IF
  1232. *
  1233. ELSE
  1234. *
  1235. * .. Full SVD ..
  1236. *
  1237. IF ( .NOT. JRACC ) THEN
  1238. *
  1239. IF ( .NOT. ALMORT ) THEN
  1240. *
  1241. * Second Preconditioning Step (QRF [with pivoting])
  1242. * Note that the composition of TRANSPOSE, QRF and TRANSPOSE is
  1243. * equivalent to an LQF CALL. Since in many libraries the QRF
  1244. * seems to be better optimized than the LQF, we do explicit
  1245. * transpose and use the QRF. This is subject to changes in an
  1246. * optimized implementation of SGEJSV.
  1247. *
  1248. DO 1968 p = 1, NR
  1249. CALL SCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
  1250. 1968 CONTINUE
  1251. *
  1252. * .. the following two loops perturb small entries to avoid
  1253. * denormals in the second QR factorization, where they are
  1254. * as good as zeros. This is done to avoid painfully slow
  1255. * computation with denormals. The relative size of the perturbation
  1256. * is a parameter that can be changed by the implementer.
  1257. * This perturbation device will be obsolete on machines with
  1258. * properly implemented arithmetic.
  1259. * To switch it off, set L2PERT=.FALSE. To remove it from the
  1260. * code, remove the action under L2PERT=.TRUE., leave the ELSE part.
  1261. * The following two loops should be blocked and fused with the
  1262. * transposed copy above.
  1263. *
  1264. IF ( L2PERT ) THEN
  1265. XSC = SQRT(SMALL)
  1266. DO 2969 q = 1, NR
  1267. TEMP1 = XSC*ABS( V(q,q) )
  1268. DO 2968 p = 1, N
  1269. IF ( ( p .GT. q ) .AND. ( ABS(V(p,q)) .LE. TEMP1 )
  1270. $ .OR. ( p .LT. q ) )
  1271. $ V(p,q) = SIGN( TEMP1, V(p,q) )
  1272. IF ( p .LT. q ) V(p,q) = - V(p,q)
  1273. 2968 CONTINUE
  1274. 2969 CONTINUE
  1275. ELSE
  1276. CALL SLASET( 'U', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
  1277. END IF
  1278. *
  1279. * Estimate the row scaled condition number of R1
  1280. * (If R1 is rectangular, N > NR, then the condition number
  1281. * of the leading NR x NR submatrix is estimated.)
  1282. *
  1283. CALL SLACPY( 'L', NR, NR, V, LDV, WORK(2*N+1), NR )
  1284. DO 3950 p = 1, NR
  1285. TEMP1 = SNRM2(NR-p+1,WORK(2*N+(p-1)*NR+p),1)
  1286. CALL SSCAL(NR-p+1,ONE/TEMP1,WORK(2*N+(p-1)*NR+p),1)
  1287. 3950 CONTINUE
  1288. CALL SPOCON('Lower',NR,WORK(2*N+1),NR,ONE,TEMP1,
  1289. $ WORK(2*N+NR*NR+1),IWORK(M+2*N+1),IERR)
  1290. CONDR1 = ONE / SQRT(TEMP1)
  1291. * .. here need a second oppinion on the condition number
  1292. * .. then assume worst case scenario
  1293. * R1 is OK for inverse <=> CONDR1 .LT. FLOAT(N)
  1294. * more conservative <=> CONDR1 .LT. SQRT(FLOAT(N))
  1295. *
  1296. COND_OK = SQRT(FLOAT(NR))
  1297. *[TP] COND_OK is a tuning parameter.
  1298. IF ( CONDR1 .LT. COND_OK ) THEN
  1299. * .. the second QRF without pivoting. Note: in an optimized
  1300. * implementation, this QRF should be implemented as the QRF
  1301. * of a lower triangular matrix.
  1302. * R1^t = Q2 * R2
  1303. CALL SGEQRF( N, NR, V, LDV, WORK(N+1), WORK(2*N+1),
  1304. $ LWORK-2*N, IERR )
  1305. *
  1306. IF ( L2PERT ) THEN
  1307. XSC = SQRT(SMALL)/EPSLN
  1308. DO 3959 p = 2, NR
  1309. DO 3958 q = 1, p - 1
  1310. TEMP1 = XSC * MIN(ABS(V(p,p)),ABS(V(q,q)))
  1311. IF ( ABS(V(q,p)) .LE. TEMP1 )
  1312. $ V(q,p) = SIGN( TEMP1, V(q,p) )
  1313. 3958 CONTINUE
  1314. 3959 CONTINUE
  1315. END IF
  1316. *
  1317. IF ( NR .NE. N )
  1318. $ CALL SLACPY( 'A', N, NR, V, LDV, WORK(2*N+1), N )
  1319. * .. save ...
  1320. *
  1321. * .. this transposed copy should be better than naive
  1322. DO 1969 p = 1, NR - 1
  1323. CALL SCOPY( NR-p, V(p,p+1), LDV, V(p+1,p), 1 )
  1324. 1969 CONTINUE
  1325. *
  1326. CONDR2 = CONDR1
  1327. *
  1328. ELSE
  1329. *
  1330. * .. ill-conditioned case: second QRF with pivoting
  1331. * Note that windowed pivoting would be equaly good
  1332. * numerically, and more run-time efficient. So, in
  1333. * an optimal implementation, the next call to SGEQP3
  1334. * should be replaced with eg. CALL SGEQPX (ACM TOMS #782)
  1335. * with properly (carefully) chosen parameters.
  1336. *
  1337. * R1^t * P2 = Q2 * R2
  1338. DO 3003 p = 1, NR
  1339. IWORK(N+p) = 0
  1340. 3003 CONTINUE
  1341. CALL SGEQP3( N, NR, V, LDV, IWORK(N+1), WORK(N+1),
  1342. $ WORK(2*N+1), LWORK-2*N, IERR )
  1343. ** CALL SGEQRF( N, NR, V, LDV, WORK(N+1), WORK(2*N+1),
  1344. ** $ LWORK-2*N, IERR )
  1345. IF ( L2PERT ) THEN
  1346. XSC = SQRT(SMALL)
  1347. DO 3969 p = 2, NR
  1348. DO 3968 q = 1, p - 1
  1349. TEMP1 = XSC * MIN(ABS(V(p,p)),ABS(V(q,q)))
  1350. IF ( ABS(V(q,p)) .LE. TEMP1 )
  1351. $ V(q,p) = SIGN( TEMP1, V(q,p) )
  1352. 3968 CONTINUE
  1353. 3969 CONTINUE
  1354. END IF
  1355. *
  1356. CALL SLACPY( 'A', N, NR, V, LDV, WORK(2*N+1), N )
  1357. *
  1358. IF ( L2PERT ) THEN
  1359. XSC = SQRT(SMALL)
  1360. DO 8970 p = 2, NR
  1361. DO 8971 q = 1, p - 1
  1362. TEMP1 = XSC * MIN(ABS(V(p,p)),ABS(V(q,q)))
  1363. V(p,q) = - SIGN( TEMP1, V(q,p) )
  1364. 8971 CONTINUE
  1365. 8970 CONTINUE
  1366. ELSE
  1367. CALL SLASET( 'L',NR-1,NR-1,ZERO,ZERO,V(2,1),LDV )
  1368. END IF
  1369. * Now, compute R2 = L3 * Q3, the LQ factorization.
  1370. CALL SGELQF( NR, NR, V, LDV, WORK(2*N+N*NR+1),
  1371. $ WORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR, IERR )
  1372. * .. and estimate the condition number
  1373. CALL SLACPY( 'L',NR,NR,V,LDV,WORK(2*N+N*NR+NR+1),NR )
  1374. DO 4950 p = 1, NR
  1375. TEMP1 = SNRM2( p, WORK(2*N+N*NR+NR+p), NR )
  1376. CALL SSCAL( p, ONE/TEMP1, WORK(2*N+N*NR+NR+p), NR )
  1377. 4950 CONTINUE
  1378. CALL SPOCON( 'L',NR,WORK(2*N+N*NR+NR+1),NR,ONE,TEMP1,
  1379. $ WORK(2*N+N*NR+NR+NR*NR+1),IWORK(M+2*N+1),IERR )
  1380. CONDR2 = ONE / SQRT(TEMP1)
  1381. *
  1382. IF ( CONDR2 .GE. COND_OK ) THEN
  1383. * .. save the Householder vectors used for Q3
  1384. * (this overwrittes the copy of R2, as it will not be
  1385. * needed in this branch, but it does not overwritte the
  1386. * Huseholder vectors of Q2.).
  1387. CALL SLACPY( 'U', NR, NR, V, LDV, WORK(2*N+1), N )
  1388. * .. and the rest of the information on Q3 is in
  1389. * WORK(2*N+N*NR+1:2*N+N*NR+N)
  1390. END IF
  1391. *
  1392. END IF
  1393. *
  1394. IF ( L2PERT ) THEN
  1395. XSC = SQRT(SMALL)
  1396. DO 4968 q = 2, NR
  1397. TEMP1 = XSC * V(q,q)
  1398. DO 4969 p = 1, q - 1
  1399. * V(p,q) = - SIGN( TEMP1, V(q,p) )
  1400. V(p,q) = - SIGN( TEMP1, V(p,q) )
  1401. 4969 CONTINUE
  1402. 4968 CONTINUE
  1403. ELSE
  1404. CALL SLASET( 'U', NR-1,NR-1, ZERO,ZERO, V(1,2), LDV )
  1405. END IF
  1406. *
  1407. * Second preconditioning finished; continue with Jacobi SVD
  1408. * The input matrix is lower trinagular.
  1409. *
  1410. * Recover the right singular vectors as solution of a well
  1411. * conditioned triangular matrix equation.
  1412. *
  1413. IF ( CONDR1 .LT. COND_OK ) THEN
  1414. *
  1415. CALL SGESVJ( 'L','U','N',NR,NR,V,LDV,SVA,NR,U,
  1416. $ LDU,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,INFO )
  1417. SCALEM = WORK(2*N+N*NR+NR+1)
  1418. NUMRANK = NINT(WORK(2*N+N*NR+NR+2))
  1419. DO 3970 p = 1, NR
  1420. CALL SCOPY( NR, V(1,p), 1, U(1,p), 1 )
  1421. CALL SSCAL( NR, SVA(p), V(1,p), 1 )
  1422. 3970 CONTINUE
  1423. * .. pick the right matrix equation and solve it
  1424. *
  1425. IF ( NR .EQ. N ) THEN
  1426. * :)) .. best case, R1 is inverted. The solution of this matrix
  1427. * equation is Q2*V2 = the product of the Jacobi rotations
  1428. * used in SGESVJ, premultiplied with the orthogonal matrix
  1429. * from the second QR factorization.
  1430. CALL STRSM( 'L','U','N','N', NR,NR,ONE, A,LDA, V,LDV )
  1431. ELSE
  1432. * .. R1 is well conditioned, but non-square. Transpose(R2)
  1433. * is inverted to get the product of the Jacobi rotations
  1434. * used in SGESVJ. The Q-factor from the second QR
  1435. * factorization is then built in explicitly.
  1436. CALL STRSM('L','U','T','N',NR,NR,ONE,WORK(2*N+1),
  1437. $ N,V,LDV)
  1438. IF ( NR .LT. N ) THEN
  1439. CALL SLASET('A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV)
  1440. CALL SLASET('A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV)
  1441. CALL SLASET('A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV)
  1442. END IF
  1443. CALL SORMQR('L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
  1444. $ V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR)
  1445. END IF
  1446. *
  1447. ELSE IF ( CONDR2 .LT. COND_OK ) THEN
  1448. *
  1449. * :) .. the input matrix A is very likely a relative of
  1450. * the Kahan matrix :)
  1451. * The matrix R2 is inverted. The solution of the matrix equation
  1452. * is Q3^T*V3 = the product of the Jacobi rotations (appplied to
  1453. * the lower triangular L3 from the LQ factorization of
  1454. * R2=L3*Q3), pre-multiplied with the transposed Q3.
  1455. CALL SGESVJ( 'L', 'U', 'N', NR, NR, V, LDV, SVA, NR, U,
  1456. $ LDU, WORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR, INFO )
  1457. SCALEM = WORK(2*N+N*NR+NR+1)
  1458. NUMRANK = NINT(WORK(2*N+N*NR+NR+2))
  1459. DO 3870 p = 1, NR
  1460. CALL SCOPY( NR, V(1,p), 1, U(1,p), 1 )
  1461. CALL SSCAL( NR, SVA(p), U(1,p), 1 )
  1462. 3870 CONTINUE
  1463. CALL STRSM('L','U','N','N',NR,NR,ONE,WORK(2*N+1),N,U,LDU)
  1464. * .. apply the permutation from the second QR factorization
  1465. DO 873 q = 1, NR
  1466. DO 872 p = 1, NR
  1467. WORK(2*N+N*NR+NR+IWORK(N+p)) = U(p,q)
  1468. 872 CONTINUE
  1469. DO 874 p = 1, NR
  1470. U(p,q) = WORK(2*N+N*NR+NR+p)
  1471. 874 CONTINUE
  1472. 873 CONTINUE
  1473. IF ( NR .LT. N ) THEN
  1474. CALL SLASET( 'A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV )
  1475. CALL SLASET( 'A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV )
  1476. CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV )
  1477. END IF
  1478. CALL SORMQR( 'L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
  1479. $ V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
  1480. ELSE
  1481. * Last line of defense.
  1482. * #:( This is a rather pathological case: no scaled condition
  1483. * improvement after two pivoted QR factorizations. Other
  1484. * possibility is that the rank revealing QR factorization
  1485. * or the condition estimator has failed, or the COND_OK
  1486. * is set very close to ONE (which is unnecessary). Normally,
  1487. * this branch should never be executed, but in rare cases of
  1488. * failure of the RRQR or condition estimator, the last line of
  1489. * defense ensures that SGEJSV completes the task.
  1490. * Compute the full SVD of L3 using SGESVJ with explicit
  1491. * accumulation of Jacobi rotations.
  1492. CALL SGESVJ( 'L', 'U', 'V', NR, NR, V, LDV, SVA, NR, U,
  1493. $ LDU, WORK(2*N+N*NR+NR+1), LWORK-2*N-N*NR-NR, INFO )
  1494. SCALEM = WORK(2*N+N*NR+NR+1)
  1495. NUMRANK = NINT(WORK(2*N+N*NR+NR+2))
  1496. IF ( NR .LT. N ) THEN
  1497. CALL SLASET( 'A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV )
  1498. CALL SLASET( 'A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV )
  1499. CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV )
  1500. END IF
  1501. CALL SORMQR( 'L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
  1502. $ V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
  1503. *
  1504. CALL SORMLQ( 'L', 'T', NR, NR, NR, WORK(2*N+1), N,
  1505. $ WORK(2*N+N*NR+1), U, LDU, WORK(2*N+N*NR+NR+1),
  1506. $ LWORK-2*N-N*NR-NR, IERR )
  1507. DO 773 q = 1, NR
  1508. DO 772 p = 1, NR
  1509. WORK(2*N+N*NR+NR+IWORK(N+p)) = U(p,q)
  1510. 772 CONTINUE
  1511. DO 774 p = 1, NR
  1512. U(p,q) = WORK(2*N+N*NR+NR+p)
  1513. 774 CONTINUE
  1514. 773 CONTINUE
  1515. *
  1516. END IF
  1517. *
  1518. * Permute the rows of V using the (column) permutation from the
  1519. * first QRF. Also, scale the columns to make them unit in
  1520. * Euclidean norm. This applies to all cases.
  1521. *
  1522. TEMP1 = SQRT(FLOAT(N)) * EPSLN
  1523. DO 1972 q = 1, N
  1524. DO 972 p = 1, N
  1525. WORK(2*N+N*NR+NR+IWORK(p)) = V(p,q)
  1526. 972 CONTINUE
  1527. DO 973 p = 1, N
  1528. V(p,q) = WORK(2*N+N*NR+NR+p)
  1529. 973 CONTINUE
  1530. XSC = ONE / SNRM2( N, V(1,q), 1 )
  1531. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  1532. $ CALL SSCAL( N, XSC, V(1,q), 1 )
  1533. 1972 CONTINUE
  1534. * At this moment, V contains the right singular vectors of A.
  1535. * Next, assemble the left singular vector matrix U (M x N).
  1536. IF ( NR .LT. M ) THEN
  1537. CALL SLASET( 'A', M-NR, NR, ZERO, ZERO, U(NR+1,1), LDU )
  1538. IF ( NR .LT. N1 ) THEN
  1539. CALL SLASET('A',NR,N1-NR,ZERO,ZERO,U(1,NR+1),LDU)
  1540. CALL SLASET('A',M-NR,N1-NR,ZERO,ONE,U(NR+1,NR+1),LDU)
  1541. END IF
  1542. END IF
  1543. *
  1544. * The Q matrix from the first QRF is built into the left singular
  1545. * matrix U. This applies to all cases.
  1546. *
  1547. CALL SORMQR( 'Left', 'No_Tr', M, N1, N, A, LDA, WORK, U,
  1548. $ LDU, WORK(N+1), LWORK-N, IERR )
  1549. * The columns of U are normalized. The cost is O(M*N) flops.
  1550. TEMP1 = SQRT(FLOAT(M)) * EPSLN
  1551. DO 1973 p = 1, NR
  1552. XSC = ONE / SNRM2( M, U(1,p), 1 )
  1553. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  1554. $ CALL SSCAL( M, XSC, U(1,p), 1 )
  1555. 1973 CONTINUE
  1556. *
  1557. * If the initial QRF is computed with row pivoting, the left
  1558. * singular vectors must be adjusted.
  1559. *
  1560. IF ( ROWPIV )
  1561. $ CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
  1562. *
  1563. ELSE
  1564. *
  1565. * .. the initial matrix A has almost orthogonal columns and
  1566. * the second QRF is not needed
  1567. *
  1568. CALL SLACPY( 'Upper', N, N, A, LDA, WORK(N+1), N )
  1569. IF ( L2PERT ) THEN
  1570. XSC = SQRT(SMALL)
  1571. DO 5970 p = 2, N
  1572. TEMP1 = XSC * WORK( N + (p-1)*N + p )
  1573. DO 5971 q = 1, p - 1
  1574. WORK(N+(q-1)*N+p)=-SIGN(TEMP1,WORK(N+(p-1)*N+q))
  1575. 5971 CONTINUE
  1576. 5970 CONTINUE
  1577. ELSE
  1578. CALL SLASET( 'Lower',N-1,N-1,ZERO,ZERO,WORK(N+2),N )
  1579. END IF
  1580. *
  1581. CALL SGESVJ( 'Upper', 'U', 'N', N, N, WORK(N+1), N, SVA,
  1582. $ N, U, LDU, WORK(N+N*N+1), LWORK-N-N*N, INFO )
  1583. *
  1584. SCALEM = WORK(N+N*N+1)
  1585. NUMRANK = NINT(WORK(N+N*N+2))
  1586. DO 6970 p = 1, N
  1587. CALL SCOPY( N, WORK(N+(p-1)*N+1), 1, U(1,p), 1 )
  1588. CALL SSCAL( N, SVA(p), WORK(N+(p-1)*N+1), 1 )
  1589. 6970 CONTINUE
  1590. *
  1591. CALL STRSM( 'Left', 'Upper', 'NoTrans', 'No UD', N, N,
  1592. $ ONE, A, LDA, WORK(N+1), N )
  1593. DO 6972 p = 1, N
  1594. CALL SCOPY( N, WORK(N+p), N, V(IWORK(p),1), LDV )
  1595. 6972 CONTINUE
  1596. TEMP1 = SQRT(FLOAT(N))*EPSLN
  1597. DO 6971 p = 1, N
  1598. XSC = ONE / SNRM2( N, V(1,p), 1 )
  1599. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  1600. $ CALL SSCAL( N, XSC, V(1,p), 1 )
  1601. 6971 CONTINUE
  1602. *
  1603. * Assemble the left singular vector matrix U (M x N).
  1604. *
  1605. IF ( N .LT. M ) THEN
  1606. CALL SLASET( 'A', M-N, N, ZERO, ZERO, U(N+1,1), LDU )
  1607. IF ( N .LT. N1 ) THEN
  1608. CALL SLASET( 'A',N, N1-N, ZERO, ZERO, U(1,N+1),LDU )
  1609. CALL SLASET( 'A',M-N,N1-N, ZERO, ONE,U(N+1,N+1),LDU )
  1610. END IF
  1611. END IF
  1612. CALL SORMQR( 'Left', 'No Tr', M, N1, N, A, LDA, WORK, U,
  1613. $ LDU, WORK(N+1), LWORK-N, IERR )
  1614. TEMP1 = SQRT(FLOAT(M))*EPSLN
  1615. DO 6973 p = 1, N1
  1616. XSC = ONE / SNRM2( M, U(1,p), 1 )
  1617. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  1618. $ CALL SSCAL( M, XSC, U(1,p), 1 )
  1619. 6973 CONTINUE
  1620. *
  1621. IF ( ROWPIV )
  1622. $ CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
  1623. *
  1624. END IF
  1625. *
  1626. * end of the >> almost orthogonal case << in the full SVD
  1627. *
  1628. ELSE
  1629. *
  1630. * This branch deploys a preconditioned Jacobi SVD with explicitly
  1631. * accumulated rotations. It is included as optional, mainly for
  1632. * experimental purposes. It does perfom well, and can also be used.
  1633. * In this implementation, this branch will be automatically activated
  1634. * if the condition number sigma_max(A) / sigma_min(A) is predicted
  1635. * to be greater than the overflow threshold. This is because the
  1636. * a posteriori computation of the singular vectors assumes robust
  1637. * implementation of BLAS and some LAPACK procedures, capable of working
  1638. * in presence of extreme values. Since that is not always the case, ...
  1639. *
  1640. DO 7968 p = 1, NR
  1641. CALL SCOPY( N-p+1, A(p,p), LDA, V(p,p), 1 )
  1642. 7968 CONTINUE
  1643. *
  1644. IF ( L2PERT ) THEN
  1645. XSC = SQRT(SMALL/EPSLN)
  1646. DO 5969 q = 1, NR
  1647. TEMP1 = XSC*ABS( V(q,q) )
  1648. DO 5968 p = 1, N
  1649. IF ( ( p .GT. q ) .AND. ( ABS(V(p,q)) .LE. TEMP1 )
  1650. $ .OR. ( p .LT. q ) )
  1651. $ V(p,q) = SIGN( TEMP1, V(p,q) )
  1652. IF ( p .LT. q ) V(p,q) = - V(p,q)
  1653. 5968 CONTINUE
  1654. 5969 CONTINUE
  1655. ELSE
  1656. CALL SLASET( 'U', NR-1, NR-1, ZERO, ZERO, V(1,2), LDV )
  1657. END IF
  1658. CALL SGEQRF( N, NR, V, LDV, WORK(N+1), WORK(2*N+1),
  1659. $ LWORK-2*N, IERR )
  1660. CALL SLACPY( 'L', N, NR, V, LDV, WORK(2*N+1), N )
  1661. *
  1662. DO 7969 p = 1, NR
  1663. CALL SCOPY( NR-p+1, V(p,p), LDV, U(p,p), 1 )
  1664. 7969 CONTINUE
  1665. IF ( L2PERT ) THEN
  1666. XSC = SQRT(SMALL/EPSLN)
  1667. DO 9970 q = 2, NR
  1668. DO 9971 p = 1, q - 1
  1669. TEMP1 = XSC * MIN(ABS(U(p,p)),ABS(U(q,q)))
  1670. U(p,q) = - SIGN( TEMP1, U(q,p) )
  1671. 9971 CONTINUE
  1672. 9970 CONTINUE
  1673. ELSE
  1674. CALL SLASET('U', NR-1, NR-1, ZERO, ZERO, U(1,2), LDU )
  1675. END IF
  1676. CALL SGESVJ( 'L', 'U', 'V', NR, NR, U, LDU, SVA,
  1677. $ N, V, LDV, WORK(2*N+N*NR+1), LWORK-2*N-N*NR, INFO )
  1678. SCALEM = WORK(2*N+N*NR+1)
  1679. NUMRANK = NINT(WORK(2*N+N*NR+2))
  1680. IF ( NR .LT. N ) THEN
  1681. CALL SLASET( 'A',N-NR,NR,ZERO,ZERO,V(NR+1,1),LDV )
  1682. CALL SLASET( 'A',NR,N-NR,ZERO,ZERO,V(1,NR+1),LDV )
  1683. CALL SLASET( 'A',N-NR,N-NR,ZERO,ONE,V(NR+1,NR+1),LDV )
  1684. END IF
  1685. CALL SORMQR( 'L','N',N,N,NR,WORK(2*N+1),N,WORK(N+1),
  1686. $ V,LDV,WORK(2*N+N*NR+NR+1),LWORK-2*N-N*NR-NR,IERR )
  1687. *
  1688. * Permute the rows of V using the (column) permutation from the
  1689. * first QRF. Also, scale the columns to make them unit in
  1690. * Euclidean norm. This applies to all cases.
  1691. *
  1692. TEMP1 = SQRT(FLOAT(N)) * EPSLN
  1693. DO 7972 q = 1, N
  1694. DO 8972 p = 1, N
  1695. WORK(2*N+N*NR+NR+IWORK(p)) = V(p,q)
  1696. 8972 CONTINUE
  1697. DO 8973 p = 1, N
  1698. V(p,q) = WORK(2*N+N*NR+NR+p)
  1699. 8973 CONTINUE
  1700. XSC = ONE / SNRM2( N, V(1,q), 1 )
  1701. IF ( (XSC .LT. (ONE-TEMP1)) .OR. (XSC .GT. (ONE+TEMP1)) )
  1702. $ CALL SSCAL( N, XSC, V(1,q), 1 )
  1703. 7972 CONTINUE
  1704. *
  1705. * At this moment, V contains the right singular vectors of A.
  1706. * Next, assemble the left singular vector matrix U (M x N).
  1707. *
  1708. IF ( NR .LT. M ) THEN
  1709. CALL SLASET( 'A', M-NR, NR, ZERO, ZERO, U(NR+1,1), LDU )
  1710. IF ( NR .LT. N1 ) THEN
  1711. CALL SLASET( 'A',NR, N1-NR, ZERO, ZERO, U(1,NR+1),LDU )
  1712. CALL SLASET( 'A',M-NR,N1-NR, ZERO, ONE,U(NR+1,NR+1),LDU )
  1713. END IF
  1714. END IF
  1715. *
  1716. CALL SORMQR( 'Left', 'No Tr', M, N1, N, A, LDA, WORK, U,
  1717. $ LDU, WORK(N+1), LWORK-N, IERR )
  1718. *
  1719. IF ( ROWPIV )
  1720. $ CALL SLASWP( N1, U, LDU, 1, M-1, IWORK(2*N+1), -1 )
  1721. *
  1722. *
  1723. END IF
  1724. IF ( TRANSP ) THEN
  1725. * .. swap U and V because the procedure worked on A^t
  1726. DO 6974 p = 1, N
  1727. CALL SSWAP( N, U(1,p), 1, V(1,p), 1 )
  1728. 6974 CONTINUE
  1729. END IF
  1730. *
  1731. END IF
  1732. * end of the full SVD
  1733. *
  1734. * Undo scaling, if necessary (and possible)
  1735. *
  1736. IF ( USCAL2 .LE. (BIG/SVA(1))*USCAL1 ) THEN
  1737. CALL SLASCL( 'G', 0, 0, USCAL1, USCAL2, NR, 1, SVA, N, IERR )
  1738. USCAL1 = ONE
  1739. USCAL2 = ONE
  1740. END IF
  1741. *
  1742. IF ( NR .LT. N ) THEN
  1743. DO 3004 p = NR+1, N
  1744. SVA(p) = ZERO
  1745. 3004 CONTINUE
  1746. END IF
  1747. *
  1748. WORK(1) = USCAL2 * SCALEM
  1749. WORK(2) = USCAL1
  1750. IF ( ERREST ) WORK(3) = SCONDA
  1751. IF ( LSVEC .AND. RSVEC ) THEN
  1752. WORK(4) = CONDR1
  1753. WORK(5) = CONDR2
  1754. END IF
  1755. IF ( L2TRAN ) THEN
  1756. WORK(6) = ENTRA
  1757. WORK(7) = ENTRAT
  1758. END IF
  1759. *
  1760. IWORK(1) = NR
  1761. IWORK(2) = NUMRANK
  1762. IWORK(3) = WARNING
  1763. *
  1764. RETURN
  1765. * ..
  1766. * .. END OF SGEJSV
  1767. * ..
  1768. END
  1769. *