C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / sgels.c
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
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #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]);}
194 #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]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #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++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #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]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static integer c__1 = 1;
516 static integer c_n1 = -1;
517 static real c_b33 = 0.f;
518 static integer c__0 = 0;
519
520 /* > \brief <b> SGELS solves overdetermined or underdetermined systems for GE matrices</b> */
521
522 /*  =========== DOCUMENTATION =========== */
523
524 /* Online html documentation available at */
525 /*            http://www.netlib.org/lapack/explore-html/ */
526
527 /* > \htmlonly */
528 /* > Download SGELS + dependencies */
529 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgels.f
530 "> */
531 /* > [TGZ]</a> */
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgels.f
533 "> */
534 /* > [ZIP]</a> */
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgels.f
536 "> */
537 /* > [TXT]</a> */
538 /* > \endhtmlonly */
539
540 /*  Definition: */
541 /*  =========== */
542
543 /*       SUBROUTINE SGELS( TRANS, M, N, NRHS, A, LDA, B, LDB, WORK, LWORK, */
544 /*                         INFO ) */
545
546 /*       CHARACTER          TRANS */
547 /*       INTEGER            INFO, LDA, LDB, LWORK, M, N, NRHS */
548 /*       REAL               A( LDA, * ), B( LDB, * ), WORK( * ) */
549
550
551 /* > \par Purpose: */
552 /*  ============= */
553 /* > */
554 /* > \verbatim */
555 /* > */
556 /* > SGELS solves overdetermined or underdetermined real linear systems */
557 /* > involving an M-by-N matrix A, or its transpose, using a QR or LQ */
558 /* > factorization of A.  It is assumed that A has full rank. */
559 /* > */
560 /* > The following options are provided: */
561 /* > */
562 /* > 1. If TRANS = 'N' and m >= n:  find the least squares solution of */
563 /* >    an overdetermined system, i.e., solve the least squares problem */
564 /* >                 minimize || B - A*X ||. */
565 /* > */
566 /* > 2. If TRANS = 'N' and m < n:  find the minimum norm solution of */
567 /* >    an underdetermined system A * X = B. */
568 /* > */
569 /* > 3. If TRANS = 'T' and m >= n:  find the minimum norm solution of */
570 /* >    an underdetermined system A**T * X = B. */
571 /* > */
572 /* > 4. If TRANS = 'T' and m < n:  find the least squares solution of */
573 /* >    an overdetermined system, i.e., solve the least squares problem */
574 /* >                 minimize || B - A**T * X ||. */
575 /* > */
576 /* > Several right hand side vectors b and solution vectors x can be */
577 /* > handled in a single call; they are stored as the columns of the */
578 /* > M-by-NRHS right hand side matrix B and the N-by-NRHS solution */
579 /* > matrix X. */
580 /* > \endverbatim */
581
582 /*  Arguments: */
583 /*  ========== */
584
585 /* > \param[in] TRANS */
586 /* > \verbatim */
587 /* >          TRANS is CHARACTER*1 */
588 /* >          = 'N': the linear system involves A; */
589 /* >          = 'T': the linear system involves A**T. */
590 /* > \endverbatim */
591 /* > */
592 /* > \param[in] M */
593 /* > \verbatim */
594 /* >          M is INTEGER */
595 /* >          The number of rows of the matrix A.  M >= 0. */
596 /* > \endverbatim */
597 /* > */
598 /* > \param[in] N */
599 /* > \verbatim */
600 /* >          N is INTEGER */
601 /* >          The number of columns of the matrix A.  N >= 0. */
602 /* > \endverbatim */
603 /* > */
604 /* > \param[in] NRHS */
605 /* > \verbatim */
606 /* >          NRHS is INTEGER */
607 /* >          The number of right hand sides, i.e., the number of */
608 /* >          columns of the matrices B and X. NRHS >=0. */
609 /* > \endverbatim */
610 /* > */
611 /* > \param[in,out] A */
612 /* > \verbatim */
613 /* >          A is REAL array, dimension (LDA,N) */
614 /* >          On entry, the M-by-N matrix A. */
615 /* >          On exit, */
616 /* >            if M >= N, A is overwritten by details of its QR */
617 /* >                       factorization as returned by SGEQRF; */
618 /* >            if M <  N, A is overwritten by details of its LQ */
619 /* >                       factorization as returned by SGELQF. */
620 /* > \endverbatim */
621 /* > */
622 /* > \param[in] LDA */
623 /* > \verbatim */
624 /* >          LDA is INTEGER */
625 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,M). */
626 /* > \endverbatim */
627 /* > */
628 /* > \param[in,out] B */
629 /* > \verbatim */
630 /* >          B is REAL array, dimension (LDB,NRHS) */
631 /* >          On entry, the matrix B of right hand side vectors, stored */
632 /* >          columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS */
633 /* >          if TRANS = 'T'. */
634 /* >          On exit, if INFO = 0, B is overwritten by the solution */
635 /* >          vectors, stored columnwise: */
636 /* >          if TRANS = 'N' and m >= n, rows 1 to n of B contain the least */
637 /* >          squares solution vectors; the residual sum of squares for the */
638 /* >          solution in each column is given by the sum of squares of */
639 /* >          elements N+1 to M in that column; */
640 /* >          if TRANS = 'N' and m < n, rows 1 to N of B contain the */
641 /* >          minimum norm solution vectors; */
642 /* >          if TRANS = 'T' and m >= n, rows 1 to M of B contain the */
643 /* >          minimum norm solution vectors; */
644 /* >          if TRANS = 'T' and m < n, rows 1 to M of B contain the */
645 /* >          least squares solution vectors; the residual sum of squares */
646 /* >          for the solution in each column is given by the sum of */
647 /* >          squares of elements M+1 to N in that column. */
648 /* > \endverbatim */
649 /* > */
650 /* > \param[in] LDB */
651 /* > \verbatim */
652 /* >          LDB is INTEGER */
653 /* >          The leading dimension of the array B. LDB >= MAX(1,M,N). */
654 /* > \endverbatim */
655 /* > */
656 /* > \param[out] WORK */
657 /* > \verbatim */
658 /* >          WORK is REAL array, dimension (MAX(1,LWORK)) */
659 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
660 /* > \endverbatim */
661 /* > */
662 /* > \param[in] LWORK */
663 /* > \verbatim */
664 /* >          LWORK is INTEGER */
665 /* >          The dimension of the array WORK. */
666 /* >          LWORK >= f2cmax( 1, MN + f2cmax( MN, NRHS ) ). */
667 /* >          For optimal performance, */
668 /* >          LWORK >= f2cmax( 1, MN + f2cmax( MN, NRHS )*NB ). */
669 /* >          where MN = f2cmin(M,N) and NB is the optimum block size. */
670 /* > */
671 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
672 /* >          only calculates the optimal size of the WORK array, returns */
673 /* >          this value as the first entry of the WORK array, and no error */
674 /* >          message related to LWORK is issued by XERBLA. */
675 /* > \endverbatim */
676 /* > */
677 /* > \param[out] INFO */
678 /* > \verbatim */
679 /* >          INFO is INTEGER */
680 /* >          = 0:  successful exit */
681 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
682 /* >          > 0:  if INFO =  i, the i-th diagonal element of the */
683 /* >                triangular factor of A is zero, so that A does not have */
684 /* >                full rank; the least squares solution could not be */
685 /* >                computed. */
686 /* > \endverbatim */
687
688 /*  Authors: */
689 /*  ======== */
690
691 /* > \author Univ. of Tennessee */
692 /* > \author Univ. of California Berkeley */
693 /* > \author Univ. of Colorado Denver */
694 /* > \author NAG Ltd. */
695
696 /* > \date December 2016 */
697
698 /* > \ingroup realGEsolve */
699
700 /*  ===================================================================== */
701 /* Subroutine */ int sgels_(char *trans, integer *m, integer *n, integer *
702         nrhs, real *a, integer *lda, real *b, integer *ldb, real *work, 
703         integer *lwork, integer *info)
704 {
705     /* System generated locals */
706     integer a_dim1, a_offset, b_dim1, b_offset, i__1, i__2;
707
708     /* Local variables */
709     real anrm, bnrm;
710     integer brow;
711     logical tpsd;
712     integer i__, j, iascl, ibscl;
713     extern logical lsame_(char *, char *);
714     integer wsize;
715     real rwork[1];
716     integer nb;
717     extern /* Subroutine */ int slabad_(real *, real *);
718     integer mn;
719     extern real slamch_(char *), slange_(char *, integer *, integer *,
720              real *, integer *, real *);
721     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
722     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
723             integer *, integer *, ftnlen, ftnlen);
724     integer scllen;
725     real bignum;
726     extern /* Subroutine */ int sgelqf_(integer *, integer *, real *, integer 
727             *, real *, real *, integer *, integer *), slascl_(char *, integer 
728             *, integer *, real *, real *, integer *, integer *, real *, 
729             integer *, integer *), sgeqrf_(integer *, integer *, real 
730             *, integer *, real *, real *, integer *, integer *), slaset_(char 
731             *, integer *, integer *, real *, real *, real *, integer *);
732     real smlnum;
733     extern /* Subroutine */ int sormlq_(char *, char *, integer *, integer *, 
734             integer *, real *, integer *, real *, real *, integer *, real *, 
735             integer *, integer *);
736     logical lquery;
737     extern /* Subroutine */ int sormqr_(char *, char *, integer *, integer *, 
738             integer *, real *, integer *, real *, real *, integer *, real *, 
739             integer *, integer *), strtrs_(char *, char *, 
740             char *, integer *, integer *, real *, integer *, real *, integer *
741             , integer *);
742
743
744 /*  -- LAPACK driver routine (version 3.7.0) -- */
745 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
746 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
747 /*     December 2016 */
748
749
750 /*  ===================================================================== */
751
752
753 /*     Test the input arguments. */
754
755     /* Parameter adjustments */
756     a_dim1 = *lda;
757     a_offset = 1 + a_dim1 * 1;
758     a -= a_offset;
759     b_dim1 = *ldb;
760     b_offset = 1 + b_dim1 * 1;
761     b -= b_offset;
762     --work;
763
764     /* Function Body */
765     *info = 0;
766     mn = f2cmin(*m,*n);
767     lquery = *lwork == -1;
768     if (! (lsame_(trans, "N") || lsame_(trans, "T"))) {
769         *info = -1;
770     } else if (*m < 0) {
771         *info = -2;
772     } else if (*n < 0) {
773         *info = -3;
774     } else if (*nrhs < 0) {
775         *info = -4;
776     } else if (*lda < f2cmax(1,*m)) {
777         *info = -6;
778     } else /* if(complicated condition) */ {
779 /* Computing MAX */
780         i__1 = f2cmax(1,*m);
781         if (*ldb < f2cmax(i__1,*n)) {
782             *info = -8;
783         } else /* if(complicated condition) */ {
784 /* Computing MAX */
785             i__1 = 1, i__2 = mn + f2cmax(mn,*nrhs);
786             if (*lwork < f2cmax(i__1,i__2) && ! lquery) {
787                 *info = -10;
788             }
789         }
790     }
791
792 /*     Figure out optimal block size */
793
794     if (*info == 0 || *info == -10) {
795
796         tpsd = TRUE_;
797         if (lsame_(trans, "N")) {
798             tpsd = FALSE_;
799         }
800
801         if (*m >= *n) {
802             nb = ilaenv_(&c__1, "SGEQRF", " ", m, n, &c_n1, &c_n1, (ftnlen)6, 
803                     (ftnlen)1);
804             if (tpsd) {
805 /* Computing MAX */
806                 i__1 = nb, i__2 = ilaenv_(&c__1, "SORMQR", "LN", m, nrhs, n, &
807                         c_n1, (ftnlen)6, (ftnlen)2);
808                 nb = f2cmax(i__1,i__2);
809             } else {
810 /* Computing MAX */
811                 i__1 = nb, i__2 = ilaenv_(&c__1, "SORMQR", "LT", m, nrhs, n, &
812                         c_n1, (ftnlen)6, (ftnlen)2);
813                 nb = f2cmax(i__1,i__2);
814             }
815         } else {
816             nb = ilaenv_(&c__1, "SGELQF", " ", m, n, &c_n1, &c_n1, (ftnlen)6, 
817                     (ftnlen)1);
818             if (tpsd) {
819 /* Computing MAX */
820                 i__1 = nb, i__2 = ilaenv_(&c__1, "SORMLQ", "LT", n, nrhs, m, &
821                         c_n1, (ftnlen)6, (ftnlen)2);
822                 nb = f2cmax(i__1,i__2);
823             } else {
824 /* Computing MAX */
825                 i__1 = nb, i__2 = ilaenv_(&c__1, "SORMLQ", "LN", n, nrhs, m, &
826                         c_n1, (ftnlen)6, (ftnlen)2);
827                 nb = f2cmax(i__1,i__2);
828             }
829         }
830
831 /* Computing MAX */
832         i__1 = 1, i__2 = mn + f2cmax(mn,*nrhs) * nb;
833         wsize = f2cmax(i__1,i__2);
834         work[1] = (real) wsize;
835
836     }
837
838     if (*info != 0) {
839         i__1 = -(*info);
840         xerbla_("SGELS ", &i__1, (ftnlen)6);
841         return 0;
842     } else if (lquery) {
843         return 0;
844     }
845
846 /*     Quick return if possible */
847
848 /* Computing MIN */
849     i__1 = f2cmin(*m,*n);
850     if (f2cmin(i__1,*nrhs) == 0) {
851         i__1 = f2cmax(*m,*n);
852         slaset_("Full", &i__1, nrhs, &c_b33, &c_b33, &b[b_offset], ldb);
853         return 0;
854     }
855
856 /*     Get machine parameters */
857
858     smlnum = slamch_("S") / slamch_("P");
859     bignum = 1.f / smlnum;
860     slabad_(&smlnum, &bignum);
861
862 /*     Scale A, B if f2cmax element outside range [SMLNUM,BIGNUM] */
863
864     anrm = slange_("M", m, n, &a[a_offset], lda, rwork);
865     iascl = 0;
866     if (anrm > 0.f && anrm < smlnum) {
867
868 /*        Scale matrix norm up to SMLNUM */
869
870         slascl_("G", &c__0, &c__0, &anrm, &smlnum, m, n, &a[a_offset], lda, 
871                 info);
872         iascl = 1;
873     } else if (anrm > bignum) {
874
875 /*        Scale matrix norm down to BIGNUM */
876
877         slascl_("G", &c__0, &c__0, &anrm, &bignum, m, n, &a[a_offset], lda, 
878                 info);
879         iascl = 2;
880     } else if (anrm == 0.f) {
881
882 /*        Matrix all zero. Return zero solution. */
883
884         i__1 = f2cmax(*m,*n);
885         slaset_("F", &i__1, nrhs, &c_b33, &c_b33, &b[b_offset], ldb);
886         goto L50;
887     }
888
889     brow = *m;
890     if (tpsd) {
891         brow = *n;
892     }
893     bnrm = slange_("M", &brow, nrhs, &b[b_offset], ldb, rwork);
894     ibscl = 0;
895     if (bnrm > 0.f && bnrm < smlnum) {
896
897 /*        Scale matrix norm up to SMLNUM */
898
899         slascl_("G", &c__0, &c__0, &bnrm, &smlnum, &brow, nrhs, &b[b_offset], 
900                 ldb, info);
901         ibscl = 1;
902     } else if (bnrm > bignum) {
903
904 /*        Scale matrix norm down to BIGNUM */
905
906         slascl_("G", &c__0, &c__0, &bnrm, &bignum, &brow, nrhs, &b[b_offset], 
907                 ldb, info);
908         ibscl = 2;
909     }
910
911     if (*m >= *n) {
912
913 /*        compute QR factorization of A */
914
915         i__1 = *lwork - mn;
916         sgeqrf_(m, n, &a[a_offset], lda, &work[1], &work[mn + 1], &i__1, info)
917                 ;
918
919 /*        workspace at least N, optimally N*NB */
920
921         if (! tpsd) {
922
923 /*           Least-Squares Problem f2cmin || A * X - B || */
924
925 /*           B(1:M,1:NRHS) := Q**T * B(1:M,1:NRHS) */
926
927             i__1 = *lwork - mn;
928             sormqr_("Left", "Transpose", m, nrhs, n, &a[a_offset], lda, &work[
929                     1], &b[b_offset], ldb, &work[mn + 1], &i__1, info);
930
931 /*           workspace at least NRHS, optimally NRHS*NB */
932
933 /*           B(1:N,1:NRHS) := inv(R) * B(1:N,1:NRHS) */
934
935             strtrs_("Upper", "No transpose", "Non-unit", n, nrhs, &a[a_offset]
936                     , lda, &b[b_offset], ldb, info);
937
938             if (*info > 0) {
939                 return 0;
940             }
941
942             scllen = *n;
943
944         } else {
945
946 /*           Underdetermined system of equations A**T * X = B */
947
948 /*           B(1:N,1:NRHS) := inv(R**T) * B(1:N,1:NRHS) */
949
950             strtrs_("Upper", "Transpose", "Non-unit", n, nrhs, &a[a_offset], 
951                     lda, &b[b_offset], ldb, info);
952
953             if (*info > 0) {
954                 return 0;
955             }
956
957 /*           B(N+1:M,1:NRHS) = ZERO */
958
959             i__1 = *nrhs;
960             for (j = 1; j <= i__1; ++j) {
961                 i__2 = *m;
962                 for (i__ = *n + 1; i__ <= i__2; ++i__) {
963                     b[i__ + j * b_dim1] = 0.f;
964 /* L10: */
965                 }
966 /* L20: */
967             }
968
969 /*           B(1:M,1:NRHS) := Q(1:N,:) * B(1:N,1:NRHS) */
970
971             i__1 = *lwork - mn;
972             sormqr_("Left", "No transpose", m, nrhs, n, &a[a_offset], lda, &
973                     work[1], &b[b_offset], ldb, &work[mn + 1], &i__1, info);
974
975 /*           workspace at least NRHS, optimally NRHS*NB */
976
977             scllen = *m;
978
979         }
980
981     } else {
982
983 /*        Compute LQ factorization of A */
984
985         i__1 = *lwork - mn;
986         sgelqf_(m, n, &a[a_offset], lda, &work[1], &work[mn + 1], &i__1, info)
987                 ;
988
989 /*        workspace at least M, optimally M*NB. */
990
991         if (! tpsd) {
992
993 /*           underdetermined system of equations A * X = B */
994
995 /*           B(1:M,1:NRHS) := inv(L) * B(1:M,1:NRHS) */
996
997             strtrs_("Lower", "No transpose", "Non-unit", m, nrhs, &a[a_offset]
998                     , lda, &b[b_offset], ldb, info);
999
1000             if (*info > 0) {
1001                 return 0;
1002             }
1003
1004 /*           B(M+1:N,1:NRHS) = 0 */
1005
1006             i__1 = *nrhs;
1007             for (j = 1; j <= i__1; ++j) {
1008                 i__2 = *n;
1009                 for (i__ = *m + 1; i__ <= i__2; ++i__) {
1010                     b[i__ + j * b_dim1] = 0.f;
1011 /* L30: */
1012                 }
1013 /* L40: */
1014             }
1015
1016 /*           B(1:N,1:NRHS) := Q(1:N,:)**T * B(1:M,1:NRHS) */
1017
1018             i__1 = *lwork - mn;
1019             sormlq_("Left", "Transpose", n, nrhs, m, &a[a_offset], lda, &work[
1020                     1], &b[b_offset], ldb, &work[mn + 1], &i__1, info);
1021
1022 /*           workspace at least NRHS, optimally NRHS*NB */
1023
1024             scllen = *n;
1025
1026         } else {
1027
1028 /*           overdetermined system f2cmin || A**T * X - B || */
1029
1030 /*           B(1:N,1:NRHS) := Q * B(1:N,1:NRHS) */
1031
1032             i__1 = *lwork - mn;
1033             sormlq_("Left", "No transpose", n, nrhs, m, &a[a_offset], lda, &
1034                     work[1], &b[b_offset], ldb, &work[mn + 1], &i__1, info);
1035
1036 /*           workspace at least NRHS, optimally NRHS*NB */
1037
1038 /*           B(1:M,1:NRHS) := inv(L**T) * B(1:M,1:NRHS) */
1039
1040             strtrs_("Lower", "Transpose", "Non-unit", m, nrhs, &a[a_offset], 
1041                     lda, &b[b_offset], ldb, info);
1042
1043             if (*info > 0) {
1044                 return 0;
1045             }
1046
1047             scllen = *m;
1048
1049         }
1050
1051     }
1052
1053 /*     Undo scaling */
1054
1055     if (iascl == 1) {
1056         slascl_("G", &c__0, &c__0, &anrm, &smlnum, &scllen, nrhs, &b[b_offset]
1057                 , ldb, info);
1058     } else if (iascl == 2) {
1059         slascl_("G", &c__0, &c__0, &anrm, &bignum, &scllen, nrhs, &b[b_offset]
1060                 , ldb, info);
1061     }
1062     if (ibscl == 1) {
1063         slascl_("G", &c__0, &c__0, &smlnum, &bnrm, &scllen, nrhs, &b[b_offset]
1064                 , ldb, info);
1065     } else if (ibscl == 2) {
1066         slascl_("G", &c__0, &c__0, &bignum, &bnrm, &scllen, nrhs, &b[b_offset]
1067                 , ldb, info);
1068     }
1069
1070 L50:
1071     work[1] = (real) wsize;
1072
1073     return 0;
1074
1075 /*     End of SGELS */
1076
1077 } /* sgels_ */
1078