C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / cgetsls.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 complex c_b1 = {0.f,0.f};
516 static integer c_n1 = -1;
517 static integer c_n2 = -2;
518 static integer c__0 = 0;
519
520 /* > \brief \b CGETSLS */
521
522 /*  Definition: */
523 /*  =========== */
524
525 /*       SUBROUTINE CGETSLS( TRANS, M, N, NRHS, A, LDA, B, LDB, */
526 /*     $                     WORK, LWORK, INFO ) */
527
528 /*       CHARACTER          TRANS */
529 /*       INTEGER            INFO, LDA, LDB, LWORK, M, N, NRHS */
530 /*       COMPLEX            A( LDA, * ), B( LDB, * ), WORK( * ) */
531
532
533 /* > \par Purpose: */
534 /*  ============= */
535 /* > */
536 /* > \verbatim */
537 /* > */
538 /* > CGETSLS solves overdetermined or underdetermined complex linear systems */
539 /* > involving an M-by-N matrix A, using a tall skinny QR or short wide LQ */
540 /* > factorization of A.  It is assumed that A has full rank. */
541 /* > */
542 /* > */
543 /* > */
544 /* > The following options are provided: */
545 /* > */
546 /* > 1. If TRANS = 'N' and m >= n:  find the least squares solution of */
547 /* >    an overdetermined system, i.e., solve the least squares problem */
548 /* >                 minimize || B - A*X ||. */
549 /* > */
550 /* > 2. If TRANS = 'N' and m < n:  find the minimum norm solution of */
551 /* >    an underdetermined system A * X = B. */
552 /* > */
553 /* > 3. If TRANS = 'C' and m >= n:  find the minimum norm solution of */
554 /* >    an undetermined system A**T * X = B. */
555 /* > */
556 /* > 4. If TRANS = 'C' and m < n:  find the least squares solution of */
557 /* >    an overdetermined system, i.e., solve the least squares problem */
558 /* >                 minimize || B - A**T * X ||. */
559 /* > */
560 /* > Several right hand side vectors b and solution vectors x can be */
561 /* > handled in a single call; they are stored as the columns of the */
562 /* > M-by-NRHS right hand side matrix B and the N-by-NRHS solution */
563 /* > matrix X. */
564 /* > \endverbatim */
565
566 /*  Arguments: */
567 /*  ========== */
568
569 /* > \param[in] TRANS */
570 /* > \verbatim */
571 /* >          TRANS is CHARACTER*1 */
572 /* >          = 'N': the linear system involves A; */
573 /* >          = 'C': the linear system involves A**H. */
574 /* > \endverbatim */
575 /* > */
576 /* > \param[in] M */
577 /* > \verbatim */
578 /* >          M is INTEGER */
579 /* >          The number of rows of the matrix A.  M >= 0. */
580 /* > \endverbatim */
581 /* > */
582 /* > \param[in] N */
583 /* > \verbatim */
584 /* >          N is INTEGER */
585 /* >          The number of columns of the matrix A.  N >= 0. */
586 /* > \endverbatim */
587 /* > */
588 /* > \param[in] NRHS */
589 /* > \verbatim */
590 /* >          NRHS is INTEGER */
591 /* >          The number of right hand sides, i.e., the number of */
592 /* >          columns of the matrices B and X. NRHS >=0. */
593 /* > \endverbatim */
594 /* > */
595 /* > \param[in,out] A */
596 /* > \verbatim */
597 /* >          A is COMPLEX array, dimension (LDA,N) */
598 /* >          On entry, the M-by-N matrix A. */
599 /* >          On exit, */
600 /* >          A is overwritten by details of its QR or LQ */
601 /* >          factorization as returned by CGEQR or CGELQ. */
602 /* > \endverbatim */
603 /* > */
604 /* > \param[in] LDA */
605 /* > \verbatim */
606 /* >          LDA is INTEGER */
607 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,M). */
608 /* > \endverbatim */
609 /* > */
610 /* > \param[in,out] B */
611 /* > \verbatim */
612 /* >          B is COMPLEX array, dimension (LDB,NRHS) */
613 /* >          On entry, the matrix B of right hand side vectors, stored */
614 /* >          columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS */
615 /* >          if TRANS = 'C'. */
616 /* >          On exit, if INFO = 0, B is overwritten by the solution */
617 /* >          vectors, stored columnwise: */
618 /* >          if TRANS = 'N' and m >= n, rows 1 to n of B contain the least */
619 /* >          squares solution vectors. */
620 /* >          if TRANS = 'N' and m < n, rows 1 to N of B contain the */
621 /* >          minimum norm solution vectors; */
622 /* >          if TRANS = 'C' and m >= n, rows 1 to M of B contain the */
623 /* >          minimum norm solution vectors; */
624 /* >          if TRANS = 'C' and m < n, rows 1 to M of B contain the */
625 /* >          least squares solution vectors. */
626 /* > \endverbatim */
627 /* > */
628 /* > \param[in] LDB */
629 /* > \verbatim */
630 /* >          LDB is INTEGER */
631 /* >          The leading dimension of the array B. LDB >= MAX(1,M,N). */
632 /* > \endverbatim */
633 /* > */
634 /* > \param[out] WORK */
635 /* > \verbatim */
636 /* >          (workspace) COMPLEX array, dimension (MAX(1,LWORK)) */
637 /* >          On exit, if INFO = 0, WORK(1) contains optimal (or either minimal */
638 /* >          or optimal, if query was assumed) LWORK. */
639 /* >          See LWORK for details. */
640 /* > \endverbatim */
641 /* > */
642 /* > \param[in] LWORK */
643 /* > \verbatim */
644 /* >          LWORK is INTEGER */
645 /* >          The dimension of the array WORK. */
646 /* >          If LWORK = -1 or -2, then a workspace query is assumed. */
647 /* >          If LWORK = -1, the routine calculates optimal size of WORK for the */
648 /* >          optimal performance and returns this value in WORK(1). */
649 /* >          If LWORK = -2, the routine calculates minimal size of WORK and */
650 /* >          returns this value in WORK(1). */
651 /* > \endverbatim */
652 /* > */
653 /* > \param[out] INFO */
654 /* > \verbatim */
655 /* >          INFO is INTEGER */
656 /* >          = 0:  successful exit */
657 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
658 /* >          > 0:  if INFO =  i, the i-th diagonal element of the */
659 /* >                triangular factor of A is zero, so that A does not have */
660 /* >                full rank; the least squares solution could not be */
661 /* >                computed. */
662 /* > \endverbatim */
663
664 /*  Authors: */
665 /*  ======== */
666
667 /* > \author Univ. of Tennessee */
668 /* > \author Univ. of California Berkeley */
669 /* > \author Univ. of Colorado Denver */
670 /* > \author NAG Ltd. */
671
672 /* > \date June 2017 */
673
674 /* > \ingroup complexGEsolve */
675
676 /*  ===================================================================== */
677 /* Subroutine */ int cgetsls_(char *trans, integer *m, integer *n, integer *
678         nrhs, complex *a, integer *lda, complex *b, integer *ldb, complex *
679         work, integer *lwork, integer *info)
680 {
681     /* System generated locals */
682     integer a_dim1, a_offset, b_dim1, b_offset, i__1, i__2, i__3;
683     real r__1;
684
685     /* Local variables */
686     real anrm, bnrm;
687     logical tran;
688     integer brow, tszm, tszo, info2, i__, j, iascl, ibscl;
689     extern /* Subroutine */ int cgelq_(integer *, integer *, complex *, 
690             integer *, complex *, integer *, complex *, integer *, integer *);
691     extern logical lsame_(char *, char *);
692     extern /* Subroutine */ int cgeqr_(integer *, integer *, complex *, 
693             integer *, complex *, integer *, complex *, integer *, integer *);
694     integer minmn, maxmn;
695     complex workq[1];
696     extern /* Subroutine */ int slabad_(real *, real *);
697     extern real clange_(char *, integer *, integer *, complex *, integer *, 
698             real *);
699     extern /* Subroutine */ int clascl_(char *, integer *, integer *, real *, 
700             real *, integer *, integer *, complex *, integer *, integer *);
701     complex tq[5];
702     extern real slamch_(char *);
703     extern /* Subroutine */ int cgemlq_(char *, char *, integer *, integer *, 
704             integer *, complex *, integer *, complex *, integer *, complex *, 
705             integer *, complex *, integer *, integer *), 
706             claset_(char *, integer *, integer *, complex *, complex *, 
707             complex *, integer *), xerbla_(char *, integer *, ftnlen),
708              cgemqr_(char *, char *, integer *, integer *, integer *, complex 
709             *, integer *, complex *, integer *, complex *, integer *, complex 
710             *, integer *, integer *);
711     integer scllen;
712     real bignum, smlnum;
713     integer wsizem, wsizeo;
714     logical lquery;
715     extern /* Subroutine */ int ctrtrs_(char *, char *, char *, integer *, 
716             integer *, complex *, integer *, complex *, integer *, integer *);
717     integer lw1, lw2, mnk;
718     real dum[1];
719     integer lwm, lwo;
720
721
722 /*  -- LAPACK driver routine (version 3.7.1) -- */
723 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
724 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
725 /*     June 2017 */
726
727
728
729 /*  ===================================================================== */
730
731
732 /*     Test the input arguments. */
733
734     /* Parameter adjustments */
735     a_dim1 = *lda;
736     a_offset = 1 + a_dim1 * 1;
737     a -= a_offset;
738     b_dim1 = *ldb;
739     b_offset = 1 + b_dim1 * 1;
740     b -= b_offset;
741     --work;
742
743     /* Function Body */
744     *info = 0;
745     minmn = f2cmin(*m,*n);
746     maxmn = f2cmax(*m,*n);
747     mnk = f2cmax(minmn,*nrhs);
748     tran = lsame_(trans, "C");
749
750     lquery = *lwork == -1 || *lwork == -2;
751     if (! (lsame_(trans, "N") || lsame_(trans, "C"))) {
752         *info = -1;
753     } else if (*m < 0) {
754         *info = -2;
755     } else if (*n < 0) {
756         *info = -3;
757     } else if (*nrhs < 0) {
758         *info = -4;
759     } else if (*lda < f2cmax(1,*m)) {
760         *info = -6;
761     } else /* if(complicated condition) */ {
762 /* Computing MAX */
763         i__1 = f2cmax(1,*m);
764         if (*ldb < f2cmax(i__1,*n)) {
765             *info = -8;
766         }
767     }
768
769     if (*info == 0) {
770
771 /*     Determine the block size and minimum LWORK */
772
773         if (*m >= *n) {
774             cgeqr_(m, n, &a[a_offset], lda, tq, &c_n1, workq, &c_n1, &info2);
775             tszo = (integer) tq[0].r;
776             lwo = (integer) workq[0].r;
777             cgemqr_("L", trans, m, nrhs, n, &a[a_offset], lda, tq, &tszo, &b[
778                     b_offset], ldb, workq, &c_n1, &info2);
779 /* Computing MAX */
780             i__1 = lwo, i__2 = (integer) workq[0].r;
781             lwo = f2cmax(i__1,i__2);
782             cgeqr_(m, n, &a[a_offset], lda, tq, &c_n2, workq, &c_n2, &info2);
783             tszm = (integer) tq[0].r;
784             lwm = (integer) workq[0].r;
785             cgemqr_("L", trans, m, nrhs, n, &a[a_offset], lda, tq, &tszm, &b[
786                     b_offset], ldb, workq, &c_n1, &info2);
787 /* Computing MAX */
788             i__1 = lwm, i__2 = (integer) workq[0].r;
789             lwm = f2cmax(i__1,i__2);
790             wsizeo = tszo + lwo;
791             wsizem = tszm + lwm;
792         } else {
793             cgelq_(m, n, &a[a_offset], lda, tq, &c_n1, workq, &c_n1, &info2);
794             tszo = (integer) tq[0].r;
795             lwo = (integer) workq[0].r;
796             cgemlq_("L", trans, n, nrhs, m, &a[a_offset], lda, tq, &tszo, &b[
797                     b_offset], ldb, workq, &c_n1, &info2);
798 /* Computing MAX */
799             i__1 = lwo, i__2 = (integer) workq[0].r;
800             lwo = f2cmax(i__1,i__2);
801             cgelq_(m, n, &a[a_offset], lda, tq, &c_n2, workq, &c_n2, &info2);
802             tszm = (integer) tq[0].r;
803             lwm = (integer) workq[0].r;
804             cgemlq_("L", trans, n, nrhs, m, &a[a_offset], lda, tq, &tszm, &b[
805                     b_offset], ldb, workq, &c_n1, &info2);
806 /* Computing MAX */
807             i__1 = lwm, i__2 = (integer) workq[0].r;
808             lwm = f2cmax(i__1,i__2);
809             wsizeo = tszo + lwo;
810             wsizem = tszm + lwm;
811         }
812
813         if (*lwork < wsizem && ! lquery) {
814             *info = -10;
815         }
816
817     }
818
819     if (*info != 0) {
820         i__1 = -(*info);
821         xerbla_("CGETSLS", &i__1, (ftnlen)7);
822         r__1 = (real) wsizeo;
823         work[1].r = r__1, work[1].i = 0.f;
824         return 0;
825     }
826     if (lquery) {
827         if (*lwork == -1) {
828             r__1 = (real) wsizeo;
829             work[1].r = r__1, work[1].i = 0.f;
830         }
831         if (*lwork == -2) {
832             r__1 = (real) wsizem;
833             work[1].r = r__1, work[1].i = 0.f;
834         }
835         return 0;
836     }
837     if (*lwork < wsizeo) {
838         lw1 = tszm;
839         lw2 = lwm;
840     } else {
841         lw1 = tszo;
842         lw2 = lwo;
843     }
844
845 /*     Quick return if possible */
846
847 /* Computing MIN */
848     i__1 = f2cmin(*m,*n);
849     if (f2cmin(i__1,*nrhs) == 0) {
850         i__1 = f2cmax(*m,*n);
851         claset_("FULL", &i__1, nrhs, &c_b1, &c_b1, &b[b_offset], ldb);
852         return 0;
853     }
854
855 /*     Get machine parameters */
856
857     smlnum = slamch_("S") / slamch_("P");
858     bignum = 1.f / smlnum;
859     slabad_(&smlnum, &bignum);
860
861 /*     Scale A, B if f2cmax element outside range [SMLNUM,BIGNUM] */
862
863     anrm = clange_("M", m, n, &a[a_offset], lda, dum);
864     iascl = 0;
865     if (anrm > 0.f && anrm < smlnum) {
866
867 /*        Scale matrix norm up to SMLNUM */
868
869         clascl_("G", &c__0, &c__0, &anrm, &smlnum, m, n, &a[a_offset], lda, 
870                 info);
871         iascl = 1;
872     } else if (anrm > bignum) {
873
874 /*        Scale matrix norm down to BIGNUM */
875
876         clascl_("G", &c__0, &c__0, &anrm, &bignum, m, n, &a[a_offset], lda, 
877                 info);
878         iascl = 2;
879     } else if (anrm == 0.f) {
880
881 /*        Matrix all zero. Return zero solution. */
882
883         claset_("F", &maxmn, nrhs, &c_b1, &c_b1, &b[b_offset], ldb)
884                 ;
885         goto L50;
886     }
887
888     brow = *m;
889     if (tran) {
890         brow = *n;
891     }
892     bnrm = clange_("M", &brow, nrhs, &b[b_offset], ldb, dum);
893     ibscl = 0;
894     if (bnrm > 0.f && bnrm < smlnum) {
895
896 /*        Scale matrix norm up to SMLNUM */
897
898         clascl_("G", &c__0, &c__0, &bnrm, &smlnum, &brow, nrhs, &b[b_offset], 
899                 ldb, info);
900         ibscl = 1;
901     } else if (bnrm > bignum) {
902
903 /*        Scale matrix norm down to BIGNUM */
904
905         clascl_("G", &c__0, &c__0, &bnrm, &bignum, &brow, nrhs, &b[b_offset], 
906                 ldb, info);
907         ibscl = 2;
908     }
909
910     if (*m >= *n) {
911
912 /*        compute QR factorization of A */
913
914         cgeqr_(m, n, &a[a_offset], lda, &work[lw2 + 1], &lw1, &work[1], &lw2, 
915                 info);
916         if (! tran) {
917
918 /*           Least-Squares Problem f2cmin || A * X - B || */
919
920 /*           B(1:M,1:NRHS) := Q**T * B(1:M,1:NRHS) */
921
922             cgemqr_("L", "C", m, nrhs, n, &a[a_offset], lda, &work[lw2 + 1], &
923                     lw1, &b[b_offset], ldb, &work[1], &lw2, info);
924
925 /*           B(1:N,1:NRHS) := inv(R) * B(1:N,1:NRHS) */
926
927             ctrtrs_("U", "N", "N", n, nrhs, &a[a_offset], lda, &b[b_offset], 
928                     ldb, info);
929             if (*info > 0) {
930                 return 0;
931             }
932             scllen = *n;
933         } else {
934
935 /*           Overdetermined system of equations A**T * X = B */
936
937 /*           B(1:N,1:NRHS) := inv(R**T) * B(1:N,1:NRHS) */
938
939             ctrtrs_("U", "C", "N", n, nrhs, &a[a_offset], lda, &b[b_offset], 
940                     ldb, info);
941
942             if (*info > 0) {
943                 return 0;
944             }
945
946 /*           B(N+1:M,1:NRHS) = CZERO */
947
948             i__1 = *nrhs;
949             for (j = 1; j <= i__1; ++j) {
950                 i__2 = *m;
951                 for (i__ = *n + 1; i__ <= i__2; ++i__) {
952                     i__3 = i__ + j * b_dim1;
953                     b[i__3].r = 0.f, b[i__3].i = 0.f;
954 /* L10: */
955                 }
956 /* L20: */
957             }
958
959 /*           B(1:M,1:NRHS) := Q(1:N,:) * B(1:N,1:NRHS) */
960
961             cgemqr_("L", "N", m, nrhs, n, &a[a_offset], lda, &work[lw2 + 1], &
962                     lw1, &b[b_offset], ldb, &work[1], &lw2, info);
963
964             scllen = *m;
965
966         }
967
968     } else {
969
970 /*        Compute LQ factorization of A */
971
972         cgelq_(m, n, &a[a_offset], lda, &work[lw2 + 1], &lw1, &work[1], &lw2, 
973                 info);
974
975 /*        workspace at least M, optimally M*NB. */
976
977         if (! tran) {
978
979 /*           underdetermined system of equations A * X = B */
980
981 /*           B(1:M,1:NRHS) := inv(L) * B(1:M,1:NRHS) */
982
983             ctrtrs_("L", "N", "N", m, nrhs, &a[a_offset], lda, &b[b_offset], 
984                     ldb, info);
985
986             if (*info > 0) {
987                 return 0;
988             }
989
990 /*           B(M+1:N,1:NRHS) = 0 */
991
992             i__1 = *nrhs;
993             for (j = 1; j <= i__1; ++j) {
994                 i__2 = *n;
995                 for (i__ = *m + 1; i__ <= i__2; ++i__) {
996                     i__3 = i__ + j * b_dim1;
997                     b[i__3].r = 0.f, b[i__3].i = 0.f;
998 /* L30: */
999                 }
1000 /* L40: */
1001             }
1002
1003 /*           B(1:N,1:NRHS) := Q(1:N,:)**T * B(1:M,1:NRHS) */
1004
1005             cgemlq_("L", "C", n, nrhs, m, &a[a_offset], lda, &work[lw2 + 1], &
1006                     lw1, &b[b_offset], ldb, &work[1], &lw2, info);
1007
1008 /*           workspace at least NRHS, optimally NRHS*NB */
1009
1010             scllen = *n;
1011
1012         } else {
1013
1014 /*           overdetermined system f2cmin || A**T * X - B || */
1015
1016 /*           B(1:N,1:NRHS) := Q * B(1:N,1:NRHS) */
1017
1018             cgemlq_("L", "N", n, nrhs, m, &a[a_offset], lda, &work[lw2 + 1], &
1019                     lw1, &b[b_offset], ldb, &work[1], &lw2, info);
1020
1021 /*           workspace at least NRHS, optimally NRHS*NB */
1022
1023 /*           B(1:M,1:NRHS) := inv(L**T) * B(1:M,1:NRHS) */
1024
1025             ctrtrs_("L", "C", "N", m, nrhs, &a[a_offset], lda, &b[b_offset], 
1026                     ldb, info);
1027
1028             if (*info > 0) {
1029                 return 0;
1030             }
1031
1032             scllen = *m;
1033
1034         }
1035
1036     }
1037
1038 /*     Undo scaling */
1039
1040     if (iascl == 1) {
1041         clascl_("G", &c__0, &c__0, &anrm, &smlnum, &scllen, nrhs, &b[b_offset]
1042                 , ldb, info);
1043     } else if (iascl == 2) {
1044         clascl_("G", &c__0, &c__0, &anrm, &bignum, &scllen, nrhs, &b[b_offset]
1045                 , ldb, info);
1046     }
1047     if (ibscl == 1) {
1048         clascl_("G", &c__0, &c__0, &smlnum, &bnrm, &scllen, nrhs, &b[b_offset]
1049                 , ldb, info);
1050     } else if (ibscl == 2) {
1051         clascl_("G", &c__0, &c__0, &bignum, &bnrm, &scllen, nrhs, &b[b_offset]
1052                 , ldb, info);
1053     }
1054
1055 L50:
1056     r__1 = (real) (tszo + lwo);
1057     work[1].r = r__1, work[1].i = 0.f;
1058     return 0;
1059
1060 /*     End of ZGETSLS */
1061
1062 } /* cgetsls_ */
1063