C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dlalsd.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static integer c__1 = 1;
516 static doublereal c_b6 = 0.;
517 static integer c__0 = 0;
518 static doublereal c_b11 = 1.;
519
520 /* > \brief \b DLALSD uses the singular value decomposition of A to solve the least squares problem. */
521
522 /*  =========== DOCUMENTATION =========== */
523
524 /* Online html documentation available at */
525 /*            http://www.netlib.org/lapack/explore-html/ */
526
527 /* > \htmlonly */
528 /* > Download DLALSD + dependencies */
529 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlalsd.
530 f"> */
531 /* > [TGZ]</a> */
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlalsd.
533 f"> */
534 /* > [ZIP]</a> */
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlalsd.
536 f"> */
537 /* > [TXT]</a> */
538 /* > \endhtmlonly */
539
540 /*  Definition: */
541 /*  =========== */
542
543 /*       SUBROUTINE DLALSD( UPLO, SMLSIZ, N, NRHS, D, E, B, LDB, RCOND, */
544 /*                          RANK, WORK, IWORK, INFO ) */
545
546 /*       CHARACTER          UPLO */
547 /*       INTEGER            INFO, LDB, N, NRHS, RANK, SMLSIZ */
548 /*       DOUBLE PRECISION   RCOND */
549 /*       INTEGER            IWORK( * ) */
550 /*       DOUBLE PRECISION   B( LDB, * ), D( * ), E( * ), WORK( * ) */
551
552
553 /* > \par Purpose: */
554 /*  ============= */
555 /* > */
556 /* > \verbatim */
557 /* > */
558 /* > DLALSD uses the singular value decomposition of A to solve the least */
559 /* > squares problem of finding X to minimize the Euclidean norm of each */
560 /* > column of A*X-B, where A is N-by-N upper bidiagonal, and X and B */
561 /* > are N-by-NRHS. The solution X overwrites B. */
562 /* > */
563 /* > The singular values of A smaller than RCOND times the largest */
564 /* > singular value are treated as zero in solving the least squares */
565 /* > problem; in this case a minimum norm solution is returned. */
566 /* > The actual singular values are returned in D in ascending order. */
567 /* > */
568 /* > This code makes very mild assumptions about floating point */
569 /* > arithmetic. It will work on machines with a guard digit in */
570 /* > add/subtract, or on those binary machines without guard digits */
571 /* > which subtract like the Cray XMP, Cray YMP, Cray C 90, or Cray 2. */
572 /* > It could conceivably fail on hexadecimal or decimal machines */
573 /* > without guard digits, but we know of none. */
574 /* > \endverbatim */
575
576 /*  Arguments: */
577 /*  ========== */
578
579 /* > \param[in] UPLO */
580 /* > \verbatim */
581 /* >          UPLO is CHARACTER*1 */
582 /* >         = 'U': D and E define an upper bidiagonal matrix. */
583 /* >         = 'L': D and E define a  lower bidiagonal matrix. */
584 /* > \endverbatim */
585 /* > */
586 /* > \param[in] SMLSIZ */
587 /* > \verbatim */
588 /* >          SMLSIZ is INTEGER */
589 /* >         The maximum size of the subproblems at the bottom of the */
590 /* >         computation tree. */
591 /* > \endverbatim */
592 /* > */
593 /* > \param[in] N */
594 /* > \verbatim */
595 /* >          N is INTEGER */
596 /* >         The dimension of the  bidiagonal matrix.  N >= 0. */
597 /* > \endverbatim */
598 /* > */
599 /* > \param[in] NRHS */
600 /* > \verbatim */
601 /* >          NRHS is INTEGER */
602 /* >         The number of columns of B. NRHS must be at least 1. */
603 /* > \endverbatim */
604 /* > */
605 /* > \param[in,out] D */
606 /* > \verbatim */
607 /* >          D is DOUBLE PRECISION array, dimension (N) */
608 /* >         On entry D contains the main diagonal of the bidiagonal */
609 /* >         matrix. On exit, if INFO = 0, D contains its singular values. */
610 /* > \endverbatim */
611 /* > */
612 /* > \param[in,out] E */
613 /* > \verbatim */
614 /* >          E is DOUBLE PRECISION array, dimension (N-1) */
615 /* >         Contains the super-diagonal entries of the bidiagonal matrix. */
616 /* >         On exit, E has been destroyed. */
617 /* > \endverbatim */
618 /* > */
619 /* > \param[in,out] B */
620 /* > \verbatim */
621 /* >          B is DOUBLE PRECISION array, dimension (LDB,NRHS) */
622 /* >         On input, B contains the right hand sides of the least */
623 /* >         squares problem. On output, B contains the solution X. */
624 /* > \endverbatim */
625 /* > */
626 /* > \param[in] LDB */
627 /* > \verbatim */
628 /* >          LDB is INTEGER */
629 /* >         The leading dimension of B in the calling subprogram. */
630 /* >         LDB must be at least f2cmax(1,N). */
631 /* > \endverbatim */
632 /* > */
633 /* > \param[in] RCOND */
634 /* > \verbatim */
635 /* >          RCOND is DOUBLE PRECISION */
636 /* >         The singular values of A less than or equal to RCOND times */
637 /* >         the largest singular value are treated as zero in solving */
638 /* >         the least squares problem. If RCOND is negative, */
639 /* >         machine precision is used instead. */
640 /* >         For example, if diag(S)*X=B were the least squares problem, */
641 /* >         where diag(S) is a diagonal matrix of singular values, the */
642 /* >         solution would be X(i) = B(i) / S(i) if S(i) is greater than */
643 /* >         RCOND*f2cmax(S), and X(i) = 0 if S(i) is less than or equal to */
644 /* >         RCOND*f2cmax(S). */
645 /* > \endverbatim */
646 /* > */
647 /* > \param[out] RANK */
648 /* > \verbatim */
649 /* >          RANK is INTEGER */
650 /* >         The number of singular values of A greater than RCOND times */
651 /* >         the largest singular value. */
652 /* > \endverbatim */
653 /* > */
654 /* > \param[out] WORK */
655 /* > \verbatim */
656 /* >          WORK is DOUBLE PRECISION array, dimension at least */
657 /* >         (9*N + 2*N*SMLSIZ + 8*N*NLVL + N*NRHS + (SMLSIZ+1)**2), */
658 /* >         where NLVL = f2cmax(0, INT(log_2 (N/(SMLSIZ+1))) + 1). */
659 /* > \endverbatim */
660 /* > */
661 /* > \param[out] IWORK */
662 /* > \verbatim */
663 /* >          IWORK is INTEGER array, dimension at least */
664 /* >         (3*N*NLVL + 11*N) */
665 /* > \endverbatim */
666 /* > */
667 /* > \param[out] INFO */
668 /* > \verbatim */
669 /* >          INFO is INTEGER */
670 /* >         = 0:  successful exit. */
671 /* >         < 0:  if INFO = -i, the i-th argument had an illegal value. */
672 /* >         > 0:  The algorithm failed to compute a singular value while */
673 /* >               working on the submatrix lying in rows and columns */
674 /* >               INFO/(N+1) through MOD(INFO,N+1). */
675 /* > \endverbatim */
676
677 /*  Authors: */
678 /*  ======== */
679
680 /* > \author Univ. of Tennessee */
681 /* > \author Univ. of California Berkeley */
682 /* > \author Univ. of Colorado Denver */
683 /* > \author NAG Ltd. */
684
685 /* > \date December 2016 */
686
687 /* > \ingroup doubleOTHERcomputational */
688
689 /* > \par Contributors: */
690 /*  ================== */
691 /* > */
692 /* >     Ming Gu and Ren-Cang Li, Computer Science Division, University of */
693 /* >       California at Berkeley, USA \n */
694 /* >     Osni Marques, LBNL/NERSC, USA \n */
695
696 /*  ===================================================================== */
697 /* Subroutine */ int dlalsd_(char *uplo, integer *smlsiz, integer *n, integer 
698         *nrhs, doublereal *d__, doublereal *e, doublereal *b, integer *ldb, 
699         doublereal *rcond, integer *rank, doublereal *work, integer *iwork, 
700         integer *info)
701 {
702     /* System generated locals */
703     integer b_dim1, b_offset, i__1, i__2;
704     doublereal d__1;
705
706     /* Local variables */
707     integer difl, difr;
708     doublereal rcnd;
709     integer perm, nsub;
710     extern /* Subroutine */ int drot_(integer *, doublereal *, integer *, 
711             doublereal *, integer *, doublereal *, doublereal *);
712     integer nlvl, sqre, bxst, c__, i__, j, k;
713     doublereal r__;
714     integer s, u;
715     extern /* Subroutine */ int dgemm_(char *, char *, integer *, integer *, 
716             integer *, doublereal *, doublereal *, integer *, doublereal *, 
717             integer *, doublereal *, doublereal *, integer *);
718     integer z__;
719     extern /* Subroutine */ int dcopy_(integer *, doublereal *, integer *, 
720             doublereal *, integer *);
721     integer poles, sizei, nsize, nwork, icmpq1, icmpq2;
722     doublereal cs;
723     extern doublereal dlamch_(char *);
724     extern /* Subroutine */ int dlasda_(integer *, integer *, integer *, 
725             integer *, doublereal *, doublereal *, doublereal *, integer *, 
726             doublereal *, integer *, doublereal *, doublereal *, doublereal *,
727              doublereal *, integer *, integer *, integer *, integer *, 
728             doublereal *, doublereal *, doublereal *, doublereal *, integer *,
729              integer *);
730     integer bx;
731     extern /* Subroutine */ int dlalsa_(integer *, integer *, integer *, 
732             integer *, doublereal *, integer *, doublereal *, integer *, 
733             doublereal *, integer *, doublereal *, integer *, doublereal *, 
734             doublereal *, doublereal *, doublereal *, integer *, integer *, 
735             integer *, integer *, doublereal *, doublereal *, doublereal *, 
736             doublereal *, integer *, integer *);
737     doublereal sn;
738     extern /* Subroutine */ int dlascl_(char *, integer *, integer *, 
739             doublereal *, doublereal *, integer *, integer *, doublereal *, 
740             integer *, integer *);
741     extern integer idamax_(integer *, doublereal *, integer *);
742     integer st;
743     extern /* Subroutine */ int dlasdq_(char *, integer *, integer *, integer 
744             *, integer *, integer *, doublereal *, doublereal *, doublereal *,
745              integer *, doublereal *, integer *, doublereal *, integer *, 
746             doublereal *, integer *);
747     integer vt;
748     extern /* Subroutine */ int dlacpy_(char *, integer *, integer *, 
749             doublereal *, integer *, doublereal *, integer *), 
750             dlartg_(doublereal *, doublereal *, doublereal *, doublereal *, 
751             doublereal *), dlaset_(char *, integer *, integer *, doublereal *,
752              doublereal *, doublereal *, integer *), xerbla_(char *, 
753             integer *, ftnlen);
754     integer givcol;
755     extern doublereal dlanst_(char *, integer *, doublereal *, doublereal *);
756     extern /* Subroutine */ int dlasrt_(char *, integer *, doublereal *, 
757             integer *);
758     doublereal orgnrm;
759     integer givnum, givptr, nm1, smlszp, st1;
760     doublereal eps;
761     integer iwk;
762     doublereal tol;
763
764
765 /*  -- LAPACK computational routine (version 3.7.0) -- */
766 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
767 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
768 /*     December 2016 */
769
770
771 /*  ===================================================================== */
772
773
774 /*     Test the input parameters. */
775
776     /* Parameter adjustments */
777     --d__;
778     --e;
779     b_dim1 = *ldb;
780     b_offset = 1 + b_dim1 * 1;
781     b -= b_offset;
782     --work;
783     --iwork;
784
785     /* Function Body */
786     *info = 0;
787
788     if (*n < 0) {
789         *info = -3;
790     } else if (*nrhs < 1) {
791         *info = -4;
792     } else if (*ldb < 1 || *ldb < *n) {
793         *info = -8;
794     }
795     if (*info != 0) {
796         i__1 = -(*info);
797         xerbla_("DLALSD", &i__1, (ftnlen)6);
798         return 0;
799     }
800
801     eps = dlamch_("Epsilon");
802
803 /*     Set up the tolerance. */
804
805     if (*rcond <= 0. || *rcond >= 1.) {
806         rcnd = eps;
807     } else {
808         rcnd = *rcond;
809     }
810
811     *rank = 0;
812
813 /*     Quick return if possible. */
814
815     if (*n == 0) {
816         return 0;
817     } else if (*n == 1) {
818         if (d__[1] == 0.) {
819             dlaset_("A", &c__1, nrhs, &c_b6, &c_b6, &b[b_offset], ldb);
820         } else {
821             *rank = 1;
822             dlascl_("G", &c__0, &c__0, &d__[1], &c_b11, &c__1, nrhs, &b[
823                     b_offset], ldb, info);
824             d__[1] = abs(d__[1]);
825         }
826         return 0;
827     }
828
829 /*     Rotate the matrix if it is lower bidiagonal. */
830
831     if (*(unsigned char *)uplo == 'L') {
832         i__1 = *n - 1;
833         for (i__ = 1; i__ <= i__1; ++i__) {
834             dlartg_(&d__[i__], &e[i__], &cs, &sn, &r__);
835             d__[i__] = r__;
836             e[i__] = sn * d__[i__ + 1];
837             d__[i__ + 1] = cs * d__[i__ + 1];
838             if (*nrhs == 1) {
839                 drot_(&c__1, &b[i__ + b_dim1], &c__1, &b[i__ + 1 + b_dim1], &
840                         c__1, &cs, &sn);
841             } else {
842                 work[(i__ << 1) - 1] = cs;
843                 work[i__ * 2] = sn;
844             }
845 /* L10: */
846         }
847         if (*nrhs > 1) {
848             i__1 = *nrhs;
849             for (i__ = 1; i__ <= i__1; ++i__) {
850                 i__2 = *n - 1;
851                 for (j = 1; j <= i__2; ++j) {
852                     cs = work[(j << 1) - 1];
853                     sn = work[j * 2];
854                     drot_(&c__1, &b[j + i__ * b_dim1], &c__1, &b[j + 1 + i__ *
855                              b_dim1], &c__1, &cs, &sn);
856 /* L20: */
857                 }
858 /* L30: */
859             }
860         }
861     }
862
863 /*     Scale. */
864
865     nm1 = *n - 1;
866     orgnrm = dlanst_("M", n, &d__[1], &e[1]);
867     if (orgnrm == 0.) {
868         dlaset_("A", n, nrhs, &c_b6, &c_b6, &b[b_offset], ldb);
869         return 0;
870     }
871
872     dlascl_("G", &c__0, &c__0, &orgnrm, &c_b11, n, &c__1, &d__[1], n, info);
873     dlascl_("G", &c__0, &c__0, &orgnrm, &c_b11, &nm1, &c__1, &e[1], &nm1, 
874             info);
875
876 /*     If N is smaller than the minimum divide size SMLSIZ, then solve */
877 /*     the problem with another solver. */
878
879     if (*n <= *smlsiz) {
880         nwork = *n * *n + 1;
881         dlaset_("A", n, n, &c_b6, &c_b11, &work[1], n);
882         dlasdq_("U", &c__0, n, n, &c__0, nrhs, &d__[1], &e[1], &work[1], n, &
883                 work[1], n, &b[b_offset], ldb, &work[nwork], info);
884         if (*info != 0) {
885             return 0;
886         }
887         tol = rcnd * (d__1 = d__[idamax_(n, &d__[1], &c__1)], abs(d__1));
888         i__1 = *n;
889         for (i__ = 1; i__ <= i__1; ++i__) {
890             if (d__[i__] <= tol) {
891                 dlaset_("A", &c__1, nrhs, &c_b6, &c_b6, &b[i__ + b_dim1], ldb);
892             } else {
893                 dlascl_("G", &c__0, &c__0, &d__[i__], &c_b11, &c__1, nrhs, &b[
894                         i__ + b_dim1], ldb, info);
895                 ++(*rank);
896             }
897 /* L40: */
898         }
899         dgemm_("T", "N", n, nrhs, n, &c_b11, &work[1], n, &b[b_offset], ldb, &
900                 c_b6, &work[nwork], n);
901         dlacpy_("A", n, nrhs, &work[nwork], n, &b[b_offset], ldb);
902
903 /*        Unscale. */
904
905         dlascl_("G", &c__0, &c__0, &c_b11, &orgnrm, n, &c__1, &d__[1], n, 
906                 info);
907         dlasrt_("D", n, &d__[1], info);
908         dlascl_("G", &c__0, &c__0, &orgnrm, &c_b11, n, nrhs, &b[b_offset], 
909                 ldb, info);
910
911         return 0;
912     }
913
914 /*     Book-keeping and setting up some constants. */
915
916     nlvl = (integer) (log((doublereal) (*n) / (doublereal) (*smlsiz + 1)) / 
917             log(2.)) + 1;
918
919     smlszp = *smlsiz + 1;
920
921     u = 1;
922     vt = *smlsiz * *n + 1;
923     difl = vt + smlszp * *n;
924     difr = difl + nlvl * *n;
925     z__ = difr + (nlvl * *n << 1);
926     c__ = z__ + nlvl * *n;
927     s = c__ + *n;
928     poles = s + *n;
929     givnum = poles + (nlvl << 1) * *n;
930     bx = givnum + (nlvl << 1) * *n;
931     nwork = bx + *n * *nrhs;
932
933     sizei = *n + 1;
934     k = sizei + *n;
935     givptr = k + *n;
936     perm = givptr + *n;
937     givcol = perm + nlvl * *n;
938     iwk = givcol + (nlvl * *n << 1);
939
940     st = 1;
941     sqre = 0;
942     icmpq1 = 1;
943     icmpq2 = 0;
944     nsub = 0;
945
946     i__1 = *n;
947     for (i__ = 1; i__ <= i__1; ++i__) {
948         if ((d__1 = d__[i__], abs(d__1)) < eps) {
949             d__[i__] = d_sign(&eps, &d__[i__]);
950         }
951 /* L50: */
952     }
953
954     i__1 = nm1;
955     for (i__ = 1; i__ <= i__1; ++i__) {
956         if ((d__1 = e[i__], abs(d__1)) < eps || i__ == nm1) {
957             ++nsub;
958             iwork[nsub] = st;
959
960 /*           Subproblem found. First determine its size and then */
961 /*           apply divide and conquer on it. */
962
963             if (i__ < nm1) {
964
965 /*              A subproblem with E(I) small for I < NM1. */
966
967                 nsize = i__ - st + 1;
968                 iwork[sizei + nsub - 1] = nsize;
969             } else if ((d__1 = e[i__], abs(d__1)) >= eps) {
970
971 /*              A subproblem with E(NM1) not too small but I = NM1. */
972
973                 nsize = *n - st + 1;
974                 iwork[sizei + nsub - 1] = nsize;
975             } else {
976
977 /*              A subproblem with E(NM1) small. This implies an */
978 /*              1-by-1 subproblem at D(N), which is not solved */
979 /*              explicitly. */
980
981                 nsize = i__ - st + 1;
982                 iwork[sizei + nsub - 1] = nsize;
983                 ++nsub;
984                 iwork[nsub] = *n;
985                 iwork[sizei + nsub - 1] = 1;
986                 dcopy_(nrhs, &b[*n + b_dim1], ldb, &work[bx + nm1], n);
987             }
988             st1 = st - 1;
989             if (nsize == 1) {
990
991 /*              This is a 1-by-1 subproblem and is not solved */
992 /*              explicitly. */
993
994                 dcopy_(nrhs, &b[st + b_dim1], ldb, &work[bx + st1], n);
995             } else if (nsize <= *smlsiz) {
996
997 /*              This is a small subproblem and is solved by DLASDQ. */
998
999                 dlaset_("A", &nsize, &nsize, &c_b6, &c_b11, &work[vt + st1], 
1000                         n);
1001                 dlasdq_("U", &c__0, &nsize, &nsize, &c__0, nrhs, &d__[st], &e[
1002                         st], &work[vt + st1], n, &work[nwork], n, &b[st + 
1003                         b_dim1], ldb, &work[nwork], info);
1004                 if (*info != 0) {
1005                     return 0;
1006                 }
1007                 dlacpy_("A", &nsize, nrhs, &b[st + b_dim1], ldb, &work[bx + 
1008                         st1], n);
1009             } else {
1010
1011 /*              A large problem. Solve it using divide and conquer. */
1012
1013                 dlasda_(&icmpq1, smlsiz, &nsize, &sqre, &d__[st], &e[st], &
1014                         work[u + st1], n, &work[vt + st1], &iwork[k + st1], &
1015                         work[difl + st1], &work[difr + st1], &work[z__ + st1],
1016                          &work[poles + st1], &iwork[givptr + st1], &iwork[
1017                         givcol + st1], n, &iwork[perm + st1], &work[givnum + 
1018                         st1], &work[c__ + st1], &work[s + st1], &work[nwork], 
1019                         &iwork[iwk], info);
1020                 if (*info != 0) {
1021                     return 0;
1022                 }
1023                 bxst = bx + st1;
1024                 dlalsa_(&icmpq2, smlsiz, &nsize, nrhs, &b[st + b_dim1], ldb, &
1025                         work[bxst], n, &work[u + st1], n, &work[vt + st1], &
1026                         iwork[k + st1], &work[difl + st1], &work[difr + st1], 
1027                         &work[z__ + st1], &work[poles + st1], &iwork[givptr + 
1028                         st1], &iwork[givcol + st1], n, &iwork[perm + st1], &
1029                         work[givnum + st1], &work[c__ + st1], &work[s + st1], 
1030                         &work[nwork], &iwork[iwk], info);
1031                 if (*info != 0) {
1032                     return 0;
1033                 }
1034             }
1035             st = i__ + 1;
1036         }
1037 /* L60: */
1038     }
1039
1040 /*     Apply the singular values and treat the tiny ones as zero. */
1041
1042     tol = rcnd * (d__1 = d__[idamax_(n, &d__[1], &c__1)], abs(d__1));
1043
1044     i__1 = *n;
1045     for (i__ = 1; i__ <= i__1; ++i__) {
1046
1047 /*        Some of the elements in D can be negative because 1-by-1 */
1048 /*        subproblems were not solved explicitly. */
1049
1050         if ((d__1 = d__[i__], abs(d__1)) <= tol) {
1051             dlaset_("A", &c__1, nrhs, &c_b6, &c_b6, &work[bx + i__ - 1], n);
1052         } else {
1053             ++(*rank);
1054             dlascl_("G", &c__0, &c__0, &d__[i__], &c_b11, &c__1, nrhs, &work[
1055                     bx + i__ - 1], n, info);
1056         }
1057         d__[i__] = (d__1 = d__[i__], abs(d__1));
1058 /* L70: */
1059     }
1060
1061 /*     Now apply back the right singular vectors. */
1062
1063     icmpq2 = 1;
1064     i__1 = nsub;
1065     for (i__ = 1; i__ <= i__1; ++i__) {
1066         st = iwork[i__];
1067         st1 = st - 1;
1068         nsize = iwork[sizei + i__ - 1];
1069         bxst = bx + st1;
1070         if (nsize == 1) {
1071             dcopy_(nrhs, &work[bxst], n, &b[st + b_dim1], ldb);
1072         } else if (nsize <= *smlsiz) {
1073             dgemm_("T", "N", &nsize, nrhs, &nsize, &c_b11, &work[vt + st1], n,
1074                      &work[bxst], n, &c_b6, &b[st + b_dim1], ldb);
1075         } else {
1076             dlalsa_(&icmpq2, smlsiz, &nsize, nrhs, &work[bxst], n, &b[st + 
1077                     b_dim1], ldb, &work[u + st1], n, &work[vt + st1], &iwork[
1078                     k + st1], &work[difl + st1], &work[difr + st1], &work[z__ 
1079                     + st1], &work[poles + st1], &iwork[givptr + st1], &iwork[
1080                     givcol + st1], n, &iwork[perm + st1], &work[givnum + st1],
1081                      &work[c__ + st1], &work[s + st1], &work[nwork], &iwork[
1082                     iwk], info);
1083             if (*info != 0) {
1084                 return 0;
1085             }
1086         }
1087 /* L80: */
1088     }
1089
1090 /*     Unscale and sort the singular values. */
1091
1092     dlascl_("G", &c__0, &c__0, &c_b11, &orgnrm, n, &c__1, &d__[1], n, info);
1093     dlasrt_("D", n, &d__[1], info);
1094     dlascl_("G", &c__0, &c__0, &orgnrm, &c_b11, n, nrhs, &b[b_offset], ldb, 
1095             info);
1096
1097     return 0;
1098
1099 /*     End of DLALSD */
1100
1101 } /* dlalsd_ */
1102