C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / zgelss.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 doublecomplex c_b1 = {0.,0.};
516 static doublecomplex c_b2 = {1.,0.};
517 static integer c__6 = 6;
518 static integer c_n1 = -1;
519 static integer c__1 = 1;
520 static integer c__0 = 0;
521 static doublereal c_b59 = 0.;
522
523 /* > \brief <b> ZGELSS solves overdetermined or underdetermined systems for GE matrices</b> */
524
525 /*  =========== DOCUMENTATION =========== */
526
527 /* Online html documentation available at */
528 /*            http://www.netlib.org/lapack/explore-html/ */
529
530 /* > \htmlonly */
531 /* > Download ZGELSS + dependencies */
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgelss.
533 f"> */
534 /* > [TGZ]</a> */
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgelss.
536 f"> */
537 /* > [ZIP]</a> */
538 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgelss.
539 f"> */
540 /* > [TXT]</a> */
541 /* > \endhtmlonly */
542
543 /*  Definition: */
544 /*  =========== */
545
546 /*       SUBROUTINE ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, */
547 /*                          WORK, LWORK, RWORK, INFO ) */
548
549 /*       INTEGER            INFO, LDA, LDB, LWORK, M, N, NRHS, RANK */
550 /*       DOUBLE PRECISION   RCOND */
551 /*       DOUBLE PRECISION   RWORK( * ), S( * ) */
552 /*       COMPLEX*16         A( LDA, * ), B( LDB, * ), WORK( * ) */
553
554
555 /* > \par Purpose: */
556 /*  ============= */
557 /* > */
558 /* > \verbatim */
559 /* > */
560 /* > ZGELSS computes the minimum norm solution to a complex linear */
561 /* > least squares problem: */
562 /* > */
563 /* > Minimize 2-norm(| b - A*x |). */
564 /* > */
565 /* > using the singular value decomposition (SVD) of A. A is an M-by-N */
566 /* > matrix which may be rank-deficient. */
567 /* > */
568 /* > Several right hand side vectors b and solution vectors x can be */
569 /* > handled in a single call; they are stored as the columns of the */
570 /* > M-by-NRHS right hand side matrix B and the N-by-NRHS solution matrix */
571 /* > X. */
572 /* > */
573 /* > The effective rank of A is determined by treating as zero those */
574 /* > singular values which are less than RCOND times the largest singular */
575 /* > value. */
576 /* > \endverbatim */
577
578 /*  Arguments: */
579 /*  ========== */
580
581 /* > \param[in] M */
582 /* > \verbatim */
583 /* >          M is INTEGER */
584 /* >          The number of rows of the matrix A. M >= 0. */
585 /* > \endverbatim */
586 /* > */
587 /* > \param[in] N */
588 /* > \verbatim */
589 /* >          N is INTEGER */
590 /* >          The number of columns of the matrix A. N >= 0. */
591 /* > \endverbatim */
592 /* > */
593 /* > \param[in] NRHS */
594 /* > \verbatim */
595 /* >          NRHS is INTEGER */
596 /* >          The number of right hand sides, i.e., the number of columns */
597 /* >          of the matrices B and X. NRHS >= 0. */
598 /* > \endverbatim */
599 /* > */
600 /* > \param[in,out] A */
601 /* > \verbatim */
602 /* >          A is COMPLEX*16 array, dimension (LDA,N) */
603 /* >          On entry, the M-by-N matrix A. */
604 /* >          On exit, the first f2cmin(m,n) rows of A are overwritten with */
605 /* >          its right singular vectors, stored rowwise. */
606 /* > \endverbatim */
607 /* > */
608 /* > \param[in] LDA */
609 /* > \verbatim */
610 /* >          LDA is INTEGER */
611 /* >          The leading dimension of the array A. LDA >= f2cmax(1,M). */
612 /* > \endverbatim */
613 /* > */
614 /* > \param[in,out] B */
615 /* > \verbatim */
616 /* >          B is COMPLEX*16 array, dimension (LDB,NRHS) */
617 /* >          On entry, the M-by-NRHS right hand side matrix B. */
618 /* >          On exit, B is overwritten by the N-by-NRHS solution matrix X. */
619 /* >          If m >= n and RANK = n, the residual sum-of-squares for */
620 /* >          the solution in the i-th column is given by the sum of */
621 /* >          squares of the modulus of elements n+1:m in that column. */
622 /* > \endverbatim */
623 /* > */
624 /* > \param[in] LDB */
625 /* > \verbatim */
626 /* >          LDB is INTEGER */
627 /* >          The leading dimension of the array B.  LDB >= f2cmax(1,M,N). */
628 /* > \endverbatim */
629 /* > */
630 /* > \param[out] S */
631 /* > \verbatim */
632 /* >          S is DOUBLE PRECISION array, dimension (f2cmin(M,N)) */
633 /* >          The singular values of A in decreasing order. */
634 /* >          The condition number of A in the 2-norm = S(1)/S(f2cmin(m,n)). */
635 /* > \endverbatim */
636 /* > */
637 /* > \param[in] RCOND */
638 /* > \verbatim */
639 /* >          RCOND is DOUBLE PRECISION */
640 /* >          RCOND is used to determine the effective rank of A. */
641 /* >          Singular values S(i) <= RCOND*S(1) are treated as zero. */
642 /* >          If RCOND < 0, machine precision is used instead. */
643 /* > \endverbatim */
644 /* > */
645 /* > \param[out] RANK */
646 /* > \verbatim */
647 /* >          RANK is INTEGER */
648 /* >          The effective rank of A, i.e., the number of singular values */
649 /* >          which are greater than RCOND*S(1). */
650 /* > \endverbatim */
651 /* > */
652 /* > \param[out] WORK */
653 /* > \verbatim */
654 /* >          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) */
655 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
656 /* > \endverbatim */
657 /* > */
658 /* > \param[in] LWORK */
659 /* > \verbatim */
660 /* >          LWORK is INTEGER */
661 /* >          The dimension of the array WORK. LWORK >= 1, and also: */
662 /* >          LWORK >=  2*f2cmin(M,N) + f2cmax(M,N,NRHS) */
663 /* >          For good performance, LWORK should generally be larger. */
664 /* > */
665 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
666 /* >          only calculates the optimal size of the WORK array, returns */
667 /* >          this value as the first entry of the WORK array, and no error */
668 /* >          message related to LWORK is issued by XERBLA. */
669 /* > \endverbatim */
670 /* > */
671 /* > \param[out] RWORK */
672 /* > \verbatim */
673 /* >          RWORK is DOUBLE PRECISION array, dimension (5*f2cmin(M,N)) */
674 /* > \endverbatim */
675 /* > */
676 /* > \param[out] INFO */
677 /* > \verbatim */
678 /* >          INFO is INTEGER */
679 /* >          = 0:  successful exit */
680 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
681 /* >          > 0:  the algorithm for computing the SVD failed to converge; */
682 /* >                if INFO = i, i off-diagonal elements of an intermediate */
683 /* >                bidiagonal form did not converge to zero. */
684 /* > \endverbatim */
685
686 /*  Authors: */
687 /*  ======== */
688
689 /* > \author Univ. of Tennessee */
690 /* > \author Univ. of California Berkeley */
691 /* > \author Univ. of Colorado Denver */
692 /* > \author NAG Ltd. */
693
694 /* > \date June 2016 */
695
696 /* > \ingroup complex16GEsolve */
697
698 /*  ===================================================================== */
699 /* Subroutine */ int zgelss_(integer *m, integer *n, integer *nrhs, 
700         doublecomplex *a, integer *lda, doublecomplex *b, integer *ldb, 
701         doublereal *s, doublereal *rcond, integer *rank, doublecomplex *work, 
702         integer *lwork, doublereal *rwork, integer *info)
703 {
704     /* System generated locals */
705     integer a_dim1, a_offset, b_dim1, b_offset, i__1, i__2, i__3;
706     doublereal d__1;
707
708     /* Local variables */
709     doublereal anrm, bnrm;
710     integer itau, lwork_zgebrd__, lwork_zgelqf__, i__, lwork_zgeqrf__, 
711             lwork_zungbr__, lwork_zunmbr__, iascl, ibscl, lwork_zunmlq__, 
712             chunk, lwork_zunmqr__;
713     doublereal sfmin;
714     integer minmn;
715     extern /* Subroutine */ int zgemm_(char *, char *, integer *, integer *, 
716             integer *, doublecomplex *, doublecomplex *, integer *, 
717             doublecomplex *, integer *, doublecomplex *, doublecomplex *, 
718             integer *);
719     integer maxmn, itaup, itauq, mnthr;
720     extern /* Subroutine */ int zgemv_(char *, integer *, integer *, 
721             doublecomplex *, doublecomplex *, integer *, doublecomplex *, 
722             integer *, doublecomplex *, doublecomplex *, integer *);
723     integer iwork;
724     extern /* Subroutine */ int zcopy_(integer *, doublecomplex *, integer *, 
725             doublecomplex *, integer *), dlabad_(doublereal *, doublereal *);
726     integer bl, ie, il;
727     extern doublereal dlamch_(char *);
728     integer mm;
729     extern /* Subroutine */ int dlascl_(char *, integer *, integer *, 
730             doublereal *, doublereal *, integer *, integer *, doublereal *, 
731             integer *, integer *), dlaset_(char *, integer *, integer 
732             *, doublereal *, doublereal *, doublereal *, integer *), 
733             xerbla_(char *, integer *, ftnlen), zgebrd_(integer *, integer *, 
734             doublecomplex *, integer *, doublereal *, doublereal *, 
735             doublecomplex *, doublecomplex *, doublecomplex *, integer *, 
736             integer *);
737     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
738             integer *, integer *, ftnlen, ftnlen);
739     extern doublereal zlange_(char *, integer *, integer *, doublecomplex *, 
740             integer *, doublereal *);
741     doublereal bignum;
742     extern /* Subroutine */ int zgelqf_(integer *, integer *, doublecomplex *,
743              integer *, doublecomplex *, doublecomplex *, integer *, integer *
744             ), zlascl_(char *, integer *, integer *, doublereal *, doublereal 
745             *, integer *, integer *, doublecomplex *, integer *, integer *), zgeqrf_(integer *, integer *, doublecomplex *, integer *,
746              doublecomplex *, doublecomplex *, integer *, integer *), zdrscl_(
747             integer *, doublereal *, doublecomplex *, integer *);
748     integer ldwork;
749     extern /* Subroutine */ int zlacpy_(char *, integer *, integer *, 
750             doublecomplex *, integer *, doublecomplex *, integer *), 
751             zlaset_(char *, integer *, integer *, doublecomplex *, 
752             doublecomplex *, doublecomplex *, integer *), zbdsqr_(
753             char *, integer *, integer *, integer *, integer *, doublereal *, 
754             doublereal *, doublecomplex *, integer *, doublecomplex *, 
755             integer *, doublecomplex *, integer *, doublereal *, integer *);
756     integer minwrk, maxwrk;
757     extern /* Subroutine */ int zungbr_(char *, integer *, integer *, integer 
758             *, doublecomplex *, integer *, doublecomplex *, doublecomplex *, 
759             integer *, integer *);
760     doublereal smlnum;
761     integer irwork;
762     extern /* Subroutine */ int zunmbr_(char *, char *, char *, integer *, 
763             integer *, integer *, doublecomplex *, integer *, doublecomplex *,
764              doublecomplex *, integer *, doublecomplex *, integer *, integer *
765             );
766     logical lquery;
767     extern /* Subroutine */ int zunmlq_(char *, char *, integer *, integer *, 
768             integer *, doublecomplex *, integer *, doublecomplex *, 
769             doublecomplex *, integer *, doublecomplex *, integer *, integer *), zunmqr_(char *, char *, integer *, integer *, 
770             integer *, doublecomplex *, integer *, doublecomplex *, 
771             doublecomplex *, integer *, doublecomplex *, integer *, integer *);
772     doublecomplex dum[1];
773     doublereal eps, thr;
774
775
776 /*  -- LAPACK driver routine (version 3.7.0) -- */
777 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
778 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
779 /*     June 2016 */
780
781
782 /*  ===================================================================== */
783
784
785 /*     Test the input arguments */
786
787     /* Parameter adjustments */
788     a_dim1 = *lda;
789     a_offset = 1 + a_dim1 * 1;
790     a -= a_offset;
791     b_dim1 = *ldb;
792     b_offset = 1 + b_dim1 * 1;
793     b -= b_offset;
794     --s;
795     --work;
796     --rwork;
797
798     /* Function Body */
799     *info = 0;
800     minmn = f2cmin(*m,*n);
801     maxmn = f2cmax(*m,*n);
802     lquery = *lwork == -1;
803     if (*m < 0) {
804         *info = -1;
805     } else if (*n < 0) {
806         *info = -2;
807     } else if (*nrhs < 0) {
808         *info = -3;
809     } else if (*lda < f2cmax(1,*m)) {
810         *info = -5;
811     } else if (*ldb < f2cmax(1,maxmn)) {
812         *info = -7;
813     }
814
815 /*     Compute workspace */
816 /*      (Note: Comments in the code beginning "Workspace:" describe the */
817 /*       minimal amount of workspace needed at that point in the code, */
818 /*       as well as the preferred amount for good performance. */
819 /*       CWorkspace refers to complex workspace, and RWorkspace refers */
820 /*       to real workspace. NB refers to the optimal block size for the */
821 /*       immediately following subroutine, as returned by ILAENV.) */
822
823     if (*info == 0) {
824         minwrk = 1;
825         maxwrk = 1;
826         if (minmn > 0) {
827             mm = *m;
828             mnthr = ilaenv_(&c__6, "ZGELSS", " ", m, n, nrhs, &c_n1, (ftnlen)
829                     6, (ftnlen)1);
830             if (*m >= *n && *m >= mnthr) {
831
832 /*              Path 1a - overdetermined, with many more rows than */
833 /*                        columns */
834
835 /*              Compute space needed for ZGEQRF */
836                 zgeqrf_(m, n, &a[a_offset], lda, dum, dum, &c_n1, info);
837                 lwork_zgeqrf__ = (integer) dum[0].r;
838 /*              Compute space needed for ZUNMQR */
839                 zunmqr_("L", "C", m, nrhs, n, &a[a_offset], lda, dum, &b[
840                         b_offset], ldb, dum, &c_n1, info);
841                 lwork_zunmqr__ = (integer) dum[0].r;
842                 mm = *n;
843 /* Computing MAX */
844                 i__1 = maxwrk, i__2 = *n + *n * ilaenv_(&c__1, "ZGEQRF", 
845                         " ", m, n, &c_n1, &c_n1, (ftnlen)6, (ftnlen)1);
846                 maxwrk = f2cmax(i__1,i__2);
847 /* Computing MAX */
848                 i__1 = maxwrk, i__2 = *n + *nrhs * ilaenv_(&c__1, "ZUNMQR", 
849                         "LC", m, nrhs, n, &c_n1, (ftnlen)6, (ftnlen)2);
850                 maxwrk = f2cmax(i__1,i__2);
851             }
852             if (*m >= *n) {
853
854 /*              Path 1 - overdetermined or exactly determined */
855
856 /*              Compute space needed for ZGEBRD */
857                 zgebrd_(&mm, n, &a[a_offset], lda, &s[1], &s[1], dum, dum, 
858                         dum, &c_n1, info);
859                 lwork_zgebrd__ = (integer) dum[0].r;
860 /*              Compute space needed for ZUNMBR */
861                 zunmbr_("Q", "L", "C", &mm, nrhs, n, &a[a_offset], lda, dum, &
862                         b[b_offset], ldb, dum, &c_n1, info);
863                 lwork_zunmbr__ = (integer) dum[0].r;
864 /*              Compute space needed for ZUNGBR */
865                 zungbr_("P", n, n, n, &a[a_offset], lda, dum, dum, &c_n1, 
866                         info);
867                 lwork_zungbr__ = (integer) dum[0].r;
868 /*              Compute total workspace needed */
869 /* Computing MAX */
870                 i__1 = maxwrk, i__2 = (*n << 1) + lwork_zgebrd__;
871                 maxwrk = f2cmax(i__1,i__2);
872 /* Computing MAX */
873                 i__1 = maxwrk, i__2 = (*n << 1) + lwork_zunmbr__;
874                 maxwrk = f2cmax(i__1,i__2);
875 /* Computing MAX */
876                 i__1 = maxwrk, i__2 = (*n << 1) + lwork_zungbr__;
877                 maxwrk = f2cmax(i__1,i__2);
878 /* Computing MAX */
879                 i__1 = maxwrk, i__2 = *n * *nrhs;
880                 maxwrk = f2cmax(i__1,i__2);
881                 minwrk = (*n << 1) + f2cmax(*nrhs,*m);
882             }
883             if (*n > *m) {
884                 minwrk = (*m << 1) + f2cmax(*nrhs,*n);
885                 if (*n >= mnthr) {
886
887 /*                 Path 2a - underdetermined, with many more columns */
888 /*                 than rows */
889
890 /*                 Compute space needed for ZGELQF */
891                     zgelqf_(m, n, &a[a_offset], lda, dum, dum, &c_n1, info);
892                     lwork_zgelqf__ = (integer) dum[0].r;
893 /*                 Compute space needed for ZGEBRD */
894                     zgebrd_(m, m, &a[a_offset], lda, &s[1], &s[1], dum, dum, 
895                             dum, &c_n1, info);
896                     lwork_zgebrd__ = (integer) dum[0].r;
897 /*                 Compute space needed for ZUNMBR */
898                     zunmbr_("Q", "L", "C", m, nrhs, n, &a[a_offset], lda, dum,
899                              &b[b_offset], ldb, dum, &c_n1, info);
900                     lwork_zunmbr__ = (integer) dum[0].r;
901 /*                 Compute space needed for ZUNGBR */
902                     zungbr_("P", m, m, m, &a[a_offset], lda, dum, dum, &c_n1, 
903                             info);
904                     lwork_zungbr__ = (integer) dum[0].r;
905 /*                 Compute space needed for ZUNMLQ */
906                     zunmlq_("L", "C", n, nrhs, m, &a[a_offset], lda, dum, &b[
907                             b_offset], ldb, dum, &c_n1, info);
908                     lwork_zunmlq__ = (integer) dum[0].r;
909 /*                 Compute total workspace needed */
910                     maxwrk = *m + lwork_zgelqf__;
911 /* Computing MAX */
912                     i__1 = maxwrk, i__2 = *m * 3 + *m * *m + lwork_zgebrd__;
913                     maxwrk = f2cmax(i__1,i__2);
914 /* Computing MAX */
915                     i__1 = maxwrk, i__2 = *m * 3 + *m * *m + lwork_zunmbr__;
916                     maxwrk = f2cmax(i__1,i__2);
917 /* Computing MAX */
918                     i__1 = maxwrk, i__2 = *m * 3 + *m * *m + lwork_zungbr__;
919                     maxwrk = f2cmax(i__1,i__2);
920                     if (*nrhs > 1) {
921 /* Computing MAX */
922                         i__1 = maxwrk, i__2 = *m * *m + *m + *m * *nrhs;
923                         maxwrk = f2cmax(i__1,i__2);
924                     } else {
925 /* Computing MAX */
926                         i__1 = maxwrk, i__2 = *m * *m + (*m << 1);
927                         maxwrk = f2cmax(i__1,i__2);
928                     }
929 /* Computing MAX */
930                     i__1 = maxwrk, i__2 = *m + lwork_zunmlq__;
931                     maxwrk = f2cmax(i__1,i__2);
932                 } else {
933
934 /*                 Path 2 - underdetermined */
935
936 /*                 Compute space needed for ZGEBRD */
937                     zgebrd_(m, n, &a[a_offset], lda, &s[1], &s[1], dum, dum, 
938                             dum, &c_n1, info);
939                     lwork_zgebrd__ = (integer) dum[0].r;
940 /*                 Compute space needed for ZUNMBR */
941                     zunmbr_("Q", "L", "C", m, nrhs, m, &a[a_offset], lda, dum,
942                              &b[b_offset], ldb, dum, &c_n1, info);
943                     lwork_zunmbr__ = (integer) dum[0].r;
944 /*                 Compute space needed for ZUNGBR */
945                     zungbr_("P", m, n, m, &a[a_offset], lda, dum, dum, &c_n1, 
946                             info);
947                     lwork_zungbr__ = (integer) dum[0].r;
948                     maxwrk = (*m << 1) + lwork_zgebrd__;
949 /* Computing MAX */
950                     i__1 = maxwrk, i__2 = (*m << 1) + lwork_zunmbr__;
951                     maxwrk = f2cmax(i__1,i__2);
952 /* Computing MAX */
953                     i__1 = maxwrk, i__2 = (*m << 1) + lwork_zungbr__;
954                     maxwrk = f2cmax(i__1,i__2);
955 /* Computing MAX */
956                     i__1 = maxwrk, i__2 = *n * *nrhs;
957                     maxwrk = f2cmax(i__1,i__2);
958                 }
959             }
960             maxwrk = f2cmax(minwrk,maxwrk);
961         }
962         work[1].r = (doublereal) maxwrk, work[1].i = 0.;
963
964         if (*lwork < minwrk && ! lquery) {
965             *info = -12;
966         }
967     }
968
969     if (*info != 0) {
970         i__1 = -(*info);
971         xerbla_("ZGELSS", &i__1, (ftnlen)6);
972         return 0;
973     } else if (lquery) {
974         return 0;
975     }
976
977 /*     Quick return if possible */
978
979     if (*m == 0 || *n == 0) {
980         *rank = 0;
981         return 0;
982     }
983
984 /*     Get machine parameters */
985
986     eps = dlamch_("P");
987     sfmin = dlamch_("S");
988     smlnum = sfmin / eps;
989     bignum = 1. / smlnum;
990     dlabad_(&smlnum, &bignum);
991
992 /*     Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
993
994     anrm = zlange_("M", m, n, &a[a_offset], lda, &rwork[1]);
995     iascl = 0;
996     if (anrm > 0. && anrm < smlnum) {
997
998 /*        Scale matrix norm up to SMLNUM */
999
1000         zlascl_("G", &c__0, &c__0, &anrm, &smlnum, m, n, &a[a_offset], lda, 
1001                 info);
1002         iascl = 1;
1003     } else if (anrm > bignum) {
1004
1005 /*        Scale matrix norm down to BIGNUM */
1006
1007         zlascl_("G", &c__0, &c__0, &anrm, &bignum, m, n, &a[a_offset], lda, 
1008                 info);
1009         iascl = 2;
1010     } else if (anrm == 0.) {
1011
1012 /*        Matrix all zero. Return zero solution. */
1013
1014         i__1 = f2cmax(*m,*n);
1015         zlaset_("F", &i__1, nrhs, &c_b1, &c_b1, &b[b_offset], ldb);
1016         dlaset_("F", &minmn, &c__1, &c_b59, &c_b59, &s[1], &minmn);
1017         *rank = 0;
1018         goto L70;
1019     }
1020
1021 /*     Scale B if f2cmax element outside range [SMLNUM,BIGNUM] */
1022
1023     bnrm = zlange_("M", m, nrhs, &b[b_offset], ldb, &rwork[1]);
1024     ibscl = 0;
1025     if (bnrm > 0. && bnrm < smlnum) {
1026
1027 /*        Scale matrix norm up to SMLNUM */
1028
1029         zlascl_("G", &c__0, &c__0, &bnrm, &smlnum, m, nrhs, &b[b_offset], ldb,
1030                  info);
1031         ibscl = 1;
1032     } else if (bnrm > bignum) {
1033
1034 /*        Scale matrix norm down to BIGNUM */
1035
1036         zlascl_("G", &c__0, &c__0, &bnrm, &bignum, m, nrhs, &b[b_offset], ldb,
1037                  info);
1038         ibscl = 2;
1039     }
1040
1041 /*     Overdetermined case */
1042
1043     if (*m >= *n) {
1044
1045 /*        Path 1 - overdetermined or exactly determined */
1046
1047         mm = *m;
1048         if (*m >= mnthr) {
1049
1050 /*           Path 1a - overdetermined, with many more rows than columns */
1051
1052             mm = *n;
1053             itau = 1;
1054             iwork = itau + *n;
1055
1056 /*           Compute A=Q*R */
1057 /*           (CWorkspace: need 2*N, prefer N+N*NB) */
1058 /*           (RWorkspace: none) */
1059
1060             i__1 = *lwork - iwork + 1;
1061             zgeqrf_(m, n, &a[a_offset], lda, &work[itau], &work[iwork], &i__1,
1062                      info);
1063
1064 /*           Multiply B by transpose(Q) */
1065 /*           (CWorkspace: need N+NRHS, prefer N+NRHS*NB) */
1066 /*           (RWorkspace: none) */
1067
1068             i__1 = *lwork - iwork + 1;
1069             zunmqr_("L", "C", m, nrhs, n, &a[a_offset], lda, &work[itau], &b[
1070                     b_offset], ldb, &work[iwork], &i__1, info);
1071
1072 /*           Zero out below R */
1073
1074             if (*n > 1) {
1075                 i__1 = *n - 1;
1076                 i__2 = *n - 1;
1077                 zlaset_("L", &i__1, &i__2, &c_b1, &c_b1, &a[a_dim1 + 2], lda);
1078             }
1079         }
1080
1081         ie = 1;
1082         itauq = 1;
1083         itaup = itauq + *n;
1084         iwork = itaup + *n;
1085
1086 /*        Bidiagonalize R in A */
1087 /*        (CWorkspace: need 2*N+MM, prefer 2*N+(MM+N)*NB) */
1088 /*        (RWorkspace: need N) */
1089
1090         i__1 = *lwork - iwork + 1;
1091         zgebrd_(&mm, n, &a[a_offset], lda, &s[1], &rwork[ie], &work[itauq], &
1092                 work[itaup], &work[iwork], &i__1, info);
1093
1094 /*        Multiply B by transpose of left bidiagonalizing vectors of R */
1095 /*        (CWorkspace: need 2*N+NRHS, prefer 2*N+NRHS*NB) */
1096 /*        (RWorkspace: none) */
1097
1098         i__1 = *lwork - iwork + 1;
1099         zunmbr_("Q", "L", "C", &mm, nrhs, n, &a[a_offset], lda, &work[itauq], 
1100                 &b[b_offset], ldb, &work[iwork], &i__1, info);
1101
1102 /*        Generate right bidiagonalizing vectors of R in A */
1103 /*        (CWorkspace: need 3*N-1, prefer 2*N+(N-1)*NB) */
1104 /*        (RWorkspace: none) */
1105
1106         i__1 = *lwork - iwork + 1;
1107         zungbr_("P", n, n, n, &a[a_offset], lda, &work[itaup], &work[iwork], &
1108                 i__1, info);
1109         irwork = ie + *n;
1110
1111 /*        Perform bidiagonal QR iteration */
1112 /*          multiply B by transpose of left singular vectors */
1113 /*          compute right singular vectors in A */
1114 /*        (CWorkspace: none) */
1115 /*        (RWorkspace: need BDSPAC) */
1116
1117         zbdsqr_("U", n, n, &c__0, nrhs, &s[1], &rwork[ie], &a[a_offset], lda, 
1118                 dum, &c__1, &b[b_offset], ldb, &rwork[irwork], info);
1119         if (*info != 0) {
1120             goto L70;
1121         }
1122
1123 /*        Multiply B by reciprocals of singular values */
1124
1125 /* Computing MAX */
1126         d__1 = *rcond * s[1];
1127         thr = f2cmax(d__1,sfmin);
1128         if (*rcond < 0.) {
1129 /* Computing MAX */
1130             d__1 = eps * s[1];
1131             thr = f2cmax(d__1,sfmin);
1132         }
1133         *rank = 0;
1134         i__1 = *n;
1135         for (i__ = 1; i__ <= i__1; ++i__) {
1136             if (s[i__] > thr) {
1137                 zdrscl_(nrhs, &s[i__], &b[i__ + b_dim1], ldb);
1138                 ++(*rank);
1139             } else {
1140                 zlaset_("F", &c__1, nrhs, &c_b1, &c_b1, &b[i__ + b_dim1], ldb);
1141             }
1142 /* L10: */
1143         }
1144
1145 /*        Multiply B by right singular vectors */
1146 /*        (CWorkspace: need N, prefer N*NRHS) */
1147 /*        (RWorkspace: none) */
1148
1149         if (*lwork >= *ldb * *nrhs && *nrhs > 1) {
1150             zgemm_("C", "N", n, nrhs, n, &c_b2, &a[a_offset], lda, &b[
1151                     b_offset], ldb, &c_b1, &work[1], ldb);
1152             zlacpy_("G", n, nrhs, &work[1], ldb, &b[b_offset], ldb)
1153                     ;
1154         } else if (*nrhs > 1) {
1155             chunk = *lwork / *n;
1156             i__1 = *nrhs;
1157             i__2 = chunk;
1158             for (i__ = 1; i__2 < 0 ? i__ >= i__1 : i__ <= i__1; i__ += i__2) {
1159 /* Computing MIN */
1160                 i__3 = *nrhs - i__ + 1;
1161                 bl = f2cmin(i__3,chunk);
1162                 zgemm_("C", "N", n, &bl, n, &c_b2, &a[a_offset], lda, &b[i__ *
1163                          b_dim1 + 1], ldb, &c_b1, &work[1], n);
1164                 zlacpy_("G", n, &bl, &work[1], n, &b[i__ * b_dim1 + 1], ldb);
1165 /* L20: */
1166             }
1167         } else {
1168             zgemv_("C", n, n, &c_b2, &a[a_offset], lda, &b[b_offset], &c__1, &
1169                     c_b1, &work[1], &c__1);
1170             zcopy_(n, &work[1], &c__1, &b[b_offset], &c__1);
1171         }
1172
1173     } else /* if(complicated condition) */ {
1174 /* Computing MAX */
1175         i__2 = f2cmax(*m,*nrhs), i__1 = *n - (*m << 1);
1176         if (*n >= mnthr && *lwork >= *m * 3 + *m * *m + f2cmax(i__2,i__1)) {
1177
1178 /*        Underdetermined case, M much less than N */
1179
1180 /*        Path 2a - underdetermined, with many more columns than rows */
1181 /*        and sufficient workspace for an efficient algorithm */
1182
1183             ldwork = *m;
1184 /* Computing MAX */
1185             i__2 = f2cmax(*m,*nrhs), i__1 = *n - (*m << 1);
1186             if (*lwork >= *m * 3 + *m * *lda + f2cmax(i__2,i__1)) {
1187                 ldwork = *lda;
1188             }
1189             itau = 1;
1190             iwork = *m + 1;
1191
1192 /*        Compute A=L*Q */
1193 /*        (CWorkspace: need 2*M, prefer M+M*NB) */
1194 /*        (RWorkspace: none) */
1195
1196             i__2 = *lwork - iwork + 1;
1197             zgelqf_(m, n, &a[a_offset], lda, &work[itau], &work[iwork], &i__2,
1198                      info);
1199             il = iwork;
1200
1201 /*        Copy L to WORK(IL), zeroing out above it */
1202
1203             zlacpy_("L", m, m, &a[a_offset], lda, &work[il], &ldwork);
1204             i__2 = *m - 1;
1205             i__1 = *m - 1;
1206             zlaset_("U", &i__2, &i__1, &c_b1, &c_b1, &work[il + ldwork], &
1207                     ldwork);
1208             ie = 1;
1209             itauq = il + ldwork * *m;
1210             itaup = itauq + *m;
1211             iwork = itaup + *m;
1212
1213 /*        Bidiagonalize L in WORK(IL) */
1214 /*        (CWorkspace: need M*M+4*M, prefer M*M+3*M+2*M*NB) */
1215 /*        (RWorkspace: need M) */
1216
1217             i__2 = *lwork - iwork + 1;
1218             zgebrd_(m, m, &work[il], &ldwork, &s[1], &rwork[ie], &work[itauq],
1219                      &work[itaup], &work[iwork], &i__2, info);
1220
1221 /*        Multiply B by transpose of left bidiagonalizing vectors of L */
1222 /*        (CWorkspace: need M*M+3*M+NRHS, prefer M*M+3*M+NRHS*NB) */
1223 /*        (RWorkspace: none) */
1224
1225             i__2 = *lwork - iwork + 1;
1226             zunmbr_("Q", "L", "C", m, nrhs, m, &work[il], &ldwork, &work[
1227                     itauq], &b[b_offset], ldb, &work[iwork], &i__2, info);
1228
1229 /*        Generate right bidiagonalizing vectors of R in WORK(IL) */
1230 /*        (CWorkspace: need M*M+4*M-1, prefer M*M+3*M+(M-1)*NB) */
1231 /*        (RWorkspace: none) */
1232
1233             i__2 = *lwork - iwork + 1;
1234             zungbr_("P", m, m, m, &work[il], &ldwork, &work[itaup], &work[
1235                     iwork], &i__2, info);
1236             irwork = ie + *m;
1237
1238 /*        Perform bidiagonal QR iteration, computing right singular */
1239 /*        vectors of L in WORK(IL) and multiplying B by transpose of */
1240 /*        left singular vectors */
1241 /*        (CWorkspace: need M*M) */
1242 /*        (RWorkspace: need BDSPAC) */
1243
1244             zbdsqr_("U", m, m, &c__0, nrhs, &s[1], &rwork[ie], &work[il], &
1245                     ldwork, &a[a_offset], lda, &b[b_offset], ldb, &rwork[
1246                     irwork], info);
1247             if (*info != 0) {
1248                 goto L70;
1249             }
1250
1251 /*        Multiply B by reciprocals of singular values */
1252
1253 /* Computing MAX */
1254             d__1 = *rcond * s[1];
1255             thr = f2cmax(d__1,sfmin);
1256             if (*rcond < 0.) {
1257 /* Computing MAX */
1258                 d__1 = eps * s[1];
1259                 thr = f2cmax(d__1,sfmin);
1260             }
1261             *rank = 0;
1262             i__2 = *m;
1263             for (i__ = 1; i__ <= i__2; ++i__) {
1264                 if (s[i__] > thr) {
1265                     zdrscl_(nrhs, &s[i__], &b[i__ + b_dim1], ldb);
1266                     ++(*rank);
1267                 } else {
1268                     zlaset_("F", &c__1, nrhs, &c_b1, &c_b1, &b[i__ + b_dim1], 
1269                             ldb);
1270                 }
1271 /* L30: */
1272             }
1273             iwork = il + *m * ldwork;
1274
1275 /*        Multiply B by right singular vectors of L in WORK(IL) */
1276 /*        (CWorkspace: need M*M+2*M, prefer M*M+M+M*NRHS) */
1277 /*        (RWorkspace: none) */
1278
1279             if (*lwork >= *ldb * *nrhs + iwork - 1 && *nrhs > 1) {
1280                 zgemm_("C", "N", m, nrhs, m, &c_b2, &work[il], &ldwork, &b[
1281                         b_offset], ldb, &c_b1, &work[iwork], ldb);
1282                 zlacpy_("G", m, nrhs, &work[iwork], ldb, &b[b_offset], ldb);
1283             } else if (*nrhs > 1) {
1284                 chunk = (*lwork - iwork + 1) / *m;
1285                 i__2 = *nrhs;
1286                 i__1 = chunk;
1287                 for (i__ = 1; i__1 < 0 ? i__ >= i__2 : i__ <= i__2; i__ += 
1288                         i__1) {
1289 /* Computing MIN */
1290                     i__3 = *nrhs - i__ + 1;
1291                     bl = f2cmin(i__3,chunk);
1292                     zgemm_("C", "N", m, &bl, m, &c_b2, &work[il], &ldwork, &b[
1293                             i__ * b_dim1 + 1], ldb, &c_b1, &work[iwork], m);
1294                     zlacpy_("G", m, &bl, &work[iwork], m, &b[i__ * b_dim1 + 1]
1295                             , ldb);
1296 /* L40: */
1297                 }
1298             } else {
1299                 zgemv_("C", m, m, &c_b2, &work[il], &ldwork, &b[b_dim1 + 1], &
1300                         c__1, &c_b1, &work[iwork], &c__1);
1301                 zcopy_(m, &work[iwork], &c__1, &b[b_dim1 + 1], &c__1);
1302             }
1303
1304 /*        Zero out below first M rows of B */
1305
1306             i__1 = *n - *m;
1307             zlaset_("F", &i__1, nrhs, &c_b1, &c_b1, &b[*m + 1 + b_dim1], ldb);
1308             iwork = itau + *m;
1309
1310 /*        Multiply transpose(Q) by B */
1311 /*        (CWorkspace: need M+NRHS, prefer M+NHRS*NB) */
1312 /*        (RWorkspace: none) */
1313
1314             i__1 = *lwork - iwork + 1;
1315             zunmlq_("L", "C", n, nrhs, m, &a[a_offset], lda, &work[itau], &b[
1316                     b_offset], ldb, &work[iwork], &i__1, info);
1317
1318         } else {
1319
1320 /*        Path 2 - remaining underdetermined cases */
1321
1322             ie = 1;
1323             itauq = 1;
1324             itaup = itauq + *m;
1325             iwork = itaup + *m;
1326
1327 /*        Bidiagonalize A */
1328 /*        (CWorkspace: need 3*M, prefer 2*M+(M+N)*NB) */
1329 /*        (RWorkspace: need N) */
1330
1331             i__1 = *lwork - iwork + 1;
1332             zgebrd_(m, n, &a[a_offset], lda, &s[1], &rwork[ie], &work[itauq], 
1333                     &work[itaup], &work[iwork], &i__1, info);
1334
1335 /*        Multiply B by transpose of left bidiagonalizing vectors */
1336 /*        (CWorkspace: need 2*M+NRHS, prefer 2*M+NRHS*NB) */
1337 /*        (RWorkspace: none) */
1338
1339             i__1 = *lwork - iwork + 1;
1340             zunmbr_("Q", "L", "C", m, nrhs, n, &a[a_offset], lda, &work[itauq]
1341                     , &b[b_offset], ldb, &work[iwork], &i__1, info);
1342
1343 /*        Generate right bidiagonalizing vectors in A */
1344 /*        (CWorkspace: need 3*M, prefer 2*M+M*NB) */
1345 /*        (RWorkspace: none) */
1346
1347             i__1 = *lwork - iwork + 1;
1348             zungbr_("P", m, n, m, &a[a_offset], lda, &work[itaup], &work[
1349                     iwork], &i__1, info);
1350             irwork = ie + *m;
1351
1352 /*        Perform bidiagonal QR iteration, */
1353 /*           computing right singular vectors of A in A and */
1354 /*           multiplying B by transpose of left singular vectors */
1355 /*        (CWorkspace: none) */
1356 /*        (RWorkspace: need BDSPAC) */
1357
1358             zbdsqr_("L", m, n, &c__0, nrhs, &s[1], &rwork[ie], &a[a_offset], 
1359                     lda, dum, &c__1, &b[b_offset], ldb, &rwork[irwork], info);
1360             if (*info != 0) {
1361                 goto L70;
1362             }
1363
1364 /*        Multiply B by reciprocals of singular values */
1365
1366 /* Computing MAX */
1367             d__1 = *rcond * s[1];
1368             thr = f2cmax(d__1,sfmin);
1369             if (*rcond < 0.) {
1370 /* Computing MAX */
1371                 d__1 = eps * s[1];
1372                 thr = f2cmax(d__1,sfmin);
1373             }
1374             *rank = 0;
1375             i__1 = *m;
1376             for (i__ = 1; i__ <= i__1; ++i__) {
1377                 if (s[i__] > thr) {
1378                     zdrscl_(nrhs, &s[i__], &b[i__ + b_dim1], ldb);
1379                     ++(*rank);
1380                 } else {
1381                     zlaset_("F", &c__1, nrhs, &c_b1, &c_b1, &b[i__ + b_dim1], 
1382                             ldb);
1383                 }
1384 /* L50: */
1385             }
1386
1387 /*        Multiply B by right singular vectors of A */
1388 /*        (CWorkspace: need N, prefer N*NRHS) */
1389 /*        (RWorkspace: none) */
1390
1391             if (*lwork >= *ldb * *nrhs && *nrhs > 1) {
1392                 zgemm_("C", "N", n, nrhs, m, &c_b2, &a[a_offset], lda, &b[
1393                         b_offset], ldb, &c_b1, &work[1], ldb);
1394                 zlacpy_("G", n, nrhs, &work[1], ldb, &b[b_offset], ldb);
1395             } else if (*nrhs > 1) {
1396                 chunk = *lwork / *n;
1397                 i__1 = *nrhs;
1398                 i__2 = chunk;
1399                 for (i__ = 1; i__2 < 0 ? i__ >= i__1 : i__ <= i__1; i__ += 
1400                         i__2) {
1401 /* Computing MIN */
1402                     i__3 = *nrhs - i__ + 1;
1403                     bl = f2cmin(i__3,chunk);
1404                     zgemm_("C", "N", n, &bl, m, &c_b2, &a[a_offset], lda, &b[
1405                             i__ * b_dim1 + 1], ldb, &c_b1, &work[1], n);
1406                     zlacpy_("F", n, &bl, &work[1], n, &b[i__ * b_dim1 + 1], 
1407                             ldb);
1408 /* L60: */
1409                 }
1410             } else {
1411                 zgemv_("C", m, n, &c_b2, &a[a_offset], lda, &b[b_offset], &
1412                         c__1, &c_b1, &work[1], &c__1);
1413                 zcopy_(n, &work[1], &c__1, &b[b_offset], &c__1);
1414             }
1415         }
1416     }
1417
1418 /*     Undo scaling */
1419
1420     if (iascl == 1) {
1421         zlascl_("G", &c__0, &c__0, &anrm, &smlnum, n, nrhs, &b[b_offset], ldb,
1422                  info);
1423         dlascl_("G", &c__0, &c__0, &smlnum, &anrm, &minmn, &c__1, &s[1], &
1424                 minmn, info);
1425     } else if (iascl == 2) {
1426         zlascl_("G", &c__0, &c__0, &anrm, &bignum, n, nrhs, &b[b_offset], ldb,
1427                  info);
1428         dlascl_("G", &c__0, &c__0, &bignum, &anrm, &minmn, &c__1, &s[1], &
1429                 minmn, info);
1430     }
1431     if (ibscl == 1) {
1432         zlascl_("G", &c__0, &c__0, &smlnum, &bnrm, n, nrhs, &b[b_offset], ldb,
1433                  info);
1434     } else if (ibscl == 2) {
1435         zlascl_("G", &c__0, &c__0, &bignum, &bnrm, n, nrhs, &b[b_offset], ldb,
1436                  info);
1437     }
1438 L70:
1439     work[1].r = (doublereal) maxwrk, work[1].i = 0.;
1440     return 0;
1441
1442 /*     End of ZGELSS */
1443
1444 } /* zgelss_ */
1445