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.

sgelsx.c 22 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755
  1. #include <math.h>
  2. #include <stdlib.h>
  3. #include <string.h>
  4. #include <stdio.h>
  5. #include <complex.h>
  6. #ifdef complex
  7. #undef complex
  8. #endif
  9. #ifdef I
  10. #undef I
  11. #endif
  12. #if defined(_WIN64)
  13. typedef long long BLASLONG;
  14. typedef unsigned long long BLASULONG;
  15. #else
  16. typedef long BLASLONG;
  17. typedef unsigned long BLASULONG;
  18. #endif
  19. #ifdef LAPACK_ILP64
  20. typedef BLASLONG blasint;
  21. #if defined(_WIN64)
  22. #define blasabs(x) llabs(x)
  23. #else
  24. #define blasabs(x) labs(x)
  25. #endif
  26. #else
  27. typedef int blasint;
  28. #define blasabs(x) abs(x)
  29. #endif
  30. typedef blasint integer;
  31. typedef unsigned int uinteger;
  32. typedef char *address;
  33. typedef short int shortint;
  34. typedef float real;
  35. typedef double doublereal;
  36. typedef struct { real r, i; } complex;
  37. typedef struct { doublereal r, i; } doublecomplex;
  38. #ifdef _MSC_VER
  39. static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
  40. static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
  41. static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
  42. static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
  43. #else
  44. static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
  45. static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
  46. static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
  47. static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
  48. #endif
  49. #define pCf(z) (*_pCf(z))
  50. #define pCd(z) (*_pCd(z))
  51. typedef char integer1;
  52. #define TRUE_ (1)
  53. #define FALSE_ (0)
  54. /* Extern is for use with -E */
  55. #ifndef Extern
  56. #define Extern extern
  57. #endif
  58. /* I/O stuff */
  59. typedef int flag;
  60. typedef int ftnlen;
  61. typedef int ftnint;
  62. /*external read, write*/
  63. typedef struct
  64. { flag cierr;
  65. ftnint ciunit;
  66. flag ciend;
  67. char *cifmt;
  68. ftnint cirec;
  69. } cilist;
  70. /*internal read, write*/
  71. typedef struct
  72. { flag icierr;
  73. char *iciunit;
  74. flag iciend;
  75. char *icifmt;
  76. ftnint icirlen;
  77. ftnint icirnum;
  78. } icilist;
  79. /*open*/
  80. typedef struct
  81. { flag oerr;
  82. ftnint ounit;
  83. char *ofnm;
  84. ftnlen ofnmlen;
  85. char *osta;
  86. char *oacc;
  87. char *ofm;
  88. ftnint orl;
  89. char *oblnk;
  90. } olist;
  91. /*close*/
  92. typedef struct
  93. { flag cerr;
  94. ftnint cunit;
  95. char *csta;
  96. } cllist;
  97. /*rewind, backspace, endfile*/
  98. typedef struct
  99. { flag aerr;
  100. ftnint aunit;
  101. } alist;
  102. /* inquire */
  103. typedef struct
  104. { flag inerr;
  105. ftnint inunit;
  106. char *infile;
  107. ftnlen infilen;
  108. ftnint *inex; /*parameters in standard's order*/
  109. ftnint *inopen;
  110. ftnint *innum;
  111. ftnint *innamed;
  112. char *inname;
  113. ftnlen innamlen;
  114. char *inacc;
  115. ftnlen inacclen;
  116. char *inseq;
  117. ftnlen inseqlen;
  118. char *indir;
  119. ftnlen indirlen;
  120. char *infmt;
  121. ftnlen infmtlen;
  122. char *inform;
  123. ftnint informlen;
  124. char *inunf;
  125. ftnlen inunflen;
  126. ftnint *inrecl;
  127. ftnint *innrec;
  128. char *inblank;
  129. ftnlen inblanklen;
  130. } inlist;
  131. #define VOID void
  132. union Multitype { /* for multiple entry points */
  133. integer1 g;
  134. shortint h;
  135. integer i;
  136. /* longint j; */
  137. real r;
  138. doublereal d;
  139. complex c;
  140. doublecomplex z;
  141. };
  142. typedef union Multitype Multitype;
  143. struct Vardesc { /* for Namelist */
  144. char *name;
  145. char *addr;
  146. ftnlen *dims;
  147. int type;
  148. };
  149. typedef struct Vardesc Vardesc;
  150. struct Namelist {
  151. char *name;
  152. Vardesc **vars;
  153. int nvars;
  154. };
  155. typedef struct Namelist Namelist;
  156. #define abs(x) ((x) >= 0 ? (x) : -(x))
  157. #define dabs(x) (fabs(x))
  158. #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
  159. #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
  160. #define dmin(a,b) (f2cmin(a,b))
  161. #define dmax(a,b) (f2cmax(a,b))
  162. #define bit_test(a,b) ((a) >> (b) & 1)
  163. #define bit_clear(a,b) ((a) & ~((uinteger)1 << (b)))
  164. #define bit_set(a,b) ((a) | ((uinteger)1 << (b)))
  165. #define abort_() { sig_die("Fortran abort routine called", 1); }
  166. #define c_abs(z) (cabsf(Cf(z)))
  167. #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
  168. #ifdef _MSC_VER
  169. #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
  170. #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
  171. #else
  172. #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
  173. #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
  174. #endif
  175. #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
  176. #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
  177. #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
  178. //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
  179. #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
  180. #define d_abs(x) (fabs(*(x)))
  181. #define d_acos(x) (acos(*(x)))
  182. #define d_asin(x) (asin(*(x)))
  183. #define d_atan(x) (atan(*(x)))
  184. #define d_atn2(x, y) (atan2(*(x),*(y)))
  185. #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
  186. #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
  187. #define d_cos(x) (cos(*(x)))
  188. #define d_cosh(x) (cosh(*(x)))
  189. #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
  190. #define d_exp(x) (exp(*(x)))
  191. #define d_imag(z) (cimag(Cd(z)))
  192. #define r_imag(z) (cimagf(Cf(z)))
  193. #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  194. #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  195. #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  196. #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  197. #define d_log(x) (log(*(x)))
  198. #define d_mod(x, y) (fmod(*(x), *(y)))
  199. #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
  200. #define d_nint(x) u_nint(*(x))
  201. #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
  202. #define d_sign(a,b) u_sign(*(a),*(b))
  203. #define r_sign(a,b) u_sign(*(a),*(b))
  204. #define d_sin(x) (sin(*(x)))
  205. #define d_sinh(x) (sinh(*(x)))
  206. #define d_sqrt(x) (sqrt(*(x)))
  207. #define d_tan(x) (tan(*(x)))
  208. #define d_tanh(x) (tanh(*(x)))
  209. #define i_abs(x) abs(*(x))
  210. #define i_dnnt(x) ((integer)u_nint(*(x)))
  211. #define i_len(s, n) (n)
  212. #define i_nint(x) ((integer)u_nint(*(x)))
  213. #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
  214. #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
  215. #define pow_si(B,E) spow_ui(*(B),*(E))
  216. #define pow_ri(B,E) spow_ui(*(B),*(E))
  217. #define pow_di(B,E) dpow_ui(*(B),*(E))
  218. #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
  219. #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
  220. #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
  221. #define s_cat(lpp, rpp, rnp, np, llp) { ftnlen i, nc, ll; char *f__rp, *lp; ll = (llp); lp = (lpp); for(i=0; i < (int)*(np); ++i) { nc = ll; if((rnp)[i] < nc) nc = (rnp)[i]; ll -= nc; f__rp = (rpp)[i]; while(--nc >= 0) *lp++ = *(f__rp)++; } while(--ll >= 0) *lp++ = ' '; }
  222. #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
  223. #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
  224. #define sig_die(s, kill) { exit(1); }
  225. #define s_stop(s, n) {exit(0);}
  226. #define z_abs(z) (cabs(Cd(z)))
  227. #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
  228. #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
  229. #define myexit_() break;
  230. #define mycycle() continue;
  231. #define myceiling(w) {ceil(w)}
  232. #define myhuge(w) {HUGE_VAL}
  233. //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
  234. #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
  235. /* procedure parameter types for -A and -C++ */
  236. /* -- translated by f2c (version 20000121).
  237. You must link the resulting object file with the libraries:
  238. -lf2c -lm (in that order)
  239. */
  240. /* Table of constant values */
  241. static integer c__0 = 0;
  242. static real c_b13 = 0.f;
  243. static integer c__2 = 2;
  244. static integer c__1 = 1;
  245. static real c_b36 = 1.f;
  246. /* > \brief <b> SGELSX solves overdetermined or underdetermined systems for GE matrices</b> */
  247. /* =========== DOCUMENTATION =========== */
  248. /* Online html documentation available at */
  249. /* http://www.netlib.org/lapack/explore-html/ */
  250. /* > \htmlonly */
  251. /* > Download SGELSX + dependencies */
  252. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgelsx.
  253. f"> */
  254. /* > [TGZ]</a> */
  255. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgelsx.
  256. f"> */
  257. /* > [ZIP]</a> */
  258. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgelsx.
  259. f"> */
  260. /* > [TXT]</a> */
  261. /* > \endhtmlonly */
  262. /* Definition: */
  263. /* =========== */
  264. /* SUBROUTINE SGELSX( M, N, NRHS, A, LDA, B, LDB, JPVT, RCOND, RANK, */
  265. /* WORK, INFO ) */
  266. /* INTEGER INFO, LDA, LDB, M, N, NRHS, RANK */
  267. /* REAL RCOND */
  268. /* INTEGER JPVT( * ) */
  269. /* REAL A( LDA, * ), B( LDB, * ), WORK( * ) */
  270. /* > \par Purpose: */
  271. /* ============= */
  272. /* > */
  273. /* > \verbatim */
  274. /* > */
  275. /* > This routine is deprecated and has been replaced by routine SGELSY. */
  276. /* > */
  277. /* > SGELSX computes the minimum-norm solution to a real linear least */
  278. /* > squares problem: */
  279. /* > minimize || A * X - B || */
  280. /* > using a complete orthogonal factorization of A. A is an M-by-N */
  281. /* > matrix which may be rank-deficient. */
  282. /* > */
  283. /* > Several right hand side vectors b and solution vectors x can be */
  284. /* > handled in a single call; they are stored as the columns of the */
  285. /* > M-by-NRHS right hand side matrix B and the N-by-NRHS solution */
  286. /* > matrix X. */
  287. /* > */
  288. /* > The routine first computes a QR factorization with column pivoting: */
  289. /* > A * P = Q * [ R11 R12 ] */
  290. /* > [ 0 R22 ] */
  291. /* > with R11 defined as the largest leading submatrix whose estimated */
  292. /* > condition number is less than 1/RCOND. The order of R11, RANK, */
  293. /* > is the effective rank of A. */
  294. /* > */
  295. /* > Then, R22 is considered to be negligible, and R12 is annihilated */
  296. /* > by orthogonal transformations from the right, arriving at the */
  297. /* > complete orthogonal factorization: */
  298. /* > A * P = Q * [ T11 0 ] * Z */
  299. /* > [ 0 0 ] */
  300. /* > The minimum-norm solution is then */
  301. /* > X = P * Z**T [ inv(T11)*Q1**T*B ] */
  302. /* > [ 0 ] */
  303. /* > where Q1 consists of the first RANK columns of Q. */
  304. /* > \endverbatim */
  305. /* Arguments: */
  306. /* ========== */
  307. /* > \param[in] M */
  308. /* > \verbatim */
  309. /* > M is INTEGER */
  310. /* > The number of rows of the matrix A. M >= 0. */
  311. /* > \endverbatim */
  312. /* > */
  313. /* > \param[in] N */
  314. /* > \verbatim */
  315. /* > N is INTEGER */
  316. /* > The number of columns of the matrix A. N >= 0. */
  317. /* > \endverbatim */
  318. /* > */
  319. /* > \param[in] NRHS */
  320. /* > \verbatim */
  321. /* > NRHS is INTEGER */
  322. /* > The number of right hand sides, i.e., the number of */
  323. /* > columns of matrices B and X. NRHS >= 0. */
  324. /* > \endverbatim */
  325. /* > */
  326. /* > \param[in,out] A */
  327. /* > \verbatim */
  328. /* > A is REAL array, dimension (LDA,N) */
  329. /* > On entry, the M-by-N matrix A. */
  330. /* > On exit, A has been overwritten by details of its */
  331. /* > complete orthogonal factorization. */
  332. /* > \endverbatim */
  333. /* > */
  334. /* > \param[in] LDA */
  335. /* > \verbatim */
  336. /* > LDA is INTEGER */
  337. /* > The leading dimension of the array A. LDA >= f2cmax(1,M). */
  338. /* > \endverbatim */
  339. /* > */
  340. /* > \param[in,out] B */
  341. /* > \verbatim */
  342. /* > B is REAL array, dimension (LDB,NRHS) */
  343. /* > On entry, the M-by-NRHS right hand side matrix B. */
  344. /* > On exit, the N-by-NRHS solution matrix X. */
  345. /* > If m >= n and RANK = n, the residual sum-of-squares for */
  346. /* > the solution in the i-th column is given by the sum of */
  347. /* > squares of elements N+1:M in that column. */
  348. /* > \endverbatim */
  349. /* > */
  350. /* > \param[in] LDB */
  351. /* > \verbatim */
  352. /* > LDB is INTEGER */
  353. /* > The leading dimension of the array B. LDB >= f2cmax(1,M,N). */
  354. /* > \endverbatim */
  355. /* > */
  356. /* > \param[in,out] JPVT */
  357. /* > \verbatim */
  358. /* > JPVT is INTEGER array, dimension (N) */
  359. /* > On entry, if JPVT(i) .ne. 0, the i-th column of A is an */
  360. /* > initial column, otherwise it is a free column. Before */
  361. /* > the QR factorization of A, all initial columns are */
  362. /* > permuted to the leading positions; only the remaining */
  363. /* > free columns are moved as a result of column pivoting */
  364. /* > during the factorization. */
  365. /* > On exit, if JPVT(i) = k, then the i-th column of A*P */
  366. /* > was the k-th column of A. */
  367. /* > \endverbatim */
  368. /* > */
  369. /* > \param[in] RCOND */
  370. /* > \verbatim */
  371. /* > RCOND is REAL */
  372. /* > RCOND is used to determine the effective rank of A, which */
  373. /* > is defined as the order of the largest leading triangular */
  374. /* > submatrix R11 in the QR factorization with pivoting of A, */
  375. /* > whose estimated condition number < 1/RCOND. */
  376. /* > \endverbatim */
  377. /* > */
  378. /* > \param[out] RANK */
  379. /* > \verbatim */
  380. /* > RANK is INTEGER */
  381. /* > The effective rank of A, i.e., the order of the submatrix */
  382. /* > R11. This is the same as the order of the submatrix T11 */
  383. /* > in the complete orthogonal factorization of A. */
  384. /* > \endverbatim */
  385. /* > */
  386. /* > \param[out] WORK */
  387. /* > \verbatim */
  388. /* > WORK is REAL array, dimension */
  389. /* > (f2cmax( f2cmin(M,N)+3*N, 2*f2cmin(M,N)+NRHS )), */
  390. /* > \endverbatim */
  391. /* > */
  392. /* > \param[out] INFO */
  393. /* > \verbatim */
  394. /* > INFO is INTEGER */
  395. /* > = 0: successful exit */
  396. /* > < 0: if INFO = -i, the i-th argument had an illegal value */
  397. /* > \endverbatim */
  398. /* Authors: */
  399. /* ======== */
  400. /* > \author Univ. of Tennessee */
  401. /* > \author Univ. of California Berkeley */
  402. /* > \author Univ. of Colorado Denver */
  403. /* > \author NAG Ltd. */
  404. /* > \date December 2016 */
  405. /* > \ingroup realGEsolve */
  406. /* ===================================================================== */
  407. /* Subroutine */ void sgelsx_(integer *m, integer *n, integer *nrhs, real *a,
  408. integer *lda, real *b, integer *ldb, integer *jpvt, real *rcond,
  409. integer *rank, real *work, integer *info)
  410. {
  411. /* System generated locals */
  412. integer a_dim1, a_offset, b_dim1, b_offset, i__1, i__2;
  413. real r__1;
  414. /* Local variables */
  415. real anrm, bnrm, smin, smax;
  416. integer i__, j, k, iascl, ibscl, ismin, ismax;
  417. real c1, c2, s1, s2, t1, t2;
  418. extern /* Subroutine */ void strsm_(char *, char *, char *, char *,
  419. integer *, integer *, real *, real *, integer *, real *, integer *
  420. ), slaic1_(integer *, integer *,
  421. real *, real *, real *, real *, real *, real *, real *), sorm2r_(
  422. char *, char *, integer *, integer *, integer *, real *, integer *
  423. , real *, real *, integer *, real *, integer *),
  424. slabad_(real *, real *);
  425. integer mn;
  426. extern real slamch_(char *), slange_(char *, integer *, integer *,
  427. real *, integer *, real *);
  428. extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
  429. real bignum;
  430. extern /* Subroutine */ void slascl_(char *, integer *, integer *, real *,
  431. real *, integer *, integer *, real *, integer *, integer *), sgeqpf_(integer *, integer *, real *, integer *, integer
  432. *, real *, real *, integer *), slaset_(char *, integer *, integer
  433. *, real *, real *, real *, integer *);
  434. real sminpr, smaxpr, smlnum;
  435. extern /* Subroutine */ void slatzm_(char *, integer *, integer *, real *,
  436. integer *, real *, real *, real *, integer *, real *),
  437. stzrqf_(integer *, integer *, real *, integer *, real *, integer *
  438. );
  439. /* -- LAPACK driver routine (version 3.7.0) -- */
  440. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  441. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  442. /* December 2016 */
  443. /* ===================================================================== */
  444. /* Parameter adjustments */
  445. a_dim1 = *lda;
  446. a_offset = 1 + a_dim1 * 1;
  447. a -= a_offset;
  448. b_dim1 = *ldb;
  449. b_offset = 1 + b_dim1 * 1;
  450. b -= b_offset;
  451. --jpvt;
  452. --work;
  453. /* Function Body */
  454. mn = f2cmin(*m,*n);
  455. ismin = mn + 1;
  456. ismax = (mn << 1) + 1;
  457. /* Test the input arguments. */
  458. *info = 0;
  459. if (*m < 0) {
  460. *info = -1;
  461. } else if (*n < 0) {
  462. *info = -2;
  463. } else if (*nrhs < 0) {
  464. *info = -3;
  465. } else if (*lda < f2cmax(1,*m)) {
  466. *info = -5;
  467. } else /* if(complicated condition) */ {
  468. /* Computing MAX */
  469. i__1 = f2cmax(1,*m);
  470. if (*ldb < f2cmax(i__1,*n)) {
  471. *info = -7;
  472. }
  473. }
  474. if (*info != 0) {
  475. i__1 = -(*info);
  476. xerbla_("SGELSX", &i__1, 6);
  477. return;
  478. }
  479. /* Quick return if possible */
  480. /* Computing MIN */
  481. i__1 = f2cmin(*m,*n);
  482. if (f2cmin(i__1,*nrhs) == 0) {
  483. *rank = 0;
  484. return;
  485. }
  486. /* Get machine parameters */
  487. smlnum = slamch_("S") / slamch_("P");
  488. bignum = 1.f / smlnum;
  489. slabad_(&smlnum, &bignum);
  490. /* Scale A, B if f2cmax elements outside range [SMLNUM,BIGNUM] */
  491. anrm = slange_("M", m, n, &a[a_offset], lda, &work[1]);
  492. iascl = 0;
  493. if (anrm > 0.f && anrm < smlnum) {
  494. /* Scale matrix norm up to SMLNUM */
  495. slascl_("G", &c__0, &c__0, &anrm, &smlnum, m, n, &a[a_offset], lda,
  496. info);
  497. iascl = 1;
  498. } else if (anrm > bignum) {
  499. /* Scale matrix norm down to BIGNUM */
  500. slascl_("G", &c__0, &c__0, &anrm, &bignum, m, n, &a[a_offset], lda,
  501. info);
  502. iascl = 2;
  503. } else if (anrm == 0.f) {
  504. /* Matrix all zero. Return zero solution. */
  505. i__1 = f2cmax(*m,*n);
  506. slaset_("F", &i__1, nrhs, &c_b13, &c_b13, &b[b_offset], ldb);
  507. *rank = 0;
  508. goto L100;
  509. }
  510. bnrm = slange_("M", m, nrhs, &b[b_offset], ldb, &work[1]);
  511. ibscl = 0;
  512. if (bnrm > 0.f && bnrm < smlnum) {
  513. /* Scale matrix norm up to SMLNUM */
  514. slascl_("G", &c__0, &c__0, &bnrm, &smlnum, m, nrhs, &b[b_offset], ldb,
  515. info);
  516. ibscl = 1;
  517. } else if (bnrm > bignum) {
  518. /* Scale matrix norm down to BIGNUM */
  519. slascl_("G", &c__0, &c__0, &bnrm, &bignum, m, nrhs, &b[b_offset], ldb,
  520. info);
  521. ibscl = 2;
  522. }
  523. /* Compute QR factorization with column pivoting of A: */
  524. /* A * P = Q * R */
  525. sgeqpf_(m, n, &a[a_offset], lda, &jpvt[1], &work[1], &work[mn + 1], info);
  526. /* workspace 3*N. Details of Householder rotations stored */
  527. /* in WORK(1:MN). */
  528. /* Determine RANK using incremental condition estimation */
  529. work[ismin] = 1.f;
  530. work[ismax] = 1.f;
  531. smax = (r__1 = a[a_dim1 + 1], abs(r__1));
  532. smin = smax;
  533. if ((r__1 = a[a_dim1 + 1], abs(r__1)) == 0.f) {
  534. *rank = 0;
  535. i__1 = f2cmax(*m,*n);
  536. slaset_("F", &i__1, nrhs, &c_b13, &c_b13, &b[b_offset], ldb);
  537. goto L100;
  538. } else {
  539. *rank = 1;
  540. }
  541. L10:
  542. if (*rank < mn) {
  543. i__ = *rank + 1;
  544. slaic1_(&c__2, rank, &work[ismin], &smin, &a[i__ * a_dim1 + 1], &a[
  545. i__ + i__ * a_dim1], &sminpr, &s1, &c1);
  546. slaic1_(&c__1, rank, &work[ismax], &smax, &a[i__ * a_dim1 + 1], &a[
  547. i__ + i__ * a_dim1], &smaxpr, &s2, &c2);
  548. if (smaxpr * *rcond <= sminpr) {
  549. i__1 = *rank;
  550. for (i__ = 1; i__ <= i__1; ++i__) {
  551. work[ismin + i__ - 1] = s1 * work[ismin + i__ - 1];
  552. work[ismax + i__ - 1] = s2 * work[ismax + i__ - 1];
  553. /* L20: */
  554. }
  555. work[ismin + *rank] = c1;
  556. work[ismax + *rank] = c2;
  557. smin = sminpr;
  558. smax = smaxpr;
  559. ++(*rank);
  560. goto L10;
  561. }
  562. }
  563. /* Logically partition R = [ R11 R12 ] */
  564. /* [ 0 R22 ] */
  565. /* where R11 = R(1:RANK,1:RANK) */
  566. /* [R11,R12] = [ T11, 0 ] * Y */
  567. if (*rank < *n) {
  568. stzrqf_(rank, n, &a[a_offset], lda, &work[mn + 1], info);
  569. }
  570. /* Details of Householder rotations stored in WORK(MN+1:2*MN) */
  571. /* B(1:M,1:NRHS) := Q**T * B(1:M,1:NRHS) */
  572. sorm2r_("Left", "Transpose", m, nrhs, &mn, &a[a_offset], lda, &work[1], &
  573. b[b_offset], ldb, &work[(mn << 1) + 1], info);
  574. /* workspace NRHS */
  575. /* B(1:RANK,1:NRHS) := inv(T11) * B(1:RANK,1:NRHS) */
  576. strsm_("Left", "Upper", "No transpose", "Non-unit", rank, nrhs, &c_b36, &
  577. a[a_offset], lda, &b[b_offset], ldb);
  578. i__1 = *n;
  579. for (i__ = *rank + 1; i__ <= i__1; ++i__) {
  580. i__2 = *nrhs;
  581. for (j = 1; j <= i__2; ++j) {
  582. b[i__ + j * b_dim1] = 0.f;
  583. /* L30: */
  584. }
  585. /* L40: */
  586. }
  587. /* B(1:N,1:NRHS) := Y**T * B(1:N,1:NRHS) */
  588. if (*rank < *n) {
  589. i__1 = *rank;
  590. for (i__ = 1; i__ <= i__1; ++i__) {
  591. i__2 = *n - *rank + 1;
  592. slatzm_("Left", &i__2, nrhs, &a[i__ + (*rank + 1) * a_dim1], lda,
  593. &work[mn + i__], &b[i__ + b_dim1], &b[*rank + 1 + b_dim1],
  594. ldb, &work[(mn << 1) + 1]);
  595. /* L50: */
  596. }
  597. }
  598. /* workspace NRHS */
  599. /* B(1:N,1:NRHS) := P * B(1:N,1:NRHS) */
  600. i__1 = *nrhs;
  601. for (j = 1; j <= i__1; ++j) {
  602. i__2 = *n;
  603. for (i__ = 1; i__ <= i__2; ++i__) {
  604. work[(mn << 1) + i__] = 1.f;
  605. /* L60: */
  606. }
  607. i__2 = *n;
  608. for (i__ = 1; i__ <= i__2; ++i__) {
  609. if (work[(mn << 1) + i__] == 1.f) {
  610. if (jpvt[i__] != i__) {
  611. k = i__;
  612. t1 = b[k + j * b_dim1];
  613. t2 = b[jpvt[k] + j * b_dim1];
  614. L70:
  615. b[jpvt[k] + j * b_dim1] = t1;
  616. work[(mn << 1) + k] = 0.f;
  617. t1 = t2;
  618. k = jpvt[k];
  619. t2 = b[jpvt[k] + j * b_dim1];
  620. if (jpvt[k] != i__) {
  621. goto L70;
  622. }
  623. b[i__ + j * b_dim1] = t1;
  624. work[(mn << 1) + k] = 0.f;
  625. }
  626. }
  627. /* L80: */
  628. }
  629. /* L90: */
  630. }
  631. /* Undo scaling */
  632. if (iascl == 1) {
  633. slascl_("G", &c__0, &c__0, &anrm, &smlnum, n, nrhs, &b[b_offset], ldb,
  634. info);
  635. slascl_("U", &c__0, &c__0, &smlnum, &anrm, rank, rank, &a[a_offset],
  636. lda, info);
  637. } else if (iascl == 2) {
  638. slascl_("G", &c__0, &c__0, &anrm, &bignum, n, nrhs, &b[b_offset], ldb,
  639. info);
  640. slascl_("U", &c__0, &c__0, &bignum, &anrm, rank, rank, &a[a_offset],
  641. lda, info);
  642. }
  643. if (ibscl == 1) {
  644. slascl_("G", &c__0, &c__0, &smlnum, &bnrm, n, nrhs, &b[b_offset], ldb,
  645. info);
  646. } else if (ibscl == 2) {
  647. slascl_("G", &c__0, &c__0, &bignum, &bnrm, n, nrhs, &b[b_offset], ldb,
  648. info);
  649. }
  650. L100:
  651. return;
  652. /* End of SGELSX */
  653. } /* sgelsx_ */