C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dlals0.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 /* Table of constant values */
513
514 static doublereal c_b5 = -1.;
515 static integer c__1 = 1;
516 static doublereal c_b11 = 1.;
517 static doublereal c_b13 = 0.;
518 static integer c__0 = 0;
519
520 /* > \brief \b DLALS0 applies back multiplying factors in solving the least squares problem using divide and c
521 onquer SVD approach. Used by sgelsd. */
522
523 /*  =========== DOCUMENTATION =========== */
524
525 /* Online html documentation available at */
526 /*            http://www.netlib.org/lapack/explore-html/ */
527
528 /* > \htmlonly */
529 /* > Download DLALS0 + dependencies */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlals0.
531 f"> */
532 /* > [TGZ]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlals0.
534 f"> */
535 /* > [ZIP]</a> */
536 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlals0.
537 f"> */
538 /* > [TXT]</a> */
539 /* > \endhtmlonly */
540
541 /*  Definition: */
542 /*  =========== */
543
544 /*       SUBROUTINE DLALS0( ICOMPQ, NL, NR, SQRE, NRHS, B, LDB, BX, LDBX, */
545 /*                          PERM, GIVPTR, GIVCOL, LDGCOL, GIVNUM, LDGNUM, */
546 /*                          POLES, DIFL, DIFR, Z, K, C, S, WORK, INFO ) */
547
548 /*       INTEGER            GIVPTR, ICOMPQ, INFO, K, LDB, LDBX, LDGCOL, */
549 /*      $                   LDGNUM, NL, NR, NRHS, SQRE */
550 /*       DOUBLE PRECISION   C, S */
551 /*       INTEGER            GIVCOL( LDGCOL, * ), PERM( * ) */
552 /*       DOUBLE PRECISION   B( LDB, * ), BX( LDBX, * ), DIFL( * ), */
553 /*      $                   DIFR( LDGNUM, * ), GIVNUM( LDGNUM, * ), */
554 /*      $                   POLES( LDGNUM, * ), WORK( * ), Z( * ) */
555
556
557 /* > \par Purpose: */
558 /*  ============= */
559 /* > */
560 /* > \verbatim */
561 /* > */
562 /* > DLALS0 applies back the multiplying factors of either the left or the */
563 /* > right singular vector matrix of a diagonal matrix appended by a row */
564 /* > to the right hand side matrix B in solving the least squares problem */
565 /* > using the divide-and-conquer SVD approach. */
566 /* > */
567 /* > For the left singular vector matrix, three types of orthogonal */
568 /* > matrices are involved: */
569 /* > */
570 /* > (1L) Givens rotations: the number of such rotations is GIVPTR; the */
571 /* >      pairs of columns/rows they were applied to are stored in GIVCOL; */
572 /* >      and the C- and S-values of these rotations are stored in GIVNUM. */
573 /* > */
574 /* > (2L) Permutation. The (NL+1)-st row of B is to be moved to the first */
575 /* >      row, and for J=2:N, PERM(J)-th row of B is to be moved to the */
576 /* >      J-th row. */
577 /* > */
578 /* > (3L) The left singular vector matrix of the remaining matrix. */
579 /* > */
580 /* > For the right singular vector matrix, four types of orthogonal */
581 /* > matrices are involved: */
582 /* > */
583 /* > (1R) The right singular vector matrix of the remaining matrix. */
584 /* > */
585 /* > (2R) If SQRE = 1, one extra Givens rotation to generate the right */
586 /* >      null space. */
587 /* > */
588 /* > (3R) The inverse transformation of (2L). */
589 /* > */
590 /* > (4R) The inverse transformation of (1L). */
591 /* > \endverbatim */
592
593 /*  Arguments: */
594 /*  ========== */
595
596 /* > \param[in] ICOMPQ */
597 /* > \verbatim */
598 /* >          ICOMPQ is INTEGER */
599 /* >         Specifies whether singular vectors are to be computed in */
600 /* >         factored form: */
601 /* >         = 0: Left singular vector matrix. */
602 /* >         = 1: Right singular vector matrix. */
603 /* > \endverbatim */
604 /* > */
605 /* > \param[in] NL */
606 /* > \verbatim */
607 /* >          NL is INTEGER */
608 /* >         The row dimension of the upper block. NL >= 1. */
609 /* > \endverbatim */
610 /* > */
611 /* > \param[in] NR */
612 /* > \verbatim */
613 /* >          NR is INTEGER */
614 /* >         The row dimension of the lower block. NR >= 1. */
615 /* > \endverbatim */
616 /* > */
617 /* > \param[in] SQRE */
618 /* > \verbatim */
619 /* >          SQRE is INTEGER */
620 /* >         = 0: the lower block is an NR-by-NR square matrix. */
621 /* >         = 1: the lower block is an NR-by-(NR+1) rectangular matrix. */
622 /* > */
623 /* >         The bidiagonal matrix has row dimension N = NL + NR + 1, */
624 /* >         and column dimension M = N + SQRE. */
625 /* > \endverbatim */
626 /* > */
627 /* > \param[in] NRHS */
628 /* > \verbatim */
629 /* >          NRHS is INTEGER */
630 /* >         The number of columns of B and BX. NRHS must be at least 1. */
631 /* > \endverbatim */
632 /* > */
633 /* > \param[in,out] B */
634 /* > \verbatim */
635 /* >          B is DOUBLE PRECISION array, dimension ( LDB, NRHS ) */
636 /* >         On input, B contains the right hand sides of the least */
637 /* >         squares problem in rows 1 through M. On output, B contains */
638 /* >         the solution X in rows 1 through N. */
639 /* > \endverbatim */
640 /* > */
641 /* > \param[in] LDB */
642 /* > \verbatim */
643 /* >          LDB is INTEGER */
644 /* >         The leading dimension of B. LDB must be at least */
645 /* >         f2cmax(1,MAX( M, N ) ). */
646 /* > \endverbatim */
647 /* > */
648 /* > \param[out] BX */
649 /* > \verbatim */
650 /* >          BX is DOUBLE PRECISION array, dimension ( LDBX, NRHS ) */
651 /* > \endverbatim */
652 /* > */
653 /* > \param[in] LDBX */
654 /* > \verbatim */
655 /* >          LDBX is INTEGER */
656 /* >         The leading dimension of BX. */
657 /* > \endverbatim */
658 /* > */
659 /* > \param[in] PERM */
660 /* > \verbatim */
661 /* >          PERM is INTEGER array, dimension ( N ) */
662 /* >         The permutations (from deflation and sorting) applied */
663 /* >         to the two blocks. */
664 /* > \endverbatim */
665 /* > */
666 /* > \param[in] GIVPTR */
667 /* > \verbatim */
668 /* >          GIVPTR is INTEGER */
669 /* >         The number of Givens rotations which took place in this */
670 /* >         subproblem. */
671 /* > \endverbatim */
672 /* > */
673 /* > \param[in] GIVCOL */
674 /* > \verbatim */
675 /* >          GIVCOL is INTEGER array, dimension ( LDGCOL, 2 ) */
676 /* >         Each pair of numbers indicates a pair of rows/columns */
677 /* >         involved in a Givens rotation. */
678 /* > \endverbatim */
679 /* > */
680 /* > \param[in] LDGCOL */
681 /* > \verbatim */
682 /* >          LDGCOL is INTEGER */
683 /* >         The leading dimension of GIVCOL, must be at least N. */
684 /* > \endverbatim */
685 /* > */
686 /* > \param[in] GIVNUM */
687 /* > \verbatim */
688 /* >          GIVNUM is DOUBLE PRECISION array, dimension ( LDGNUM, 2 ) */
689 /* >         Each number indicates the C or S value used in the */
690 /* >         corresponding Givens rotation. */
691 /* > \endverbatim */
692 /* > */
693 /* > \param[in] LDGNUM */
694 /* > \verbatim */
695 /* >          LDGNUM is INTEGER */
696 /* >         The leading dimension of arrays DIFR, POLES and */
697 /* >         GIVNUM, must be at least K. */
698 /* > \endverbatim */
699 /* > */
700 /* > \param[in] POLES */
701 /* > \verbatim */
702 /* >          POLES is DOUBLE PRECISION array, dimension ( LDGNUM, 2 ) */
703 /* >         On entry, POLES(1:K, 1) contains the new singular */
704 /* >         values obtained from solving the secular equation, and */
705 /* >         POLES(1:K, 2) is an array containing the poles in the secular */
706 /* >         equation. */
707 /* > \endverbatim */
708 /* > */
709 /* > \param[in] DIFL */
710 /* > \verbatim */
711 /* >          DIFL is DOUBLE PRECISION array, dimension ( K ). */
712 /* >         On entry, DIFL(I) is the distance between I-th updated */
713 /* >         (undeflated) singular value and the I-th (undeflated) old */
714 /* >         singular value. */
715 /* > \endverbatim */
716 /* > */
717 /* > \param[in] DIFR */
718 /* > \verbatim */
719 /* >          DIFR is DOUBLE PRECISION array, dimension ( LDGNUM, 2 ). */
720 /* >         On entry, DIFR(I, 1) contains the distances between I-th */
721 /* >         updated (undeflated) singular value and the I+1-th */
722 /* >         (undeflated) old singular value. And DIFR(I, 2) is the */
723 /* >         normalizing factor for the I-th right singular vector. */
724 /* > \endverbatim */
725 /* > */
726 /* > \param[in] Z */
727 /* > \verbatim */
728 /* >          Z is DOUBLE PRECISION array, dimension ( K ) */
729 /* >         Contain the components of the deflation-adjusted updating row */
730 /* >         vector. */
731 /* > \endverbatim */
732 /* > */
733 /* > \param[in] K */
734 /* > \verbatim */
735 /* >          K is INTEGER */
736 /* >         Contains the dimension of the non-deflated matrix, */
737 /* >         This is the order of the related secular equation. 1 <= K <=N. */
738 /* > \endverbatim */
739 /* > */
740 /* > \param[in] C */
741 /* > \verbatim */
742 /* >          C is DOUBLE PRECISION */
743 /* >         C contains garbage if SQRE =0 and the C-value of a Givens */
744 /* >         rotation related to the right null space if SQRE = 1. */
745 /* > \endverbatim */
746 /* > */
747 /* > \param[in] S */
748 /* > \verbatim */
749 /* >          S is DOUBLE PRECISION */
750 /* >         S contains garbage if SQRE =0 and the S-value of a Givens */
751 /* >         rotation related to the right null space if SQRE = 1. */
752 /* > \endverbatim */
753 /* > */
754 /* > \param[out] WORK */
755 /* > \verbatim */
756 /* >          WORK is DOUBLE PRECISION array, dimension ( K ) */
757 /* > \endverbatim */
758 /* > */
759 /* > \param[out] INFO */
760 /* > \verbatim */
761 /* >          INFO is INTEGER */
762 /* >          = 0:  successful exit. */
763 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
764 /* > \endverbatim */
765
766 /*  Authors: */
767 /*  ======== */
768
769 /* > \author Univ. of Tennessee */
770 /* > \author Univ. of California Berkeley */
771 /* > \author Univ. of Colorado Denver */
772 /* > \author NAG Ltd. */
773
774 /* > \date December 2016 */
775
776 /* > \ingroup doubleOTHERcomputational */
777
778 /* > \par Contributors: */
779 /*  ================== */
780 /* > */
781 /* >     Ming Gu and Ren-Cang Li, Computer Science Division, University of */
782 /* >       California at Berkeley, USA \n */
783 /* >     Osni Marques, LBNL/NERSC, USA \n */
784
785 /*  ===================================================================== */
786 /* Subroutine */ int dlals0_(integer *icompq, integer *nl, integer *nr, 
787         integer *sqre, integer *nrhs, doublereal *b, integer *ldb, doublereal 
788         *bx, integer *ldbx, integer *perm, integer *givptr, integer *givcol, 
789         integer *ldgcol, doublereal *givnum, integer *ldgnum, doublereal *
790         poles, doublereal *difl, doublereal *difr, doublereal *z__, integer *
791         k, doublereal *c__, doublereal *s, doublereal *work, integer *info)
792 {
793     /* System generated locals */
794     integer givcol_dim1, givcol_offset, b_dim1, b_offset, bx_dim1, bx_offset, 
795             difr_dim1, difr_offset, givnum_dim1, givnum_offset, poles_dim1, 
796             poles_offset, i__1, i__2;
797     doublereal d__1;
798
799     /* Local variables */
800     doublereal temp;
801     extern /* Subroutine */ int drot_(integer *, doublereal *, integer *, 
802             doublereal *, integer *, doublereal *, doublereal *);
803     extern doublereal dnrm2_(integer *, doublereal *, integer *);
804     integer i__, j, m, n;
805     extern /* Subroutine */ int dscal_(integer *, doublereal *, doublereal *, 
806             integer *);
807     doublereal diflj, difrj, dsigj;
808     extern /* Subroutine */ int dgemv_(char *, integer *, integer *, 
809             doublereal *, doublereal *, integer *, doublereal *, integer *, 
810             doublereal *, doublereal *, integer *), dcopy_(integer *, 
811             doublereal *, integer *, doublereal *, integer *);
812     extern doublereal dlamc3_(doublereal *, doublereal *);
813     doublereal dj;
814     extern /* Subroutine */ int dlascl_(char *, integer *, integer *, 
815             doublereal *, doublereal *, integer *, integer *, doublereal *, 
816             integer *, integer *), dlacpy_(char *, integer *, integer 
817             *, doublereal *, integer *, doublereal *, integer *), 
818             xerbla_(char *, integer *, ftnlen);
819     doublereal dsigjp;
820     integer nlp1;
821
822
823 /*  -- LAPACK computational routine (version 3.7.0) -- */
824 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
825 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
826 /*     December 2016 */
827
828
829 /*  ===================================================================== */
830
831
832 /*     Test the input parameters. */
833
834     /* Parameter adjustments */
835     b_dim1 = *ldb;
836     b_offset = 1 + b_dim1 * 1;
837     b -= b_offset;
838     bx_dim1 = *ldbx;
839     bx_offset = 1 + bx_dim1 * 1;
840     bx -= bx_offset;
841     --perm;
842     givcol_dim1 = *ldgcol;
843     givcol_offset = 1 + givcol_dim1 * 1;
844     givcol -= givcol_offset;
845     difr_dim1 = *ldgnum;
846     difr_offset = 1 + difr_dim1 * 1;
847     difr -= difr_offset;
848     poles_dim1 = *ldgnum;
849     poles_offset = 1 + poles_dim1 * 1;
850     poles -= poles_offset;
851     givnum_dim1 = *ldgnum;
852     givnum_offset = 1 + givnum_dim1 * 1;
853     givnum -= givnum_offset;
854     --difl;
855     --z__;
856     --work;
857
858     /* Function Body */
859     *info = 0;
860     n = *nl + *nr + 1;
861
862     if (*icompq < 0 || *icompq > 1) {
863         *info = -1;
864     } else if (*nl < 1) {
865         *info = -2;
866     } else if (*nr < 1) {
867         *info = -3;
868     } else if (*sqre < 0 || *sqre > 1) {
869         *info = -4;
870     } else if (*nrhs < 1) {
871         *info = -5;
872     } else if (*ldb < n) {
873         *info = -7;
874     } else if (*ldbx < n) {
875         *info = -9;
876     } else if (*givptr < 0) {
877         *info = -11;
878     } else if (*ldgcol < n) {
879         *info = -13;
880     } else if (*ldgnum < n) {
881         *info = -15;
882     //} else if (*k < 1) {
883     } else if (*k < 0) {
884         *info = -20;
885     }
886     if (*info != 0) {
887         i__1 = -(*info);
888         xerbla_("DLALS0", &i__1, (ftnlen)6);
889         return 0;
890     }
891
892     m = n + *sqre;
893     nlp1 = *nl + 1;
894
895     if (*icompq == 0) {
896
897 /*        Apply back orthogonal transformations from the left. */
898
899 /*        Step (1L): apply back the Givens rotations performed. */
900
901         i__1 = *givptr;
902         for (i__ = 1; i__ <= i__1; ++i__) {
903             drot_(nrhs, &b[givcol[i__ + (givcol_dim1 << 1)] + b_dim1], ldb, &
904                     b[givcol[i__ + givcol_dim1] + b_dim1], ldb, &givnum[i__ + 
905                     (givnum_dim1 << 1)], &givnum[i__ + givnum_dim1]);
906 /* L10: */
907         }
908
909 /*        Step (2L): permute rows of B. */
910
911         dcopy_(nrhs, &b[nlp1 + b_dim1], ldb, &bx[bx_dim1 + 1], ldbx);
912         i__1 = n;
913         for (i__ = 2; i__ <= i__1; ++i__) {
914             dcopy_(nrhs, &b[perm[i__] + b_dim1], ldb, &bx[i__ + bx_dim1], 
915                     ldbx);
916 /* L20: */
917         }
918
919 /*        Step (3L): apply the inverse of the left singular vector */
920 /*        matrix to BX. */
921
922         if (*k == 1) {
923             dcopy_(nrhs, &bx[bx_offset], ldbx, &b[b_offset], ldb);
924             if (z__[1] < 0.) {
925                 dscal_(nrhs, &c_b5, &b[b_offset], ldb);
926             }
927         } else {
928             i__1 = *k;
929             for (j = 1; j <= i__1; ++j) {
930                 diflj = difl[j];
931                 dj = poles[j + poles_dim1];
932                 dsigj = -poles[j + (poles_dim1 << 1)];
933                 if (j < *k) {
934                     difrj = -difr[j + difr_dim1];
935                     dsigjp = -poles[j + 1 + (poles_dim1 << 1)];
936                 }
937                 if (z__[j] == 0. || poles[j + (poles_dim1 << 1)] == 0.) {
938                     work[j] = 0.;
939                 } else {
940                     work[j] = -poles[j + (poles_dim1 << 1)] * z__[j] / diflj /
941                              (poles[j + (poles_dim1 << 1)] + dj);
942                 }
943                 i__2 = j - 1;
944                 for (i__ = 1; i__ <= i__2; ++i__) {
945                     if (z__[i__] == 0. || poles[i__ + (poles_dim1 << 1)] == 
946                             0.) {
947                         work[i__] = 0.;
948                     } else {
949                         work[i__] = poles[i__ + (poles_dim1 << 1)] * z__[i__] 
950                                 / (dlamc3_(&poles[i__ + (poles_dim1 << 1)], &
951                                 dsigj) - diflj) / (poles[i__ + (poles_dim1 << 
952                                 1)] + dj);
953                     }
954 /* L30: */
955                 }
956                 i__2 = *k;
957                 for (i__ = j + 1; i__ <= i__2; ++i__) {
958                     if (z__[i__] == 0. || poles[i__ + (poles_dim1 << 1)] == 
959                             0.) {
960                         work[i__] = 0.;
961                     } else {
962                         work[i__] = poles[i__ + (poles_dim1 << 1)] * z__[i__] 
963                                 / (dlamc3_(&poles[i__ + (poles_dim1 << 1)], &
964                                 dsigjp) + difrj) / (poles[i__ + (poles_dim1 <<
965                                  1)] + dj);
966                     }
967 /* L40: */
968                 }
969                 work[1] = -1.;
970                 temp = dnrm2_(k, &work[1], &c__1);
971                 dgemv_("T", k, nrhs, &c_b11, &bx[bx_offset], ldbx, &work[1], &
972                         c__1, &c_b13, &b[j + b_dim1], ldb);
973                 dlascl_("G", &c__0, &c__0, &temp, &c_b11, &c__1, nrhs, &b[j + 
974                         b_dim1], ldb, info);
975 /* L50: */
976             }
977         }
978
979 /*        Move the deflated rows of BX to B also. */
980
981         if (*k < f2cmax(m,n)) {
982             i__1 = n - *k;
983             dlacpy_("A", &i__1, nrhs, &bx[*k + 1 + bx_dim1], ldbx, &b[*k + 1 
984                     + b_dim1], ldb);
985         }
986     } else {
987
988 /*        Apply back the right orthogonal transformations. */
989
990 /*        Step (1R): apply back the new right singular vector matrix */
991 /*        to B. */
992
993         if (*k == 1) {
994             dcopy_(nrhs, &b[b_offset], ldb, &bx[bx_offset], ldbx);
995         } else {
996             i__1 = *k;
997             for (j = 1; j <= i__1; ++j) {
998                 dsigj = poles[j + (poles_dim1 << 1)];
999                 if (z__[j] == 0.) {
1000                     work[j] = 0.;
1001                 } else {
1002                     work[j] = -z__[j] / difl[j] / (dsigj + poles[j + 
1003                             poles_dim1]) / difr[j + (difr_dim1 << 1)];
1004                 }
1005                 i__2 = j - 1;
1006                 for (i__ = 1; i__ <= i__2; ++i__) {
1007                     if (z__[j] == 0.) {
1008                         work[i__] = 0.;
1009                     } else {
1010                         d__1 = -poles[i__ + 1 + (poles_dim1 << 1)];
1011                         work[i__] = z__[j] / (dlamc3_(&dsigj, &d__1) - difr[
1012                                 i__ + difr_dim1]) / (dsigj + poles[i__ + 
1013                                 poles_dim1]) / difr[i__ + (difr_dim1 << 1)];
1014                     }
1015 /* L60: */
1016                 }
1017                 i__2 = *k;
1018                 for (i__ = j + 1; i__ <= i__2; ++i__) {
1019                     if (z__[j] == 0.) {
1020                         work[i__] = 0.;
1021                     } else {
1022                         d__1 = -poles[i__ + (poles_dim1 << 1)];
1023                         work[i__] = z__[j] / (dlamc3_(&dsigj, &d__1) - difl[
1024                                 i__]) / (dsigj + poles[i__ + poles_dim1]) / 
1025                                 difr[i__ + (difr_dim1 << 1)];
1026                     }
1027 /* L70: */
1028                 }
1029                 dgemv_("T", k, nrhs, &c_b11, &b[b_offset], ldb, &work[1], &
1030                         c__1, &c_b13, &bx[j + bx_dim1], ldbx);
1031 /* L80: */
1032             }
1033         }
1034
1035 /*        Step (2R): if SQRE = 1, apply back the rotation that is */
1036 /*        related to the right null space of the subproblem. */
1037
1038         if (*sqre == 1) {
1039             dcopy_(nrhs, &b[m + b_dim1], ldb, &bx[m + bx_dim1], ldbx);
1040             drot_(nrhs, &bx[bx_dim1 + 1], ldbx, &bx[m + bx_dim1], ldbx, c__, 
1041                     s);
1042         }
1043         if (*k < f2cmax(m,n)) {
1044             i__1 = n - *k;
1045             dlacpy_("A", &i__1, nrhs, &b[*k + 1 + b_dim1], ldb, &bx[*k + 1 + 
1046                     bx_dim1], ldbx);
1047         }
1048
1049 /*        Step (3R): permute rows of B. */
1050
1051         dcopy_(nrhs, &bx[bx_dim1 + 1], ldbx, &b[nlp1 + b_dim1], ldb);
1052         if (*sqre == 1) {
1053             dcopy_(nrhs, &bx[m + bx_dim1], ldbx, &b[m + b_dim1], ldb);
1054         }
1055         i__1 = n;
1056         for (i__ = 2; i__ <= i__1; ++i__) {
1057             dcopy_(nrhs, &bx[i__ + bx_dim1], ldbx, &b[perm[i__] + b_dim1], 
1058                     ldb);
1059 /* L90: */
1060         }
1061
1062 /*        Step (4R): apply back the Givens rotations performed. */
1063
1064         for (i__ = *givptr; i__ >= 1; --i__) {
1065             d__1 = -givnum[i__ + givnum_dim1];
1066             drot_(nrhs, &b[givcol[i__ + (givcol_dim1 << 1)] + b_dim1], ldb, &
1067                     b[givcol[i__ + givcol_dim1] + b_dim1], ldb, &givnum[i__ + 
1068                     (givnum_dim1 << 1)], &d__1);
1069 /* L100: */
1070         }
1071     }
1072
1073     return 0;
1074
1075 /*     End of DLALS0 */
1076
1077 } /* dlals0_ */
1078