dgesvj.c 54 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797
  1. /* dgesvj.f -- translated by f2c (version 20061008).
  2. You must link the resulting object file with libf2c:
  3. on Microsoft Windows system, link with libf2c.lib;
  4. on Linux or Unix systems, link with .../path/to/libf2c.a -lm
  5. or, if you install libf2c.a in a standard place, with -lf2c -lm
  6. -- in that order, at the end of the command line, as in
  7. cc *.o -lf2c -lm
  8. Source for libf2c is in /netlib/f2c/libf2c.zip, e.g.,
  9. http://www.netlib.org/f2c/libf2c.zip
  10. */
  11. #include "f2c.h"
  12. #include "blaswrap.h"
  13. /* Table of constant values */
  14. static doublereal c_b17 = 0.;
  15. static doublereal c_b18 = 1.;
  16. static integer c__1 = 1;
  17. static integer c__0 = 0;
  18. static integer c__2 = 2;
  19. /* Subroutine */ int dgesvj_(char *joba, char *jobu, char *jobv, integer *m,
  20. integer *n, doublereal *a, integer *lda, doublereal *sva, integer *mv,
  21. doublereal *v, integer *ldv, doublereal *work, integer *lwork,
  22. integer *info)
  23. {
  24. /* System generated locals */
  25. integer a_dim1, a_offset, v_dim1, v_offset, i__1, i__2, i__3, i__4, i__5;
  26. doublereal d__1, d__2;
  27. /* Builtin functions */
  28. double sqrt(doublereal), d_sign(doublereal *, doublereal *);
  29. /* Local variables */
  30. doublereal bigtheta;
  31. integer pskipped, i__, p, q;
  32. doublereal t;
  33. integer n2, n4;
  34. doublereal rootsfmin;
  35. integer n34;
  36. doublereal cs, sn;
  37. integer ir1, jbc;
  38. doublereal big;
  39. integer kbl, igl, ibr, jgl, nbl;
  40. doublereal tol;
  41. integer mvl;
  42. doublereal aapp, aapq, aaqq;
  43. extern doublereal ddot_(integer *, doublereal *, integer *, doublereal *,
  44. integer *);
  45. doublereal ctol;
  46. integer ierr;
  47. doublereal aapp0;
  48. extern doublereal dnrm2_(integer *, doublereal *, integer *);
  49. doublereal temp1;
  50. extern /* Subroutine */ int dscal_(integer *, doublereal *, doublereal *,
  51. integer *);
  52. doublereal scale, large, apoaq, aqoap;
  53. extern logical lsame_(char *, char *);
  54. doublereal theta, small, sfmin;
  55. logical lsvec;
  56. extern /* Subroutine */ int dcopy_(integer *, doublereal *, integer *,
  57. doublereal *, integer *);
  58. doublereal fastr[5];
  59. extern /* Subroutine */ int dswap_(integer *, doublereal *, integer *,
  60. doublereal *, integer *);
  61. logical applv, rsvec;
  62. extern /* Subroutine */ int daxpy_(integer *, doublereal *, doublereal *,
  63. integer *, doublereal *, integer *);
  64. logical uctol;
  65. extern /* Subroutine */ int drotm_(integer *, doublereal *, integer *,
  66. doublereal *, integer *, doublereal *);
  67. logical lower, upper, rotok;
  68. extern /* Subroutine */ int dgsvj0_(char *, integer *, integer *,
  69. doublereal *, integer *, doublereal *, doublereal *, integer *,
  70. doublereal *, integer *, doublereal *, doublereal *, doublereal *,
  71. integer *, doublereal *, integer *, integer *), dgsvj1_(
  72. char *, integer *, integer *, integer *, doublereal *, integer *,
  73. doublereal *, doublereal *, integer *, doublereal *, integer *,
  74. doublereal *, doublereal *, doublereal *, integer *, doublereal *,
  75. integer *, integer *);
  76. extern doublereal dlamch_(char *);
  77. extern /* Subroutine */ int dlascl_(char *, integer *, integer *,
  78. doublereal *, doublereal *, integer *, integer *, doublereal *,
  79. integer *, integer *);
  80. extern integer idamax_(integer *, doublereal *, integer *);
  81. extern /* Subroutine */ int dlaset_(char *, integer *, integer *,
  82. doublereal *, doublereal *, doublereal *, integer *),
  83. xerbla_(char *, integer *);
  84. integer ijblsk, swband, blskip;
  85. doublereal mxaapq;
  86. extern /* Subroutine */ int dlassq_(integer *, doublereal *, integer *,
  87. doublereal *, doublereal *);
  88. doublereal thsign, mxsinj;
  89. integer emptsw, notrot, iswrot, lkahead;
  90. logical goscale, noscale;
  91. doublereal rootbig, epsilon, rooteps;
  92. integer rowskip;
  93. doublereal roottol;
  94. /* -- LAPACK routine (version 3.2) -- */
  95. /* -- Contributed by Zlatko Drmac of the University of Zagreb and -- */
  96. /* -- Kresimir Veselic of the Fernuniversitaet Hagen -- */
  97. /* -- November 2008 -- */
  98. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  99. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  100. /* This routine is also part of SIGMA (version 1.23, October 23. 2008.) */
  101. /* SIGMA is a library of algorithms for highly accurate algorithms for */
  102. /* computation of SVD, PSVD, QSVD, (H,K)-SVD, and for solution of the */
  103. /* eigenvalue problems Hx = lambda M x, H M x = lambda x with H, M > 0. */
  104. /* -#- Scalar Arguments -#- */
  105. /* -#- Array Arguments -#- */
  106. /* .. */
  107. /* Purpose */
  108. /* ~~~~~~~ */
  109. /* DGESVJ computes the singular value decomposition (SVD) of a real */
  110. /* M-by-N matrix A, where M >= N. The SVD of A is written as */
  111. /* [++] [xx] [x0] [xx] */
  112. /* A = U * SIGMA * V^t, [++] = [xx] * [ox] * [xx] */
  113. /* [++] [xx] */
  114. /* where SIGMA is an N-by-N diagonal matrix, U is an M-by-N orthonormal */
  115. /* matrix, and V is an N-by-N orthogonal matrix. The diagonal elements */
  116. /* of SIGMA are the singular values of A. The columns of U and V are the */
  117. /* left and the right singular vectors of A, respectively. */
  118. /* Further Details */
  119. /* ~~~~~~~~~~~~~~~ */
  120. /* The orthogonal N-by-N matrix V is obtained as a product of Jacobi plane */
  121. /* rotations. The rotations are implemented as fast scaled rotations of */
  122. /* Anda and Park [1]. In the case of underflow of the Jacobi angle, a */
  123. /* modified Jacobi transformation of Drmac [4] is used. Pivot strategy uses */
  124. /* column interchanges of de Rijk [2]. The relative accuracy of the computed */
  125. /* singular values and the accuracy of the computed singular vectors (in */
  126. /* angle metric) is as guaranteed by the theory of Demmel and Veselic [3]. */
  127. /* The condition number that determines the accuracy in the full rank case */
  128. /* is essentially min_{D=diag} kappa(A*D), where kappa(.) is the */
  129. /* spectral condition number. The best performance of this Jacobi SVD */
  130. /* procedure is achieved if used in an accelerated version of Drmac and */
  131. /* Veselic [5,6], and it is the kernel routine in the SIGMA library [7]. */
  132. /* Some tunning parameters (marked with [TP]) are available for the */
  133. /* implementer. */
  134. /* The computational range for the nonzero singular values is the machine */
  135. /* number interval ( UNDERFLOW , OVERFLOW ). In extreme cases, even */
  136. /* denormalized singular values can be computed with the corresponding */
  137. /* gradual loss of accurate digits. */
  138. /* Contributors */
  139. /* ~~~~~~~~~~~~ */
  140. /* Zlatko Drmac (Zagreb, Croatia) and Kresimir Veselic (Hagen, Germany) */
  141. /* References */
  142. /* ~~~~~~~~~~ */
  143. /* [1] A. A. Anda and H. Park: Fast plane rotations with dynamic scaling. */
  144. /* SIAM J. matrix Anal. Appl., Vol. 15 (1994), pp. 162-174. */
  145. /* [2] P. P. M. De Rijk: A one-sided Jacobi algorithm for computing the */
  146. /* singular value decomposition on a vector computer. */
  147. /* SIAM J. Sci. Stat. Comp., Vol. 10 (1998), pp. 359-371. */
  148. /* [3] J. Demmel and K. Veselic: Jacobi method is more accurate than QR. */
  149. /* [4] Z. Drmac: Implementation of Jacobi rotations for accurate singular */
  150. /* value computation in floating point arithmetic. */
  151. /* SIAM J. Sci. Comp., Vol. 18 (1997), pp. 1200-1222. */
  152. /* [5] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm I. */
  153. /* SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1322-1342. */
  154. /* LAPACK Working note 169. */
  155. /* [6] Z. Drmac and K. Veselic: New fast and accurate Jacobi SVD algorithm II. */
  156. /* SIAM J. Matrix Anal. Appl. Vol. 35, No. 2 (2008), pp. 1343-1362. */
  157. /* LAPACK Working note 170. */
  158. /* [7] Z. Drmac: SIGMA - mathematical software library for accurate SVD, PSV, */
  159. /* QSVD, (H,K)-SVD computations. */
  160. /* Department of Mathematics, University of Zagreb, 2008. */
  161. /* Bugs, Examples and Comments */
  162. /* ~~~~~~~~~~~~~~~~~~~~~~~~~~~ */
  163. /* Please report all bugs and send interesting test examples and comments to */
  164. /* drmac@math.hr. Thank you. */
  165. /* Arguments */
  166. /* ~~~~~~~~~ */
  167. /* JOBA (input) CHARACTER* 1 */
  168. /* Specifies the structure of A. */
  169. /* = 'L': The input matrix A is lower triangular; */
  170. /* = 'U': The input matrix A is upper triangular; */
  171. /* = 'G': The input matrix A is general M-by-N matrix, M >= N. */
  172. /* JOBU (input) CHARACTER*1 */
  173. /* Specifies whether to compute the left singular vectors */
  174. /* (columns of U): */
  175. /* = 'U': The left singular vectors corresponding to the nonzero */
  176. /* singular values are computed and returned in the leading */
  177. /* columns of A. See more details in the description of A. */
  178. /* The default numerical orthogonality threshold is set to */
  179. /* approximately TOL=CTOL*EPS, CTOL=DSQRT(M), EPS=DLAMCH('E'). */
  180. /* = 'C': Analogous to JOBU='U', except that user can control the */
  181. /* level of numerical orthogonality of the computed left */
  182. /* singular vectors. TOL can be set to TOL = CTOL*EPS, where */
  183. /* CTOL is given on input in the array WORK. */
  184. /* No CTOL smaller than ONE is allowed. CTOL greater */
  185. /* than 1 / EPS is meaningless. The option 'C' */
  186. /* can be used if M*EPS is satisfactory orthogonality */
  187. /* of the computed left singular vectors, so CTOL=M could */
  188. /* save few sweeps of Jacobi rotations. */
  189. /* See the descriptions of A and WORK(1). */
  190. /* = 'N': The matrix U is not computed. However, see the */
  191. /* description of A. */
  192. /* JOBV (input) CHARACTER*1 */
  193. /* Specifies whether to compute the right singular vectors, that */
  194. /* is, the matrix V: */
  195. /* = 'V' : the matrix V is computed and returned in the array V */
  196. /* = 'A' : the Jacobi rotations are applied to the MV-by-N */
  197. /* array V. In other words, the right singular vector */
  198. /* matrix V is not computed explicitly, instead it is */
  199. /* applied to an MV-by-N matrix initially stored in the */
  200. /* first MV rows of V. */
  201. /* = 'N' : the matrix V is not computed and the array V is not */
  202. /* referenced */
  203. /* M (input) INTEGER */
  204. /* The number of rows of the input matrix A. M >= 0. */
  205. /* N (input) INTEGER */
  206. /* The number of columns of the input matrix A. */
  207. /* M >= N >= 0. */
  208. /* A (input/output) REAL array, dimension (LDA,N) */
  209. /* On entry, the M-by-N matrix A. */
  210. /* On exit, */
  211. /* If JOBU .EQ. 'U' .OR. JOBU .EQ. 'C': */
  212. /* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ */
  213. /* If INFO .EQ. 0, */
  214. /* ~~~~~~~~~~~~~~~ */
  215. /* RANKA orthonormal columns of U are returned in the */
  216. /* leading RANKA columns of the array A. Here RANKA <= N */
  217. /* is the number of computed singular values of A that are */
  218. /* above the underflow threshold DLAMCH('S'). The singular */
  219. /* vectors corresponding to underflowed or zero singular */
  220. /* values are not computed. The value of RANKA is returned */
  221. /* in the array WORK as RANKA=NINT(WORK(2)). Also see the */
  222. /* descriptions of SVA and WORK. The computed columns of U */
  223. /* are mutually numerically orthogonal up to approximately */
  224. /* TOL=DSQRT(M)*EPS (default); or TOL=CTOL*EPS (JOBU.EQ.'C'), */
  225. /* see the description of JOBU. */
  226. /* If INFO .GT. 0, */
  227. /* ~~~~~~~~~~~~~~~ */
  228. /* the procedure DGESVJ did not converge in the given number */
  229. /* of iterations (sweeps). In that case, the computed */
  230. /* columns of U may not be orthogonal up to TOL. The output */
  231. /* U (stored in A), SIGMA (given by the computed singular */
  232. /* values in SVA(1:N)) and V is still a decomposition of the */
  233. /* input matrix A in the sense that the residual */
  234. /* ||A-SCALE*U*SIGMA*V^T||_2 / ||A||_2 is small. */
  235. /* If JOBU .EQ. 'N': */
  236. /* ~~~~~~~~~~~~~~~~~ */
  237. /* If INFO .EQ. 0 */
  238. /* ~~~~~~~~~~~~~~ */
  239. /* Note that the left singular vectors are 'for free' in the */
  240. /* one-sided Jacobi SVD algorithm. However, if only the */
  241. /* singular values are needed, the level of numerical */
  242. /* orthogonality of U is not an issue and iterations are */
  243. /* stopped when the columns of the iterated matrix are */
  244. /* numerically orthogonal up to approximately M*EPS. Thus, */
  245. /* on exit, A contains the columns of U scaled with the */
  246. /* corresponding singular values. */
  247. /* If INFO .GT. 0, */
  248. /* ~~~~~~~~~~~~~~~ */
  249. /* the procedure DGESVJ did not converge in the given number */
  250. /* of iterations (sweeps). */
  251. /* LDA (input) INTEGER */
  252. /* The leading dimension of the array A. LDA >= max(1,M). */
  253. /* SVA (workspace/output) REAL array, dimension (N) */
  254. /* On exit, */
  255. /* If INFO .EQ. 0, */
  256. /* ~~~~~~~~~~~~~~~ */
  257. /* depending on the value SCALE = WORK(1), we have: */
  258. /* If SCALE .EQ. ONE: */
  259. /* ~~~~~~~~~~~~~~~~~~ */
  260. /* SVA(1:N) contains the computed singular values of A. */
  261. /* During the computation SVA contains the Euclidean column */
  262. /* norms of the iterated matrices in the array A. */
  263. /* If SCALE .NE. ONE: */
  264. /* ~~~~~~~~~~~~~~~~~~ */
  265. /* The singular values of A are SCALE*SVA(1:N), and this */
  266. /* factored representation is due to the fact that some of the */
  267. /* singular values of A might underflow or overflow. */
  268. /* If INFO .GT. 0, */
  269. /* ~~~~~~~~~~~~~~~ */
  270. /* the procedure DGESVJ did not converge in the given number of */
  271. /* iterations (sweeps) and SCALE*SVA(1:N) may not be accurate. */
  272. /* MV (input) INTEGER */
  273. /* If JOBV .EQ. 'A', then the product of Jacobi rotations in DGESVJ */
  274. /* is applied to the first MV rows of V. See the description of JOBV. */
  275. /* V (input/output) REAL array, dimension (LDV,N) */
  276. /* If JOBV = 'V', then V contains on exit the N-by-N matrix of */
  277. /* the right singular vectors; */
  278. /* If JOBV = 'A', then V contains the product of the computed right */
  279. /* singular vector matrix and the initial matrix in */
  280. /* the array V. */
  281. /* If JOBV = 'N', then V is not referenced. */
  282. /* LDV (input) INTEGER */
  283. /* The leading dimension of the array V, LDV .GE. 1. */
  284. /* If JOBV .EQ. 'V', then LDV .GE. max(1,N). */
  285. /* If JOBV .EQ. 'A', then LDV .GE. max(1,MV) . */
  286. /* WORK (input/workspace/output) REAL array, dimension max(4,M+N). */
  287. /* On entry, */
  288. /* If JOBU .EQ. 'C', */
  289. /* ~~~~~~~~~~~~~~~~~ */
  290. /* WORK(1) = CTOL, where CTOL defines the threshold for convergence. */
  291. /* The process stops if all columns of A are mutually */
  292. /* orthogonal up to CTOL*EPS, EPS=DLAMCH('E'). */
  293. /* It is required that CTOL >= ONE, i.e. it is not */
  294. /* allowed to force the routine to obtain orthogonality */
  295. /* below EPSILON. */
  296. /* On exit, */
  297. /* WORK(1) = SCALE is the scaling factor such that SCALE*SVA(1:N) */
  298. /* are the computed singular vcalues of A. */
  299. /* (See description of SVA().) */
  300. /* WORK(2) = NINT(WORK(2)) is the number of the computed nonzero */
  301. /* singular values. */
  302. /* WORK(3) = NINT(WORK(3)) is the number of the computed singular */
  303. /* values that are larger than the underflow threshold. */
  304. /* WORK(4) = NINT(WORK(4)) is the number of sweeps of Jacobi */
  305. /* rotations needed for numerical convergence. */
  306. /* WORK(5) = max_{i.NE.j} |COS(A(:,i),A(:,j))| in the last sweep. */
  307. /* This is useful information in cases when DGESVJ did */
  308. /* not converge, as it can be used to estimate whether */
  309. /* the output is stil useful and for post festum analysis. */
  310. /* WORK(6) = the largest absolute value over all sines of the */
  311. /* Jacobi rotation angles in the last sweep. It can be */
  312. /* useful for a post festum analysis. */
  313. /* LWORK length of WORK, WORK >= MAX(6,M+N) */
  314. /* INFO (output) INTEGER */
  315. /* = 0 : successful exit. */
  316. /* < 0 : if INFO = -i, then the i-th argument had an illegal value */
  317. /* > 0 : DGESVJ did not converge in the maximal allowed number (30) */
  318. /* of sweeps. The output may still be useful. See the */
  319. /* description of WORK. */
  320. /* Local Parameters */
  321. /* Local Scalars */
  322. /* Local Arrays */
  323. /* Intrinsic Functions */
  324. /* External Functions */
  325. /* .. from BLAS */
  326. /* .. from LAPACK */
  327. /* External Subroutines */
  328. /* .. from BLAS */
  329. /* .. from LAPACK */
  330. /* Test the input arguments */
  331. /* Parameter adjustments */
  332. --sva;
  333. a_dim1 = *lda;
  334. a_offset = 1 + a_dim1;
  335. a -= a_offset;
  336. v_dim1 = *ldv;
  337. v_offset = 1 + v_dim1;
  338. v -= v_offset;
  339. --work;
  340. /* Function Body */
  341. lsvec = lsame_(jobu, "U");
  342. uctol = lsame_(jobu, "C");
  343. rsvec = lsame_(jobv, "V");
  344. applv = lsame_(jobv, "A");
  345. upper = lsame_(joba, "U");
  346. lower = lsame_(joba, "L");
  347. if (! (upper || lower || lsame_(joba, "G"))) {
  348. *info = -1;
  349. } else if (! (lsvec || uctol || lsame_(jobu, "N")))
  350. {
  351. *info = -2;
  352. } else if (! (rsvec || applv || lsame_(jobv, "N")))
  353. {
  354. *info = -3;
  355. } else if (*m < 0) {
  356. *info = -4;
  357. } else if (*n < 0 || *n > *m) {
  358. *info = -5;
  359. } else if (*lda < *m) {
  360. *info = -7;
  361. } else if (*mv < 0) {
  362. *info = -9;
  363. } else if (rsvec && *ldv < *n || applv && *ldv < *mv) {
  364. *info = -11;
  365. } else if (uctol && work[1] <= 1.) {
  366. *info = -12;
  367. } else /* if(complicated condition) */ {
  368. /* Computing MAX */
  369. i__1 = *m + *n;
  370. if (*lwork < max(i__1,6)) {
  371. *info = -13;
  372. } else {
  373. *info = 0;
  374. }
  375. }
  376. /* #:( */
  377. if (*info != 0) {
  378. i__1 = -(*info);
  379. xerbla_("DGESVJ", &i__1);
  380. return 0;
  381. }
  382. /* #:) Quick return for void matrix */
  383. if (*m == 0 || *n == 0) {
  384. return 0;
  385. }
  386. /* Set numerical parameters */
  387. /* The stopping criterion for Jacobi rotations is */
  388. /* max_{i<>j}|A(:,i)^T * A(:,j)|/(||A(:,i)||*||A(:,j)||) < CTOL*EPS */
  389. /* where EPS is the round-off and CTOL is defined as follows: */
  390. if (uctol) {
  391. /* ... user controlled */
  392. ctol = work[1];
  393. } else {
  394. /* ... default */
  395. if (lsvec || rsvec || applv) {
  396. ctol = sqrt((doublereal) (*m));
  397. } else {
  398. ctol = (doublereal) (*m);
  399. }
  400. }
  401. /* ... and the machine dependent parameters are */
  402. /* [!] (Make sure that DLAMCH() works properly on the target machine.) */
  403. epsilon = dlamch_("Epsilon");
  404. rooteps = sqrt(epsilon);
  405. sfmin = dlamch_("SafeMinimum");
  406. rootsfmin = sqrt(sfmin);
  407. small = sfmin / epsilon;
  408. big = dlamch_("Overflow");
  409. /* BIG = ONE / SFMIN */
  410. rootbig = 1. / rootsfmin;
  411. large = big / sqrt((doublereal) (*m * *n));
  412. bigtheta = 1. / rooteps;
  413. tol = ctol * epsilon;
  414. roottol = sqrt(tol);
  415. if ((doublereal) (*m) * epsilon >= 1.) {
  416. *info = -5;
  417. i__1 = -(*info);
  418. xerbla_("DGESVJ", &i__1);
  419. return 0;
  420. }
  421. /* Initialize the right singular vector matrix. */
  422. if (rsvec) {
  423. mvl = *n;
  424. dlaset_("A", &mvl, n, &c_b17, &c_b18, &v[v_offset], ldv);
  425. } else if (applv) {
  426. mvl = *mv;
  427. }
  428. rsvec = rsvec || applv;
  429. /* Initialize SVA( 1:N ) = ( ||A e_i||_2, i = 1:N ) */
  430. /* (!) If necessary, scale A to protect the largest singular value */
  431. /* from overflow. It is possible that saving the largest singular */
  432. /* value destroys the information about the small ones. */
  433. /* This initial scaling is almost minimal in the sense that the */
  434. /* goal is to make sure that no column norm overflows, and that */
  435. /* DSQRT(N)*max_i SVA(i) does not overflow. If INFinite entries */
  436. /* in A are detected, the procedure returns with INFO=-6. */
  437. scale = 1. / sqrt((doublereal) (*m) * (doublereal) (*n));
  438. noscale = TRUE_;
  439. goscale = TRUE_;
  440. if (lower) {
  441. /* the input matrix is M-by-N lower triangular (trapezoidal) */
  442. i__1 = *n;
  443. for (p = 1; p <= i__1; ++p) {
  444. aapp = 0.;
  445. aaqq = 0.;
  446. i__2 = *m - p + 1;
  447. dlassq_(&i__2, &a[p + p * a_dim1], &c__1, &aapp, &aaqq);
  448. if (aapp > big) {
  449. *info = -6;
  450. i__2 = -(*info);
  451. xerbla_("DGESVJ", &i__2);
  452. return 0;
  453. }
  454. aaqq = sqrt(aaqq);
  455. if (aapp < big / aaqq && noscale) {
  456. sva[p] = aapp * aaqq;
  457. } else {
  458. noscale = FALSE_;
  459. sva[p] = aapp * (aaqq * scale);
  460. if (goscale) {
  461. goscale = FALSE_;
  462. i__2 = p - 1;
  463. for (q = 1; q <= i__2; ++q) {
  464. sva[q] *= scale;
  465. /* L1873: */
  466. }
  467. }
  468. }
  469. /* L1874: */
  470. }
  471. } else if (upper) {
  472. /* the input matrix is M-by-N upper triangular (trapezoidal) */
  473. i__1 = *n;
  474. for (p = 1; p <= i__1; ++p) {
  475. aapp = 0.;
  476. aaqq = 0.;
  477. dlassq_(&p, &a[p * a_dim1 + 1], &c__1, &aapp, &aaqq);
  478. if (aapp > big) {
  479. *info = -6;
  480. i__2 = -(*info);
  481. xerbla_("DGESVJ", &i__2);
  482. return 0;
  483. }
  484. aaqq = sqrt(aaqq);
  485. if (aapp < big / aaqq && noscale) {
  486. sva[p] = aapp * aaqq;
  487. } else {
  488. noscale = FALSE_;
  489. sva[p] = aapp * (aaqq * scale);
  490. if (goscale) {
  491. goscale = FALSE_;
  492. i__2 = p - 1;
  493. for (q = 1; q <= i__2; ++q) {
  494. sva[q] *= scale;
  495. /* L2873: */
  496. }
  497. }
  498. }
  499. /* L2874: */
  500. }
  501. } else {
  502. /* the input matrix is M-by-N general dense */
  503. i__1 = *n;
  504. for (p = 1; p <= i__1; ++p) {
  505. aapp = 0.;
  506. aaqq = 0.;
  507. dlassq_(m, &a[p * a_dim1 + 1], &c__1, &aapp, &aaqq);
  508. if (aapp > big) {
  509. *info = -6;
  510. i__2 = -(*info);
  511. xerbla_("DGESVJ", &i__2);
  512. return 0;
  513. }
  514. aaqq = sqrt(aaqq);
  515. if (aapp < big / aaqq && noscale) {
  516. sva[p] = aapp * aaqq;
  517. } else {
  518. noscale = FALSE_;
  519. sva[p] = aapp * (aaqq * scale);
  520. if (goscale) {
  521. goscale = FALSE_;
  522. i__2 = p - 1;
  523. for (q = 1; q <= i__2; ++q) {
  524. sva[q] *= scale;
  525. /* L3873: */
  526. }
  527. }
  528. }
  529. /* L3874: */
  530. }
  531. }
  532. if (noscale) {
  533. scale = 1.;
  534. }
  535. /* Move the smaller part of the spectrum from the underflow threshold */
  536. /* (!) Start by determining the position of the nonzero entries of the */
  537. /* array SVA() relative to ( SFMIN, BIG ). */
  538. aapp = 0.;
  539. aaqq = big;
  540. i__1 = *n;
  541. for (p = 1; p <= i__1; ++p) {
  542. if (sva[p] != 0.) {
  543. /* Computing MIN */
  544. d__1 = aaqq, d__2 = sva[p];
  545. aaqq = min(d__1,d__2);
  546. }
  547. /* Computing MAX */
  548. d__1 = aapp, d__2 = sva[p];
  549. aapp = max(d__1,d__2);
  550. /* L4781: */
  551. }
  552. /* #:) Quick return for zero matrix */
  553. if (aapp == 0.) {
  554. if (lsvec) {
  555. dlaset_("G", m, n, &c_b17, &c_b18, &a[a_offset], lda);
  556. }
  557. work[1] = 1.;
  558. work[2] = 0.;
  559. work[3] = 0.;
  560. work[4] = 0.;
  561. work[5] = 0.;
  562. work[6] = 0.;
  563. return 0;
  564. }
  565. /* #:) Quick return for one-column matrix */
  566. if (*n == 1) {
  567. if (lsvec) {
  568. dlascl_("G", &c__0, &c__0, &sva[1], &scale, m, &c__1, &a[a_dim1 +
  569. 1], lda, &ierr);
  570. }
  571. work[1] = 1. / scale;
  572. if (sva[1] >= sfmin) {
  573. work[2] = 1.;
  574. } else {
  575. work[2] = 0.;
  576. }
  577. work[3] = 0.;
  578. work[4] = 0.;
  579. work[5] = 0.;
  580. work[6] = 0.;
  581. return 0;
  582. }
  583. /* Protect small singular values from underflow, and try to */
  584. /* avoid underflows/overflows in computing Jacobi rotations. */
  585. sn = sqrt(sfmin / epsilon);
  586. temp1 = sqrt(big / (doublereal) (*n));
  587. if (aapp <= sn || aaqq >= temp1 || sn <= aaqq && aapp <= temp1) {
  588. /* Computing MIN */
  589. d__1 = big, d__2 = temp1 / aapp;
  590. temp1 = min(d__1,d__2);
  591. /* AAQQ = AAQQ*TEMP1 */
  592. /* AAPP = AAPP*TEMP1 */
  593. } else if (aaqq <= sn && aapp <= temp1) {
  594. /* Computing MIN */
  595. d__1 = sn / aaqq, d__2 = big / (aapp * sqrt((doublereal) (*n)));
  596. temp1 = min(d__1,d__2);
  597. /* AAQQ = AAQQ*TEMP1 */
  598. /* AAPP = AAPP*TEMP1 */
  599. } else if (aaqq >= sn && aapp >= temp1) {
  600. /* Computing MAX */
  601. d__1 = sn / aaqq, d__2 = temp1 / aapp;
  602. temp1 = max(d__1,d__2);
  603. /* AAQQ = AAQQ*TEMP1 */
  604. /* AAPP = AAPP*TEMP1 */
  605. } else if (aaqq <= sn && aapp >= temp1) {
  606. /* Computing MIN */
  607. d__1 = sn / aaqq, d__2 = big / (sqrt((doublereal) (*n)) * aapp);
  608. temp1 = min(d__1,d__2);
  609. /* AAQQ = AAQQ*TEMP1 */
  610. /* AAPP = AAPP*TEMP1 */
  611. } else {
  612. temp1 = 1.;
  613. }
  614. /* Scale, if necessary */
  615. if (temp1 != 1.) {
  616. dlascl_("G", &c__0, &c__0, &c_b18, &temp1, n, &c__1, &sva[1], n, &
  617. ierr);
  618. }
  619. scale = temp1 * scale;
  620. if (scale != 1.) {
  621. dlascl_(joba, &c__0, &c__0, &c_b18, &scale, m, n, &a[a_offset], lda, &
  622. ierr);
  623. scale = 1. / scale;
  624. }
  625. /* Row-cyclic Jacobi SVD algorithm with column pivoting */
  626. emptsw = *n * (*n - 1) / 2;
  627. notrot = 0;
  628. fastr[0] = 0.;
  629. /* A is represented in factored form A = A * diag(WORK), where diag(WORK) */
  630. /* is initialized to identity. WORK is updated during fast scaled */
  631. /* rotations. */
  632. i__1 = *n;
  633. for (q = 1; q <= i__1; ++q) {
  634. work[q] = 1.;
  635. /* L1868: */
  636. }
  637. swband = 3;
  638. /* [TP] SWBAND is a tuning parameter [TP]. It is meaningful and effective */
  639. /* if DGESVJ is used as a computational routine in the preconditioned */
  640. /* Jacobi SVD algorithm DGESVJ. For sweeps i=1:SWBAND the procedure */
  641. /* works on pivots inside a band-like region around the diagonal. */
  642. /* The boundaries are determined dynamically, based on the number of */
  643. /* pivots above a threshold. */
  644. kbl = min(8,*n);
  645. /* [TP] KBL is a tuning parameter that defines the tile size in the */
  646. /* tiling of the p-q loops of pivot pairs. In general, an optimal */
  647. /* value of KBL depends on the matrix dimensions and on the */
  648. /* parameters of the computer's memory. */
  649. nbl = *n / kbl;
  650. if (nbl * kbl != *n) {
  651. ++nbl;
  652. }
  653. /* Computing 2nd power */
  654. i__1 = kbl;
  655. blskip = i__1 * i__1;
  656. /* [TP] BLKSKIP is a tuning parameter that depends on SWBAND and KBL. */
  657. rowskip = min(5,kbl);
  658. /* [TP] ROWSKIP is a tuning parameter. */
  659. lkahead = 1;
  660. /* [TP] LKAHEAD is a tuning parameter. */
  661. /* Quasi block transformations, using the lower (upper) triangular */
  662. /* structure of the input matrix. The quasi-block-cycling usually */
  663. /* invokes cubic convergence. Big part of this cycle is done inside */
  664. /* canonical subspaces of dimensions less than M. */
  665. /* Computing MAX */
  666. i__1 = 64, i__2 = kbl << 2;
  667. if ((lower || upper) && *n > max(i__1,i__2)) {
  668. /* [TP] The number of partition levels and the actual partition are */
  669. /* tuning parameters. */
  670. n4 = *n / 4;
  671. n2 = *n / 2;
  672. n34 = n4 * 3;
  673. if (applv) {
  674. q = 0;
  675. } else {
  676. q = 1;
  677. }
  678. if (lower) {
  679. /* This works very well on lower triangular matrices, in particular */
  680. /* in the framework of the preconditioned Jacobi SVD (xGEJSV). */
  681. /* The idea is simple: */
  682. /* [+ 0 0 0] Note that Jacobi transformations of [0 0] */
  683. /* [+ + 0 0] [0 0] */
  684. /* [+ + x 0] actually work on [x 0] [x 0] */
  685. /* [+ + x x] [x x]. [x x] */
  686. i__1 = *m - n34;
  687. i__2 = *n - n34;
  688. i__3 = *lwork - *n;
  689. dgsvj0_(jobv, &i__1, &i__2, &a[n34 + 1 + (n34 + 1) * a_dim1], lda,
  690. &work[n34 + 1], &sva[n34 + 1], &mvl, &v[n34 * q + 1 + (
  691. n34 + 1) * v_dim1], ldv, &epsilon, &sfmin, &tol, &c__2, &
  692. work[*n + 1], &i__3, &ierr);
  693. i__1 = *m - n2;
  694. i__2 = n34 - n2;
  695. i__3 = *lwork - *n;
  696. dgsvj0_(jobv, &i__1, &i__2, &a[n2 + 1 + (n2 + 1) * a_dim1], lda, &
  697. work[n2 + 1], &sva[n2 + 1], &mvl, &v[n2 * q + 1 + (n2 + 1)
  698. * v_dim1], ldv, &epsilon, &sfmin, &tol, &c__2, &work[*n
  699. + 1], &i__3, &ierr);
  700. i__1 = *m - n2;
  701. i__2 = *n - n2;
  702. i__3 = *lwork - *n;
  703. dgsvj1_(jobv, &i__1, &i__2, &n4, &a[n2 + 1 + (n2 + 1) * a_dim1],
  704. lda, &work[n2 + 1], &sva[n2 + 1], &mvl, &v[n2 * q + 1 + (
  705. n2 + 1) * v_dim1], ldv, &epsilon, &sfmin, &tol, &c__1, &
  706. work[*n + 1], &i__3, &ierr);
  707. i__1 = *m - n4;
  708. i__2 = n2 - n4;
  709. i__3 = *lwork - *n;
  710. dgsvj0_(jobv, &i__1, &i__2, &a[n4 + 1 + (n4 + 1) * a_dim1], lda, &
  711. work[n4 + 1], &sva[n4 + 1], &mvl, &v[n4 * q + 1 + (n4 + 1)
  712. * v_dim1], ldv, &epsilon, &sfmin, &tol, &c__1, &work[*n
  713. + 1], &i__3, &ierr);
  714. i__1 = *lwork - *n;
  715. dgsvj0_(jobv, m, &n4, &a[a_offset], lda, &work[1], &sva[1], &mvl,
  716. &v[v_offset], ldv, &epsilon, &sfmin, &tol, &c__1, &work[*
  717. n + 1], &i__1, &ierr);
  718. i__1 = *lwork - *n;
  719. dgsvj1_(jobv, m, &n2, &n4, &a[a_offset], lda, &work[1], &sva[1], &
  720. mvl, &v[v_offset], ldv, &epsilon, &sfmin, &tol, &c__1, &
  721. work[*n + 1], &i__1, &ierr);
  722. } else if (upper) {
  723. i__1 = *lwork - *n;
  724. dgsvj0_(jobv, &n4, &n4, &a[a_offset], lda, &work[1], &sva[1], &
  725. mvl, &v[v_offset], ldv, &epsilon, &sfmin, &tol, &c__2, &
  726. work[*n + 1], &i__1, &ierr);
  727. i__1 = *lwork - *n;
  728. dgsvj0_(jobv, &n2, &n4, &a[(n4 + 1) * a_dim1 + 1], lda, &work[n4
  729. + 1], &sva[n4 + 1], &mvl, &v[n4 * q + 1 + (n4 + 1) *
  730. v_dim1], ldv, &epsilon, &sfmin, &tol, &c__1, &work[*n + 1]
  731. , &i__1, &ierr);
  732. i__1 = *lwork - *n;
  733. dgsvj1_(jobv, &n2, &n2, &n4, &a[a_offset], lda, &work[1], &sva[1],
  734. &mvl, &v[v_offset], ldv, &epsilon, &sfmin, &tol, &c__1, &
  735. work[*n + 1], &i__1, &ierr);
  736. i__1 = n2 + n4;
  737. i__2 = *lwork - *n;
  738. dgsvj0_(jobv, &i__1, &n4, &a[(n2 + 1) * a_dim1 + 1], lda, &work[
  739. n2 + 1], &sva[n2 + 1], &mvl, &v[n2 * q + 1 + (n2 + 1) *
  740. v_dim1], ldv, &epsilon, &sfmin, &tol, &c__1, &work[*n + 1]
  741. , &i__2, &ierr);
  742. }
  743. }
  744. /* -#- Row-cyclic pivot strategy with de Rijk's pivoting -#- */
  745. for (i__ = 1; i__ <= 30; ++i__) {
  746. /* .. go go go ... */
  747. mxaapq = 0.;
  748. mxsinj = 0.;
  749. iswrot = 0;
  750. notrot = 0;
  751. pskipped = 0;
  752. /* Each sweep is unrolled using KBL-by-KBL tiles over the pivot pairs */
  753. /* 1 <= p < q <= N. This is the first step toward a blocked implementation */
  754. /* of the rotations. New implementation, based on block transformations, */
  755. /* is under development. */
  756. i__1 = nbl;
  757. for (ibr = 1; ibr <= i__1; ++ibr) {
  758. igl = (ibr - 1) * kbl + 1;
  759. /* Computing MIN */
  760. i__3 = lkahead, i__4 = nbl - ibr;
  761. i__2 = min(i__3,i__4);
  762. for (ir1 = 0; ir1 <= i__2; ++ir1) {
  763. igl += ir1 * kbl;
  764. /* Computing MIN */
  765. i__4 = igl + kbl - 1, i__5 = *n - 1;
  766. i__3 = min(i__4,i__5);
  767. for (p = igl; p <= i__3; ++p) {
  768. /* .. de Rijk's pivoting */
  769. i__4 = *n - p + 1;
  770. q = idamax_(&i__4, &sva[p], &c__1) + p - 1;
  771. if (p != q) {
  772. dswap_(m, &a[p * a_dim1 + 1], &c__1, &a[q * a_dim1 +
  773. 1], &c__1);
  774. if (rsvec) {
  775. dswap_(&mvl, &v[p * v_dim1 + 1], &c__1, &v[q *
  776. v_dim1 + 1], &c__1);
  777. }
  778. temp1 = sva[p];
  779. sva[p] = sva[q];
  780. sva[q] = temp1;
  781. temp1 = work[p];
  782. work[p] = work[q];
  783. work[q] = temp1;
  784. }
  785. if (ir1 == 0) {
  786. /* Column norms are periodically updated by explicit */
  787. /* norm computation. */
  788. /* Caveat: */
  789. /* Unfortunately, some BLAS implementations compute DNRM2(M,A(1,p),1) */
  790. /* as DSQRT(DDOT(M,A(1,p),1,A(1,p),1)), which may cause the result to */
  791. /* overflow for ||A(:,p)||_2 > DSQRT(overflow_threshold), and to */
  792. /* underflow for ||A(:,p)||_2 < DSQRT(underflow_threshold). */
  793. /* Hence, DNRM2 cannot be trusted, not even in the case when */
  794. /* the true norm is far from the under(over)flow boundaries. */
  795. /* If properly implemented DNRM2 is available, the IF-THEN-ELSE */
  796. /* below should read "AAPP = DNRM2( M, A(1,p), 1 ) * WORK(p)". */
  797. if (sva[p] < rootbig && sva[p] > rootsfmin) {
  798. sva[p] = dnrm2_(m, &a[p * a_dim1 + 1], &c__1) *
  799. work[p];
  800. } else {
  801. temp1 = 0.;
  802. aapp = 0.;
  803. dlassq_(m, &a[p * a_dim1 + 1], &c__1, &temp1, &
  804. aapp);
  805. sva[p] = temp1 * sqrt(aapp) * work[p];
  806. }
  807. aapp = sva[p];
  808. } else {
  809. aapp = sva[p];
  810. }
  811. if (aapp > 0.) {
  812. pskipped = 0;
  813. /* Computing MIN */
  814. i__5 = igl + kbl - 1;
  815. i__4 = min(i__5,*n);
  816. for (q = p + 1; q <= i__4; ++q) {
  817. aaqq = sva[q];
  818. if (aaqq > 0.) {
  819. aapp0 = aapp;
  820. if (aaqq >= 1.) {
  821. rotok = small * aapp <= aaqq;
  822. if (aapp < big / aaqq) {
  823. aapq = ddot_(m, &a[p * a_dim1 + 1], &
  824. c__1, &a[q * a_dim1 + 1], &
  825. c__1) * work[p] * work[q] /
  826. aaqq / aapp;
  827. } else {
  828. dcopy_(m, &a[p * a_dim1 + 1], &c__1, &
  829. work[*n + 1], &c__1);
  830. dlascl_("G", &c__0, &c__0, &aapp, &
  831. work[p], m, &c__1, &work[*n +
  832. 1], lda, &ierr);
  833. aapq = ddot_(m, &work[*n + 1], &c__1,
  834. &a[q * a_dim1 + 1], &c__1) *
  835. work[q] / aaqq;
  836. }
  837. } else {
  838. rotok = aapp <= aaqq / small;
  839. if (aapp > small / aaqq) {
  840. aapq = ddot_(m, &a[p * a_dim1 + 1], &
  841. c__1, &a[q * a_dim1 + 1], &
  842. c__1) * work[p] * work[q] /
  843. aaqq / aapp;
  844. } else {
  845. dcopy_(m, &a[q * a_dim1 + 1], &c__1, &
  846. work[*n + 1], &c__1);
  847. dlascl_("G", &c__0, &c__0, &aaqq, &
  848. work[q], m, &c__1, &work[*n +
  849. 1], lda, &ierr);
  850. aapq = ddot_(m, &work[*n + 1], &c__1,
  851. &a[p * a_dim1 + 1], &c__1) *
  852. work[p] / aapp;
  853. }
  854. }
  855. /* Computing MAX */
  856. d__1 = mxaapq, d__2 = abs(aapq);
  857. mxaapq = max(d__1,d__2);
  858. /* TO rotate or NOT to rotate, THAT is the question ... */
  859. if (abs(aapq) > tol) {
  860. /* .. rotate */
  861. /* [RTD] ROTATED = ROTATED + ONE */
  862. if (ir1 == 0) {
  863. notrot = 0;
  864. pskipped = 0;
  865. ++iswrot;
  866. }
  867. if (rotok) {
  868. aqoap = aaqq / aapp;
  869. apoaq = aapp / aaqq;
  870. theta = (d__1 = aqoap - apoaq, abs(
  871. d__1)) * -.5 / aapq;
  872. if (abs(theta) > bigtheta) {
  873. t = .5 / theta;
  874. fastr[2] = t * work[p] / work[q];
  875. fastr[3] = -t * work[q] / work[p];
  876. drotm_(m, &a[p * a_dim1 + 1], &
  877. c__1, &a[q * a_dim1 + 1],
  878. &c__1, fastr);
  879. if (rsvec) {
  880. drotm_(&mvl, &v[p * v_dim1 + 1], &c__1, &v[q *
  881. v_dim1 + 1], &c__1, fastr);
  882. }
  883. /* Computing MAX */
  884. d__1 = 0., d__2 = t * apoaq *
  885. aapq + 1.;
  886. sva[q] = aaqq * sqrt((max(d__1,
  887. d__2)));
  888. aapp *= sqrt(1. - t * aqoap *
  889. aapq);
  890. /* Computing MAX */
  891. d__1 = mxsinj, d__2 = abs(t);
  892. mxsinj = max(d__1,d__2);
  893. } else {
  894. /* .. choose correct signum for THETA and rotate */
  895. thsign = -d_sign(&c_b18, &aapq);
  896. t = 1. / (theta + thsign * sqrt(
  897. theta * theta + 1.));
  898. cs = sqrt(1. / (t * t + 1.));
  899. sn = t * cs;
  900. /* Computing MAX */
  901. d__1 = mxsinj, d__2 = abs(sn);
  902. mxsinj = max(d__1,d__2);
  903. /* Computing MAX */
  904. d__1 = 0., d__2 = t * apoaq *
  905. aapq + 1.;
  906. sva[q] = aaqq * sqrt((max(d__1,
  907. d__2)));
  908. /* Computing MAX */
  909. d__1 = 0., d__2 = 1. - t * aqoap *
  910. aapq;
  911. aapp *= sqrt((max(d__1,d__2)));
  912. apoaq = work[p] / work[q];
  913. aqoap = work[q] / work[p];
  914. if (work[p] >= 1.) {
  915. if (work[q] >= 1.) {
  916. fastr[2] = t * apoaq;
  917. fastr[3] = -t * aqoap;
  918. work[p] *= cs;
  919. work[q] *= cs;
  920. drotm_(m, &a[p * a_dim1 + 1], &c__1, &a[q *
  921. a_dim1 + 1], &c__1, fastr);
  922. if (rsvec) {
  923. drotm_(&mvl, &v[p * v_dim1 + 1], &c__1, &v[
  924. q * v_dim1 + 1], &c__1, fastr);
  925. }
  926. } else {
  927. d__1 = -t * aqoap;
  928. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1, &a[
  929. p * a_dim1 + 1], &c__1);
  930. d__1 = cs * sn * apoaq;
  931. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1, &a[
  932. q * a_dim1 + 1], &c__1);
  933. work[p] *= cs;
  934. work[q] /= cs;
  935. if (rsvec) {
  936. d__1 = -t * aqoap;
  937. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1], &
  938. c__1, &v[p * v_dim1 + 1], &c__1);
  939. d__1 = cs * sn * apoaq;
  940. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1], &
  941. c__1, &v[q * v_dim1 + 1], &c__1);
  942. }
  943. }
  944. } else {
  945. if (work[q] >= 1.) {
  946. d__1 = t * apoaq;
  947. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1, &a[
  948. q * a_dim1 + 1], &c__1);
  949. d__1 = -cs * sn * aqoap;
  950. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1, &a[
  951. p * a_dim1 + 1], &c__1);
  952. work[p] /= cs;
  953. work[q] *= cs;
  954. if (rsvec) {
  955. d__1 = t * apoaq;
  956. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1], &
  957. c__1, &v[q * v_dim1 + 1], &c__1);
  958. d__1 = -cs * sn * aqoap;
  959. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1], &
  960. c__1, &v[p * v_dim1 + 1], &c__1);
  961. }
  962. } else {
  963. if (work[p] >= work[q]) {
  964. d__1 = -t * aqoap;
  965. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1,
  966. &a[p * a_dim1 + 1], &c__1);
  967. d__1 = cs * sn * apoaq;
  968. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1,
  969. &a[q * a_dim1 + 1], &c__1);
  970. work[p] *= cs;
  971. work[q] /= cs;
  972. if (rsvec) {
  973. d__1 = -t * aqoap;
  974. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1],
  975. &c__1, &v[p * v_dim1 + 1], &
  976. c__1);
  977. d__1 = cs * sn * apoaq;
  978. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1],
  979. &c__1, &v[q * v_dim1 + 1], &
  980. c__1);
  981. }
  982. } else {
  983. d__1 = t * apoaq;
  984. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1,
  985. &a[q * a_dim1 + 1], &c__1);
  986. d__1 = -cs * sn * aqoap;
  987. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1,
  988. &a[p * a_dim1 + 1], &c__1);
  989. work[p] /= cs;
  990. work[q] *= cs;
  991. if (rsvec) {
  992. d__1 = t * apoaq;
  993. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1],
  994. &c__1, &v[q * v_dim1 + 1], &
  995. c__1);
  996. d__1 = -cs * sn * aqoap;
  997. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1],
  998. &c__1, &v[p * v_dim1 + 1], &
  999. c__1);
  1000. }
  1001. }
  1002. }
  1003. }
  1004. }
  1005. } else {
  1006. /* .. have to use modified Gram-Schmidt like transformation */
  1007. dcopy_(m, &a[p * a_dim1 + 1], &c__1, &
  1008. work[*n + 1], &c__1);
  1009. dlascl_("G", &c__0, &c__0, &aapp, &
  1010. c_b18, m, &c__1, &work[*n + 1]
  1011. , lda, &ierr);
  1012. dlascl_("G", &c__0, &c__0, &aaqq, &
  1013. c_b18, m, &c__1, &a[q *
  1014. a_dim1 + 1], lda, &ierr);
  1015. temp1 = -aapq * work[p] / work[q];
  1016. daxpy_(m, &temp1, &work[*n + 1], &
  1017. c__1, &a[q * a_dim1 + 1], &
  1018. c__1);
  1019. dlascl_("G", &c__0, &c__0, &c_b18, &
  1020. aaqq, m, &c__1, &a[q * a_dim1
  1021. + 1], lda, &ierr);
  1022. /* Computing MAX */
  1023. d__1 = 0., d__2 = 1. - aapq * aapq;
  1024. sva[q] = aaqq * sqrt((max(d__1,d__2)))
  1025. ;
  1026. mxsinj = max(mxsinj,sfmin);
  1027. }
  1028. /* END IF ROTOK THEN ... ELSE */
  1029. /* In the case of cancellation in updating SVA(q), SVA(p) */
  1030. /* recompute SVA(q), SVA(p). */
  1031. /* Computing 2nd power */
  1032. d__1 = sva[q] / aaqq;
  1033. if (d__1 * d__1 <= rooteps) {
  1034. if (aaqq < rootbig && aaqq >
  1035. rootsfmin) {
  1036. sva[q] = dnrm2_(m, &a[q * a_dim1
  1037. + 1], &c__1) * work[q];
  1038. } else {
  1039. t = 0.;
  1040. aaqq = 0.;
  1041. dlassq_(m, &a[q * a_dim1 + 1], &
  1042. c__1, &t, &aaqq);
  1043. sva[q] = t * sqrt(aaqq) * work[q];
  1044. }
  1045. }
  1046. if (aapp / aapp0 <= rooteps) {
  1047. if (aapp < rootbig && aapp >
  1048. rootsfmin) {
  1049. aapp = dnrm2_(m, &a[p * a_dim1 +
  1050. 1], &c__1) * work[p];
  1051. } else {
  1052. t = 0.;
  1053. aapp = 0.;
  1054. dlassq_(m, &a[p * a_dim1 + 1], &
  1055. c__1, &t, &aapp);
  1056. aapp = t * sqrt(aapp) * work[p];
  1057. }
  1058. sva[p] = aapp;
  1059. }
  1060. } else {
  1061. /* A(:,p) and A(:,q) already numerically orthogonal */
  1062. if (ir1 == 0) {
  1063. ++notrot;
  1064. }
  1065. /* [RTD] SKIPPED = SKIPPED + 1 */
  1066. ++pskipped;
  1067. }
  1068. } else {
  1069. /* A(:,q) is zero column */
  1070. if (ir1 == 0) {
  1071. ++notrot;
  1072. }
  1073. ++pskipped;
  1074. }
  1075. if (i__ <= swband && pskipped > rowskip) {
  1076. if (ir1 == 0) {
  1077. aapp = -aapp;
  1078. }
  1079. notrot = 0;
  1080. goto L2103;
  1081. }
  1082. /* L2002: */
  1083. }
  1084. /* END q-LOOP */
  1085. L2103:
  1086. /* bailed out of q-loop */
  1087. sva[p] = aapp;
  1088. } else {
  1089. sva[p] = aapp;
  1090. if (ir1 == 0 && aapp == 0.) {
  1091. /* Computing MIN */
  1092. i__4 = igl + kbl - 1;
  1093. notrot = notrot + min(i__4,*n) - p;
  1094. }
  1095. }
  1096. /* L2001: */
  1097. }
  1098. /* end of the p-loop */
  1099. /* end of doing the block ( ibr, ibr ) */
  1100. /* L1002: */
  1101. }
  1102. /* end of ir1-loop */
  1103. /* ... go to the off diagonal blocks */
  1104. igl = (ibr - 1) * kbl + 1;
  1105. i__2 = nbl;
  1106. for (jbc = ibr + 1; jbc <= i__2; ++jbc) {
  1107. jgl = (jbc - 1) * kbl + 1;
  1108. /* doing the block at ( ibr, jbc ) */
  1109. ijblsk = 0;
  1110. /* Computing MIN */
  1111. i__4 = igl + kbl - 1;
  1112. i__3 = min(i__4,*n);
  1113. for (p = igl; p <= i__3; ++p) {
  1114. aapp = sva[p];
  1115. if (aapp > 0.) {
  1116. pskipped = 0;
  1117. /* Computing MIN */
  1118. i__5 = jgl + kbl - 1;
  1119. i__4 = min(i__5,*n);
  1120. for (q = jgl; q <= i__4; ++q) {
  1121. aaqq = sva[q];
  1122. if (aaqq > 0.) {
  1123. aapp0 = aapp;
  1124. /* -#- M x 2 Jacobi SVD -#- */
  1125. /* Safe Gram matrix computation */
  1126. if (aaqq >= 1.) {
  1127. if (aapp >= aaqq) {
  1128. rotok = small * aapp <= aaqq;
  1129. } else {
  1130. rotok = small * aaqq <= aapp;
  1131. }
  1132. if (aapp < big / aaqq) {
  1133. aapq = ddot_(m, &a[p * a_dim1 + 1], &
  1134. c__1, &a[q * a_dim1 + 1], &
  1135. c__1) * work[p] * work[q] /
  1136. aaqq / aapp;
  1137. } else {
  1138. dcopy_(m, &a[p * a_dim1 + 1], &c__1, &
  1139. work[*n + 1], &c__1);
  1140. dlascl_("G", &c__0, &c__0, &aapp, &
  1141. work[p], m, &c__1, &work[*n +
  1142. 1], lda, &ierr);
  1143. aapq = ddot_(m, &work[*n + 1], &c__1,
  1144. &a[q * a_dim1 + 1], &c__1) *
  1145. work[q] / aaqq;
  1146. }
  1147. } else {
  1148. if (aapp >= aaqq) {
  1149. rotok = aapp <= aaqq / small;
  1150. } else {
  1151. rotok = aaqq <= aapp / small;
  1152. }
  1153. if (aapp > small / aaqq) {
  1154. aapq = ddot_(m, &a[p * a_dim1 + 1], &
  1155. c__1, &a[q * a_dim1 + 1], &
  1156. c__1) * work[p] * work[q] /
  1157. aaqq / aapp;
  1158. } else {
  1159. dcopy_(m, &a[q * a_dim1 + 1], &c__1, &
  1160. work[*n + 1], &c__1);
  1161. dlascl_("G", &c__0, &c__0, &aaqq, &
  1162. work[q], m, &c__1, &work[*n +
  1163. 1], lda, &ierr);
  1164. aapq = ddot_(m, &work[*n + 1], &c__1,
  1165. &a[p * a_dim1 + 1], &c__1) *
  1166. work[p] / aapp;
  1167. }
  1168. }
  1169. /* Computing MAX */
  1170. d__1 = mxaapq, d__2 = abs(aapq);
  1171. mxaapq = max(d__1,d__2);
  1172. /* TO rotate or NOT to rotate, THAT is the question ... */
  1173. if (abs(aapq) > tol) {
  1174. notrot = 0;
  1175. /* [RTD] ROTATED = ROTATED + 1 */
  1176. pskipped = 0;
  1177. ++iswrot;
  1178. if (rotok) {
  1179. aqoap = aaqq / aapp;
  1180. apoaq = aapp / aaqq;
  1181. theta = (d__1 = aqoap - apoaq, abs(
  1182. d__1)) * -.5 / aapq;
  1183. if (aaqq > aapp0) {
  1184. theta = -theta;
  1185. }
  1186. if (abs(theta) > bigtheta) {
  1187. t = .5 / theta;
  1188. fastr[2] = t * work[p] / work[q];
  1189. fastr[3] = -t * work[q] / work[p];
  1190. drotm_(m, &a[p * a_dim1 + 1], &
  1191. c__1, &a[q * a_dim1 + 1],
  1192. &c__1, fastr);
  1193. if (rsvec) {
  1194. drotm_(&mvl, &v[p * v_dim1 + 1], &c__1, &v[q *
  1195. v_dim1 + 1], &c__1, fastr);
  1196. }
  1197. /* Computing MAX */
  1198. d__1 = 0., d__2 = t * apoaq *
  1199. aapq + 1.;
  1200. sva[q] = aaqq * sqrt((max(d__1,
  1201. d__2)));
  1202. /* Computing MAX */
  1203. d__1 = 0., d__2 = 1. - t * aqoap *
  1204. aapq;
  1205. aapp *= sqrt((max(d__1,d__2)));
  1206. /* Computing MAX */
  1207. d__1 = mxsinj, d__2 = abs(t);
  1208. mxsinj = max(d__1,d__2);
  1209. } else {
  1210. /* .. choose correct signum for THETA and rotate */
  1211. thsign = -d_sign(&c_b18, &aapq);
  1212. if (aaqq > aapp0) {
  1213. thsign = -thsign;
  1214. }
  1215. t = 1. / (theta + thsign * sqrt(
  1216. theta * theta + 1.));
  1217. cs = sqrt(1. / (t * t + 1.));
  1218. sn = t * cs;
  1219. /* Computing MAX */
  1220. d__1 = mxsinj, d__2 = abs(sn);
  1221. mxsinj = max(d__1,d__2);
  1222. /* Computing MAX */
  1223. d__1 = 0., d__2 = t * apoaq *
  1224. aapq + 1.;
  1225. sva[q] = aaqq * sqrt((max(d__1,
  1226. d__2)));
  1227. aapp *= sqrt(1. - t * aqoap *
  1228. aapq);
  1229. apoaq = work[p] / work[q];
  1230. aqoap = work[q] / work[p];
  1231. if (work[p] >= 1.) {
  1232. if (work[q] >= 1.) {
  1233. fastr[2] = t * apoaq;
  1234. fastr[3] = -t * aqoap;
  1235. work[p] *= cs;
  1236. work[q] *= cs;
  1237. drotm_(m, &a[p * a_dim1 + 1], &c__1, &a[q *
  1238. a_dim1 + 1], &c__1, fastr);
  1239. if (rsvec) {
  1240. drotm_(&mvl, &v[p * v_dim1 + 1], &c__1, &v[
  1241. q * v_dim1 + 1], &c__1, fastr);
  1242. }
  1243. } else {
  1244. d__1 = -t * aqoap;
  1245. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1, &a[
  1246. p * a_dim1 + 1], &c__1);
  1247. d__1 = cs * sn * apoaq;
  1248. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1, &a[
  1249. q * a_dim1 + 1], &c__1);
  1250. if (rsvec) {
  1251. d__1 = -t * aqoap;
  1252. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1], &
  1253. c__1, &v[p * v_dim1 + 1], &c__1);
  1254. d__1 = cs * sn * apoaq;
  1255. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1], &
  1256. c__1, &v[q * v_dim1 + 1], &c__1);
  1257. }
  1258. work[p] *= cs;
  1259. work[q] /= cs;
  1260. }
  1261. } else {
  1262. if (work[q] >= 1.) {
  1263. d__1 = t * apoaq;
  1264. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1, &a[
  1265. q * a_dim1 + 1], &c__1);
  1266. d__1 = -cs * sn * aqoap;
  1267. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1, &a[
  1268. p * a_dim1 + 1], &c__1);
  1269. if (rsvec) {
  1270. d__1 = t * apoaq;
  1271. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1], &
  1272. c__1, &v[q * v_dim1 + 1], &c__1);
  1273. d__1 = -cs * sn * aqoap;
  1274. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1], &
  1275. c__1, &v[p * v_dim1 + 1], &c__1);
  1276. }
  1277. work[p] /= cs;
  1278. work[q] *= cs;
  1279. } else {
  1280. if (work[p] >= work[q]) {
  1281. d__1 = -t * aqoap;
  1282. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1,
  1283. &a[p * a_dim1 + 1], &c__1);
  1284. d__1 = cs * sn * apoaq;
  1285. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1,
  1286. &a[q * a_dim1 + 1], &c__1);
  1287. work[p] *= cs;
  1288. work[q] /= cs;
  1289. if (rsvec) {
  1290. d__1 = -t * aqoap;
  1291. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1],
  1292. &c__1, &v[p * v_dim1 + 1], &
  1293. c__1);
  1294. d__1 = cs * sn * apoaq;
  1295. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1],
  1296. &c__1, &v[q * v_dim1 + 1], &
  1297. c__1);
  1298. }
  1299. } else {
  1300. d__1 = t * apoaq;
  1301. daxpy_(m, &d__1, &a[p * a_dim1 + 1], &c__1,
  1302. &a[q * a_dim1 + 1], &c__1);
  1303. d__1 = -cs * sn * aqoap;
  1304. daxpy_(m, &d__1, &a[q * a_dim1 + 1], &c__1,
  1305. &a[p * a_dim1 + 1], &c__1);
  1306. work[p] /= cs;
  1307. work[q] *= cs;
  1308. if (rsvec) {
  1309. d__1 = t * apoaq;
  1310. daxpy_(&mvl, &d__1, &v[p * v_dim1 + 1],
  1311. &c__1, &v[q * v_dim1 + 1], &
  1312. c__1);
  1313. d__1 = -cs * sn * aqoap;
  1314. daxpy_(&mvl, &d__1, &v[q * v_dim1 + 1],
  1315. &c__1, &v[p * v_dim1 + 1], &
  1316. c__1);
  1317. }
  1318. }
  1319. }
  1320. }
  1321. }
  1322. } else {
  1323. if (aapp > aaqq) {
  1324. dcopy_(m, &a[p * a_dim1 + 1], &
  1325. c__1, &work[*n + 1], &
  1326. c__1);
  1327. dlascl_("G", &c__0, &c__0, &aapp,
  1328. &c_b18, m, &c__1, &work[*
  1329. n + 1], lda, &ierr);
  1330. dlascl_("G", &c__0, &c__0, &aaqq,
  1331. &c_b18, m, &c__1, &a[q *
  1332. a_dim1 + 1], lda, &ierr);
  1333. temp1 = -aapq * work[p] / work[q];
  1334. daxpy_(m, &temp1, &work[*n + 1], &
  1335. c__1, &a[q * a_dim1 + 1],
  1336. &c__1);
  1337. dlascl_("G", &c__0, &c__0, &c_b18,
  1338. &aaqq, m, &c__1, &a[q *
  1339. a_dim1 + 1], lda, &ierr);
  1340. /* Computing MAX */
  1341. d__1 = 0., d__2 = 1. - aapq *
  1342. aapq;
  1343. sva[q] = aaqq * sqrt((max(d__1,
  1344. d__2)));
  1345. mxsinj = max(mxsinj,sfmin);
  1346. } else {
  1347. dcopy_(m, &a[q * a_dim1 + 1], &
  1348. c__1, &work[*n + 1], &
  1349. c__1);
  1350. dlascl_("G", &c__0, &c__0, &aaqq,
  1351. &c_b18, m, &c__1, &work[*
  1352. n + 1], lda, &ierr);
  1353. dlascl_("G", &c__0, &c__0, &aapp,
  1354. &c_b18, m, &c__1, &a[p *
  1355. a_dim1 + 1], lda, &ierr);
  1356. temp1 = -aapq * work[q] / work[p];
  1357. daxpy_(m, &temp1, &work[*n + 1], &
  1358. c__1, &a[p * a_dim1 + 1],
  1359. &c__1);
  1360. dlascl_("G", &c__0, &c__0, &c_b18,
  1361. &aapp, m, &c__1, &a[p *
  1362. a_dim1 + 1], lda, &ierr);
  1363. /* Computing MAX */
  1364. d__1 = 0., d__2 = 1. - aapq *
  1365. aapq;
  1366. sva[p] = aapp * sqrt((max(d__1,
  1367. d__2)));
  1368. mxsinj = max(mxsinj,sfmin);
  1369. }
  1370. }
  1371. /* END IF ROTOK THEN ... ELSE */
  1372. /* In the case of cancellation in updating SVA(q) */
  1373. /* .. recompute SVA(q) */
  1374. /* Computing 2nd power */
  1375. d__1 = sva[q] / aaqq;
  1376. if (d__1 * d__1 <= rooteps) {
  1377. if (aaqq < rootbig && aaqq >
  1378. rootsfmin) {
  1379. sva[q] = dnrm2_(m, &a[q * a_dim1
  1380. + 1], &c__1) * work[q];
  1381. } else {
  1382. t = 0.;
  1383. aaqq = 0.;
  1384. dlassq_(m, &a[q * a_dim1 + 1], &
  1385. c__1, &t, &aaqq);
  1386. sva[q] = t * sqrt(aaqq) * work[q];
  1387. }
  1388. }
  1389. /* Computing 2nd power */
  1390. d__1 = aapp / aapp0;
  1391. if (d__1 * d__1 <= rooteps) {
  1392. if (aapp < rootbig && aapp >
  1393. rootsfmin) {
  1394. aapp = dnrm2_(m, &a[p * a_dim1 +
  1395. 1], &c__1) * work[p];
  1396. } else {
  1397. t = 0.;
  1398. aapp = 0.;
  1399. dlassq_(m, &a[p * a_dim1 + 1], &
  1400. c__1, &t, &aapp);
  1401. aapp = t * sqrt(aapp) * work[p];
  1402. }
  1403. sva[p] = aapp;
  1404. }
  1405. /* end of OK rotation */
  1406. } else {
  1407. ++notrot;
  1408. /* [RTD] SKIPPED = SKIPPED + 1 */
  1409. ++pskipped;
  1410. ++ijblsk;
  1411. }
  1412. } else {
  1413. ++notrot;
  1414. ++pskipped;
  1415. ++ijblsk;
  1416. }
  1417. if (i__ <= swband && ijblsk >= blskip) {
  1418. sva[p] = aapp;
  1419. notrot = 0;
  1420. goto L2011;
  1421. }
  1422. if (i__ <= swband && pskipped > rowskip) {
  1423. aapp = -aapp;
  1424. notrot = 0;
  1425. goto L2203;
  1426. }
  1427. /* L2200: */
  1428. }
  1429. /* end of the q-loop */
  1430. L2203:
  1431. sva[p] = aapp;
  1432. } else {
  1433. if (aapp == 0.) {
  1434. /* Computing MIN */
  1435. i__4 = jgl + kbl - 1;
  1436. notrot = notrot + min(i__4,*n) - jgl + 1;
  1437. }
  1438. if (aapp < 0.) {
  1439. notrot = 0;
  1440. }
  1441. }
  1442. /* L2100: */
  1443. }
  1444. /* end of the p-loop */
  1445. /* L2010: */
  1446. }
  1447. /* end of the jbc-loop */
  1448. L2011:
  1449. /* 2011 bailed out of the jbc-loop */
  1450. /* Computing MIN */
  1451. i__3 = igl + kbl - 1;
  1452. i__2 = min(i__3,*n);
  1453. for (p = igl; p <= i__2; ++p) {
  1454. sva[p] = (d__1 = sva[p], abs(d__1));
  1455. /* L2012: */
  1456. }
  1457. /* ** */
  1458. /* L2000: */
  1459. }
  1460. /* 2000 :: end of the ibr-loop */
  1461. /* .. update SVA(N) */
  1462. if (sva[*n] < rootbig && sva[*n] > rootsfmin) {
  1463. sva[*n] = dnrm2_(m, &a[*n * a_dim1 + 1], &c__1) * work[*n];
  1464. } else {
  1465. t = 0.;
  1466. aapp = 0.;
  1467. dlassq_(m, &a[*n * a_dim1 + 1], &c__1, &t, &aapp);
  1468. sva[*n] = t * sqrt(aapp) * work[*n];
  1469. }
  1470. /* Additional steering devices */
  1471. if (i__ < swband && (mxaapq <= roottol || iswrot <= *n)) {
  1472. swband = i__;
  1473. }
  1474. if (i__ > swband + 1 && mxaapq < sqrt((doublereal) (*n)) * tol && (
  1475. doublereal) (*n) * mxaapq * mxsinj < tol) {
  1476. goto L1994;
  1477. }
  1478. if (notrot >= emptsw) {
  1479. goto L1994;
  1480. }
  1481. /* L1993: */
  1482. }
  1483. /* end i=1:NSWEEP loop */
  1484. /* #:( Reaching this point means that the procedure has not converged. */
  1485. *info = 29;
  1486. goto L1995;
  1487. L1994:
  1488. /* #:) Reaching this point means numerical convergence after the i-th */
  1489. /* sweep. */
  1490. *info = 0;
  1491. /* #:) INFO = 0 confirms successful iterations. */
  1492. L1995:
  1493. /* Sort the singular values and find how many are above */
  1494. /* the underflow threshold. */
  1495. n2 = 0;
  1496. n4 = 0;
  1497. i__1 = *n - 1;
  1498. for (p = 1; p <= i__1; ++p) {
  1499. i__2 = *n - p + 1;
  1500. q = idamax_(&i__2, &sva[p], &c__1) + p - 1;
  1501. if (p != q) {
  1502. temp1 = sva[p];
  1503. sva[p] = sva[q];
  1504. sva[q] = temp1;
  1505. temp1 = work[p];
  1506. work[p] = work[q];
  1507. work[q] = temp1;
  1508. dswap_(m, &a[p * a_dim1 + 1], &c__1, &a[q * a_dim1 + 1], &c__1);
  1509. if (rsvec) {
  1510. dswap_(&mvl, &v[p * v_dim1 + 1], &c__1, &v[q * v_dim1 + 1], &
  1511. c__1);
  1512. }
  1513. }
  1514. if (sva[p] != 0.) {
  1515. ++n4;
  1516. if (sva[p] * scale > sfmin) {
  1517. ++n2;
  1518. }
  1519. }
  1520. /* L5991: */
  1521. }
  1522. if (sva[*n] != 0.) {
  1523. ++n4;
  1524. if (sva[*n] * scale > sfmin) {
  1525. ++n2;
  1526. }
  1527. }
  1528. /* Normalize the left singular vectors. */
  1529. if (lsvec || uctol) {
  1530. i__1 = n2;
  1531. for (p = 1; p <= i__1; ++p) {
  1532. d__1 = work[p] / sva[p];
  1533. dscal_(m, &d__1, &a[p * a_dim1 + 1], &c__1);
  1534. /* L1998: */
  1535. }
  1536. }
  1537. /* Scale the product of Jacobi rotations (assemble the fast rotations). */
  1538. if (rsvec) {
  1539. if (applv) {
  1540. i__1 = *n;
  1541. for (p = 1; p <= i__1; ++p) {
  1542. dscal_(&mvl, &work[p], &v[p * v_dim1 + 1], &c__1);
  1543. /* L2398: */
  1544. }
  1545. } else {
  1546. i__1 = *n;
  1547. for (p = 1; p <= i__1; ++p) {
  1548. temp1 = 1. / dnrm2_(&mvl, &v[p * v_dim1 + 1], &c__1);
  1549. dscal_(&mvl, &temp1, &v[p * v_dim1 + 1], &c__1);
  1550. /* L2399: */
  1551. }
  1552. }
  1553. }
  1554. /* Undo scaling, if necessary (and possible). */
  1555. if (scale > 1. && sva[1] < big / scale || scale < 1. && sva[n2] > sfmin /
  1556. scale) {
  1557. i__1 = *n;
  1558. for (p = 1; p <= i__1; ++p) {
  1559. sva[p] = scale * sva[p];
  1560. /* L2400: */
  1561. }
  1562. scale = 1.;
  1563. }
  1564. work[1] = scale;
  1565. /* The singular values of A are SCALE*SVA(1:N). If SCALE.NE.ONE */
  1566. /* then some of the singular values may overflow or underflow and */
  1567. /* the spectrum is given in this factored representation. */
  1568. work[2] = (doublereal) n4;
  1569. /* N4 is the number of computed nonzero singular values of A. */
  1570. work[3] = (doublereal) n2;
  1571. /* N2 is the number of singular values of A greater than SFMIN. */
  1572. /* If N2<N, SVA(N2:N) contains ZEROS and/or denormalized numbers */
  1573. /* that may carry some information. */
  1574. work[4] = (doublereal) i__;
  1575. /* i is the index of the last sweep before declaring convergence. */
  1576. work[5] = mxaapq;
  1577. /* MXAAPQ is the largest absolute value of scaled pivots in the */
  1578. /* last sweep */
  1579. work[6] = mxsinj;
  1580. /* MXSINJ is the largest absolute value of the sines of Jacobi angles */
  1581. /* in the last sweep */
  1582. return 0;
  1583. /* .. */
  1584. /* .. END OF DGESVJ */
  1585. /* .. */
  1586. } /* dgesvj_ */