C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / clalsa.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 real c_b9 = 1.f;
516 static real c_b10 = 0.f;
517 static integer c__2 = 2;
518
519 /* > \brief \b CLALSA computes the SVD of the coefficient matrix in compact form. Used by sgelsd. */
520
521 /*  =========== DOCUMENTATION =========== */
522
523 /* Online html documentation available at */
524 /*            http://www.netlib.org/lapack/explore-html/ */
525
526 /* > \htmlonly */
527 /* > Download CLALSA + dependencies */
528 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/clalsa.
529 f"> */
530 /* > [TGZ]</a> */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/clalsa.
532 f"> */
533 /* > [ZIP]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/clalsa.
535 f"> */
536 /* > [TXT]</a> */
537 /* > \endhtmlonly */
538
539 /*  Definition: */
540 /*  =========== */
541
542 /*       SUBROUTINE CLALSA( ICOMPQ, SMLSIZ, N, NRHS, B, LDB, BX, LDBX, U, */
543 /*                          LDU, VT, K, DIFL, DIFR, Z, POLES, GIVPTR, */
544 /*                          GIVCOL, LDGCOL, PERM, GIVNUM, C, S, RWORK, */
545 /*                          IWORK, INFO ) */
546
547 /*       INTEGER            ICOMPQ, INFO, LDB, LDBX, LDGCOL, LDU, N, NRHS, */
548 /*      $                   SMLSIZ */
549 /*       INTEGER            GIVCOL( LDGCOL, * ), GIVPTR( * ), IWORK( * ), */
550 /*      $                   K( * ), PERM( LDGCOL, * ) */
551 /*       REAL               C( * ), DIFL( LDU, * ), DIFR( LDU, * ), */
552 /*      $                   GIVNUM( LDU, * ), POLES( LDU, * ), RWORK( * ), */
553 /*      $                   S( * ), U( LDU, * ), VT( LDU, * ), Z( LDU, * ) */
554 /*       COMPLEX            B( LDB, * ), BX( LDBX, * ) */
555
556
557 /* > \par Purpose: */
558 /*  ============= */
559 /* > */
560 /* > \verbatim */
561 /* > */
562 /* > CLALSA is an itermediate step in solving the least squares problem */
563 /* > by computing the SVD of the coefficient matrix in compact form (The */
564 /* > singular vectors are computed as products of simple orthorgonal */
565 /* > matrices.). */
566 /* > */
567 /* > If ICOMPQ = 0, CLALSA applies the inverse of the left singular vector */
568 /* > matrix of an upper bidiagonal matrix to the right hand side; and if */
569 /* > ICOMPQ = 1, CLALSA applies the right singular vector matrix to the */
570 /* > right hand side. The singular vector matrices were generated in */
571 /* > compact form by CLALSA. */
572 /* > \endverbatim */
573
574 /*  Arguments: */
575 /*  ========== */
576
577 /* > \param[in] ICOMPQ */
578 /* > \verbatim */
579 /* >          ICOMPQ is INTEGER */
580 /* >         Specifies whether the left or the right singular vector */
581 /* >         matrix is involved. */
582 /* >         = 0: Left singular vector matrix */
583 /* >         = 1: Right singular vector 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 row and column dimensions of the upper bidiagonal matrix. */
597 /* > \endverbatim */
598 /* > */
599 /* > \param[in] NRHS */
600 /* > \verbatim */
601 /* >          NRHS is INTEGER */
602 /* >         The number of columns of B and BX. NRHS must be at least 1. */
603 /* > \endverbatim */
604 /* > */
605 /* > \param[in,out] B */
606 /* > \verbatim */
607 /* >          B is COMPLEX array, dimension ( LDB, NRHS ) */
608 /* >         On input, B contains the right hand sides of the least */
609 /* >         squares problem in rows 1 through M. */
610 /* >         On output, B contains the solution X in rows 1 through N. */
611 /* > \endverbatim */
612 /* > */
613 /* > \param[in] LDB */
614 /* > \verbatim */
615 /* >          LDB is INTEGER */
616 /* >         The leading dimension of B in the calling subprogram. */
617 /* >         LDB must be at least f2cmax(1,MAX( M, N ) ). */
618 /* > \endverbatim */
619 /* > */
620 /* > \param[out] BX */
621 /* > \verbatim */
622 /* >          BX is COMPLEX array, dimension ( LDBX, NRHS ) */
623 /* >         On exit, the result of applying the left or right singular */
624 /* >         vector matrix to B. */
625 /* > \endverbatim */
626 /* > */
627 /* > \param[in] LDBX */
628 /* > \verbatim */
629 /* >          LDBX is INTEGER */
630 /* >         The leading dimension of BX. */
631 /* > \endverbatim */
632 /* > */
633 /* > \param[in] U */
634 /* > \verbatim */
635 /* >          U is REAL array, dimension ( LDU, SMLSIZ ). */
636 /* >         On entry, U contains the left singular vector matrices of all */
637 /* >         subproblems at the bottom level. */
638 /* > \endverbatim */
639 /* > */
640 /* > \param[in] LDU */
641 /* > \verbatim */
642 /* >          LDU is INTEGER, LDU = > N. */
643 /* >         The leading dimension of arrays U, VT, DIFL, DIFR, */
644 /* >         POLES, GIVNUM, and Z. */
645 /* > \endverbatim */
646 /* > */
647 /* > \param[in] VT */
648 /* > \verbatim */
649 /* >          VT is REAL array, dimension ( LDU, SMLSIZ+1 ). */
650 /* >         On entry, VT**H contains the right singular vector matrices of */
651 /* >         all subproblems at the bottom level. */
652 /* > \endverbatim */
653 /* > */
654 /* > \param[in] K */
655 /* > \verbatim */
656 /* >          K is INTEGER array, dimension ( N ). */
657 /* > \endverbatim */
658 /* > */
659 /* > \param[in] DIFL */
660 /* > \verbatim */
661 /* >          DIFL is REAL array, dimension ( LDU, NLVL ). */
662 /* >         where NLVL = INT(log_2 (N/(SMLSIZ+1))) + 1. */
663 /* > \endverbatim */
664 /* > */
665 /* > \param[in] DIFR */
666 /* > \verbatim */
667 /* >          DIFR is REAL array, dimension ( LDU, 2 * NLVL ). */
668 /* >         On entry, DIFL(*, I) and DIFR(*, 2 * I -1) record */
669 /* >         distances between singular values on the I-th level and */
670 /* >         singular values on the (I -1)-th level, and DIFR(*, 2 * I) */
671 /* >         record the normalizing factors of the right singular vectors */
672 /* >         matrices of subproblems on I-th level. */
673 /* > \endverbatim */
674 /* > */
675 /* > \param[in] Z */
676 /* > \verbatim */
677 /* >          Z is REAL array, dimension ( LDU, NLVL ). */
678 /* >         On entry, Z(1, I) contains the components of the deflation- */
679 /* >         adjusted updating row vector for subproblems on the I-th */
680 /* >         level. */
681 /* > \endverbatim */
682 /* > */
683 /* > \param[in] POLES */
684 /* > \verbatim */
685 /* >          POLES is REAL array, dimension ( LDU, 2 * NLVL ). */
686 /* >         On entry, POLES(*, 2 * I -1: 2 * I) contains the new and old */
687 /* >         singular values involved in the secular equations on the I-th */
688 /* >         level. */
689 /* > \endverbatim */
690 /* > */
691 /* > \param[in] GIVPTR */
692 /* > \verbatim */
693 /* >          GIVPTR is INTEGER array, dimension ( N ). */
694 /* >         On entry, GIVPTR( I ) records the number of Givens */
695 /* >         rotations performed on the I-th problem on the computation */
696 /* >         tree. */
697 /* > \endverbatim */
698 /* > */
699 /* > \param[in] GIVCOL */
700 /* > \verbatim */
701 /* >          GIVCOL is INTEGER array, dimension ( LDGCOL, 2 * NLVL ). */
702 /* >         On entry, for each I, GIVCOL(*, 2 * I - 1: 2 * I) records the */
703 /* >         locations of Givens rotations performed on the I-th level on */
704 /* >         the computation tree. */
705 /* > \endverbatim */
706 /* > */
707 /* > \param[in] LDGCOL */
708 /* > \verbatim */
709 /* >          LDGCOL is INTEGER, LDGCOL = > N. */
710 /* >         The leading dimension of arrays GIVCOL and PERM. */
711 /* > \endverbatim */
712 /* > */
713 /* > \param[in] PERM */
714 /* > \verbatim */
715 /* >          PERM is INTEGER array, dimension ( LDGCOL, NLVL ). */
716 /* >         On entry, PERM(*, I) records permutations done on the I-th */
717 /* >         level of the computation tree. */
718 /* > \endverbatim */
719 /* > */
720 /* > \param[in] GIVNUM */
721 /* > \verbatim */
722 /* >          GIVNUM is REAL array, dimension ( LDU, 2 * NLVL ). */
723 /* >         On entry, GIVNUM(*, 2 *I -1 : 2 * I) records the C- and S- */
724 /* >         values of Givens rotations performed on the I-th level on the */
725 /* >         computation tree. */
726 /* > \endverbatim */
727 /* > */
728 /* > \param[in] C */
729 /* > \verbatim */
730 /* >          C is REAL array, dimension ( N ). */
731 /* >         On entry, if the I-th subproblem is not square, */
732 /* >         C( I ) contains the C-value of a Givens rotation related to */
733 /* >         the right null space of the I-th subproblem. */
734 /* > \endverbatim */
735 /* > */
736 /* > \param[in] S */
737 /* > \verbatim */
738 /* >          S is REAL array, dimension ( N ). */
739 /* >         On entry, if the I-th subproblem is not square, */
740 /* >         S( I ) contains the S-value of a Givens rotation related to */
741 /* >         the right null space of the I-th subproblem. */
742 /* > \endverbatim */
743 /* > */
744 /* > \param[out] RWORK */
745 /* > \verbatim */
746 /* >          RWORK is REAL array, dimension at least */
747 /* >         MAX( (SMLSZ+1)*NRHS*3, N*(1+NRHS) + 2*NRHS ). */
748 /* > \endverbatim */
749 /* > */
750 /* > \param[out] IWORK */
751 /* > \verbatim */
752 /* >          IWORK is INTEGER array, dimension (3*N) */
753 /* > \endverbatim */
754 /* > */
755 /* > \param[out] INFO */
756 /* > \verbatim */
757 /* >          INFO is INTEGER */
758 /* >          = 0:  successful exit. */
759 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
760 /* > \endverbatim */
761
762 /*  Authors: */
763 /*  ======== */
764
765 /* > \author Univ. of Tennessee */
766 /* > \author Univ. of California Berkeley */
767 /* > \author Univ. of Colorado Denver */
768 /* > \author NAG Ltd. */
769
770 /* > \date June 2017 */
771
772 /* > \ingroup complexOTHERcomputational */
773
774 /* > \par Contributors: */
775 /*  ================== */
776 /* > */
777 /* >     Ming Gu and Ren-Cang Li, Computer Science Division, University of */
778 /* >       California at Berkeley, USA \n */
779 /* >     Osni Marques, LBNL/NERSC, USA \n */
780
781 /*  ===================================================================== */
782 /* Subroutine */ int clalsa_(integer *icompq, integer *smlsiz, integer *n, 
783         integer *nrhs, complex *b, integer *ldb, complex *bx, integer *ldbx, 
784         real *u, integer *ldu, real *vt, integer *k, real *difl, real *difr, 
785         real *z__, real *poles, integer *givptr, integer *givcol, integer *
786         ldgcol, integer *perm, real *givnum, real *c__, real *s, real *rwork, 
787         integer *iwork, integer *info)
788 {
789     /* System generated locals */
790     integer givcol_dim1, givcol_offset, perm_dim1, perm_offset, difl_dim1, 
791             difl_offset, difr_dim1, difr_offset, givnum_dim1, givnum_offset, 
792             poles_dim1, poles_offset, u_dim1, u_offset, vt_dim1, vt_offset, 
793             z_dim1, z_offset, b_dim1, b_offset, bx_dim1, bx_offset, i__1, 
794             i__2, i__3, i__4, i__5, i__6;
795     complex q__1;
796
797     /* Local variables */
798     integer jcol, nlvl, sqre, jrow, i__, j, jimag, jreal, inode, ndiml;
799     extern /* Subroutine */ int sgemm_(char *, char *, integer *, integer *, 
800             integer *, real *, real *, integer *, real *, integer *, real *, 
801             real *, integer *);
802     integer ndimr;
803     extern /* Subroutine */ int ccopy_(integer *, complex *, integer *, 
804             complex *, integer *);
805     integer i1;
806     extern /* Subroutine */ int clals0_(integer *, integer *, integer *, 
807             integer *, integer *, complex *, integer *, complex *, integer *, 
808             integer *, integer *, integer *, integer *, real *, integer *, 
809             real *, real *, real *, real *, integer *, real *, real *, real *,
810              integer *);
811     integer ic, lf, nd, ll, nl, nr;
812     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen), slasdt_(
813             integer *, integer *, integer *, integer *, integer *, integer *, 
814             integer *);
815     integer im1, nlf, nrf, lvl, ndb1, nlp1, lvl2, nrp1;
816
817
818 /*  -- LAPACK computational routine (version 3.7.1) -- */
819 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
820 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
821 /*     June 2017 */
822
823
824 /*  ===================================================================== */
825
826
827 /*     Test the input parameters. */
828
829     /* Parameter adjustments */
830     b_dim1 = *ldb;
831     b_offset = 1 + b_dim1 * 1;
832     b -= b_offset;
833     bx_dim1 = *ldbx;
834     bx_offset = 1 + bx_dim1 * 1;
835     bx -= bx_offset;
836     givnum_dim1 = *ldu;
837     givnum_offset = 1 + givnum_dim1 * 1;
838     givnum -= givnum_offset;
839     poles_dim1 = *ldu;
840     poles_offset = 1 + poles_dim1 * 1;
841     poles -= poles_offset;
842     z_dim1 = *ldu;
843     z_offset = 1 + z_dim1 * 1;
844     z__ -= z_offset;
845     difr_dim1 = *ldu;
846     difr_offset = 1 + difr_dim1 * 1;
847     difr -= difr_offset;
848     difl_dim1 = *ldu;
849     difl_offset = 1 + difl_dim1 * 1;
850     difl -= difl_offset;
851     vt_dim1 = *ldu;
852     vt_offset = 1 + vt_dim1 * 1;
853     vt -= vt_offset;
854     u_dim1 = *ldu;
855     u_offset = 1 + u_dim1 * 1;
856     u -= u_offset;
857     --k;
858     --givptr;
859     perm_dim1 = *ldgcol;
860     perm_offset = 1 + perm_dim1 * 1;
861     perm -= perm_offset;
862     givcol_dim1 = *ldgcol;
863     givcol_offset = 1 + givcol_dim1 * 1;
864     givcol -= givcol_offset;
865     --c__;
866     --s;
867     --rwork;
868     --iwork;
869
870     /* Function Body */
871     *info = 0;
872
873     if (*icompq < 0 || *icompq > 1) {
874         *info = -1;
875     } else if (*smlsiz < 3) {
876         *info = -2;
877     } else if (*n < *smlsiz) {
878         *info = -3;
879     } else if (*nrhs < 1) {
880         *info = -4;
881     } else if (*ldb < *n) {
882         *info = -6;
883     } else if (*ldbx < *n) {
884         *info = -8;
885     } else if (*ldu < *n) {
886         *info = -10;
887     } else if (*ldgcol < *n) {
888         *info = -19;
889     }
890     if (*info != 0) {
891         i__1 = -(*info);
892         xerbla_("CLALSA", &i__1, (ftnlen)6);
893         return 0;
894     }
895
896 /*     Book-keeping and  setting up the computation tree. */
897
898     inode = 1;
899     ndiml = inode + *n;
900     ndimr = ndiml + *n;
901
902     slasdt_(n, &nlvl, &nd, &iwork[inode], &iwork[ndiml], &iwork[ndimr], 
903             smlsiz);
904
905 /*     The following code applies back the left singular vector factors. */
906 /*     For applying back the right singular vector factors, go to 170. */
907
908     if (*icompq == 1) {
909         goto L170;
910     }
911
912 /*     The nodes on the bottom level of the tree were solved */
913 /*     by SLASDQ. The corresponding left and right singular vector */
914 /*     matrices are in explicit form. First apply back the left */
915 /*     singular vector matrices. */
916
917     ndb1 = (nd + 1) / 2;
918     i__1 = nd;
919     for (i__ = ndb1; i__ <= i__1; ++i__) {
920
921 /*        IC : center row of each node */
922 /*        NL : number of rows of left  subproblem */
923 /*        NR : number of rows of right subproblem */
924 /*        NLF: starting row of the left   subproblem */
925 /*        NRF: starting row of the right  subproblem */
926
927         i1 = i__ - 1;
928         ic = iwork[inode + i1];
929         nl = iwork[ndiml + i1];
930         nr = iwork[ndimr + i1];
931         nlf = ic - nl;
932         nrf = ic + 1;
933
934 /*        Since B and BX are complex, the following call to SGEMM */
935 /*        is performed in two steps (real and imaginary parts). */
936
937 /*        CALL SGEMM( 'T', 'N', NL, NRHS, NL, ONE, U( NLF, 1 ), LDU, */
938 /*     $               B( NLF, 1 ), LDB, ZERO, BX( NLF, 1 ), LDBX ) */
939
940         j = nl * *nrhs << 1;
941         i__2 = *nrhs;
942         for (jcol = 1; jcol <= i__2; ++jcol) {
943             i__3 = nlf + nl - 1;
944             for (jrow = nlf; jrow <= i__3; ++jrow) {
945                 ++j;
946                 i__4 = jrow + jcol * b_dim1;
947                 rwork[j] = b[i__4].r;
948 /* L10: */
949             }
950 /* L20: */
951         }
952         sgemm_("T", "N", &nl, nrhs, &nl, &c_b9, &u[nlf + u_dim1], ldu, &rwork[
953                 (nl * *nrhs << 1) + 1], &nl, &c_b10, &rwork[1], &nl);
954         j = nl * *nrhs << 1;
955         i__2 = *nrhs;
956         for (jcol = 1; jcol <= i__2; ++jcol) {
957             i__3 = nlf + nl - 1;
958             for (jrow = nlf; jrow <= i__3; ++jrow) {
959                 ++j;
960                 rwork[j] = r_imag(&b[jrow + jcol * b_dim1]);
961 /* L30: */
962             }
963 /* L40: */
964         }
965         sgemm_("T", "N", &nl, nrhs, &nl, &c_b9, &u[nlf + u_dim1], ldu, &rwork[
966                 (nl * *nrhs << 1) + 1], &nl, &c_b10, &rwork[nl * *nrhs + 1], &
967                 nl);
968         jreal = 0;
969         jimag = nl * *nrhs;
970         i__2 = *nrhs;
971         for (jcol = 1; jcol <= i__2; ++jcol) {
972             i__3 = nlf + nl - 1;
973             for (jrow = nlf; jrow <= i__3; ++jrow) {
974                 ++jreal;
975                 ++jimag;
976                 i__4 = jrow + jcol * bx_dim1;
977                 i__5 = jreal;
978                 i__6 = jimag;
979                 q__1.r = rwork[i__5], q__1.i = rwork[i__6];
980                 bx[i__4].r = q__1.r, bx[i__4].i = q__1.i;
981 /* L50: */
982             }
983 /* L60: */
984         }
985
986 /*        Since B and BX are complex, the following call to SGEMM */
987 /*        is performed in two steps (real and imaginary parts). */
988
989 /*        CALL SGEMM( 'T', 'N', NR, NRHS, NR, ONE, U( NRF, 1 ), LDU, */
990 /*    $               B( NRF, 1 ), LDB, ZERO, BX( NRF, 1 ), LDBX ) */
991
992         j = nr * *nrhs << 1;
993         i__2 = *nrhs;
994         for (jcol = 1; jcol <= i__2; ++jcol) {
995             i__3 = nrf + nr - 1;
996             for (jrow = nrf; jrow <= i__3; ++jrow) {
997                 ++j;
998                 i__4 = jrow + jcol * b_dim1;
999                 rwork[j] = b[i__4].r;
1000 /* L70: */
1001             }
1002 /* L80: */
1003         }
1004         sgemm_("T", "N", &nr, nrhs, &nr, &c_b9, &u[nrf + u_dim1], ldu, &rwork[
1005                 (nr * *nrhs << 1) + 1], &nr, &c_b10, &rwork[1], &nr);
1006         j = nr * *nrhs << 1;
1007         i__2 = *nrhs;
1008         for (jcol = 1; jcol <= i__2; ++jcol) {
1009             i__3 = nrf + nr - 1;
1010             for (jrow = nrf; jrow <= i__3; ++jrow) {
1011                 ++j;
1012                 rwork[j] = r_imag(&b[jrow + jcol * b_dim1]);
1013 /* L90: */
1014             }
1015 /* L100: */
1016         }
1017         sgemm_("T", "N", &nr, nrhs, &nr, &c_b9, &u[nrf + u_dim1], ldu, &rwork[
1018                 (nr * *nrhs << 1) + 1], &nr, &c_b10, &rwork[nr * *nrhs + 1], &
1019                 nr);
1020         jreal = 0;
1021         jimag = nr * *nrhs;
1022         i__2 = *nrhs;
1023         for (jcol = 1; jcol <= i__2; ++jcol) {
1024             i__3 = nrf + nr - 1;
1025             for (jrow = nrf; jrow <= i__3; ++jrow) {
1026                 ++jreal;
1027                 ++jimag;
1028                 i__4 = jrow + jcol * bx_dim1;
1029                 i__5 = jreal;
1030                 i__6 = jimag;
1031                 q__1.r = rwork[i__5], q__1.i = rwork[i__6];
1032                 bx[i__4].r = q__1.r, bx[i__4].i = q__1.i;
1033 /* L110: */
1034             }
1035 /* L120: */
1036         }
1037
1038 /* L130: */
1039     }
1040
1041 /*     Next copy the rows of B that correspond to unchanged rows */
1042 /*     in the bidiagonal matrix to BX. */
1043
1044     i__1 = nd;
1045     for (i__ = 1; i__ <= i__1; ++i__) {
1046         ic = iwork[inode + i__ - 1];
1047         ccopy_(nrhs, &b[ic + b_dim1], ldb, &bx[ic + bx_dim1], ldbx);
1048 /* L140: */
1049     }
1050
1051 /*     Finally go through the left singular vector matrices of all */
1052 /*     the other subproblems bottom-up on the tree. */
1053
1054     j = pow_ii(&c__2, &nlvl);
1055     sqre = 0;
1056
1057     for (lvl = nlvl; lvl >= 1; --lvl) {
1058         lvl2 = (lvl << 1) - 1;
1059
1060 /*        find the first node LF and last node LL on */
1061 /*        the current level LVL */
1062
1063         if (lvl == 1) {
1064             lf = 1;
1065             ll = 1;
1066         } else {
1067             i__1 = lvl - 1;
1068             lf = pow_ii(&c__2, &i__1);
1069             ll = (lf << 1) - 1;
1070         }
1071         i__1 = ll;
1072         for (i__ = lf; i__ <= i__1; ++i__) {
1073             im1 = i__ - 1;
1074             ic = iwork[inode + im1];
1075             nl = iwork[ndiml + im1];
1076             nr = iwork[ndimr + im1];
1077             nlf = ic - nl;
1078             nrf = ic + 1;
1079             --j;
1080             clals0_(icompq, &nl, &nr, &sqre, nrhs, &bx[nlf + bx_dim1], ldbx, &
1081                     b[nlf + b_dim1], ldb, &perm[nlf + lvl * perm_dim1], &
1082                     givptr[j], &givcol[nlf + lvl2 * givcol_dim1], ldgcol, &
1083                     givnum[nlf + lvl2 * givnum_dim1], ldu, &poles[nlf + lvl2 *
1084                      poles_dim1], &difl[nlf + lvl * difl_dim1], &difr[nlf + 
1085                     lvl2 * difr_dim1], &z__[nlf + lvl * z_dim1], &k[j], &c__[
1086                     j], &s[j], &rwork[1], info);
1087 /* L150: */
1088         }
1089 /* L160: */
1090     }
1091     goto L330;
1092
1093 /*     ICOMPQ = 1: applying back the right singular vector factors. */
1094
1095 L170:
1096
1097 /*     First now go through the right singular vector matrices of all */
1098 /*     the tree nodes top-down. */
1099
1100     j = 0;
1101     i__1 = nlvl;
1102     for (lvl = 1; lvl <= i__1; ++lvl) {
1103         lvl2 = (lvl << 1) - 1;
1104
1105 /*        Find the first node LF and last node LL on */
1106 /*        the current level LVL. */
1107
1108         if (lvl == 1) {
1109             lf = 1;
1110             ll = 1;
1111         } else {
1112             i__2 = lvl - 1;
1113             lf = pow_ii(&c__2, &i__2);
1114             ll = (lf << 1) - 1;
1115         }
1116         i__2 = lf;
1117         for (i__ = ll; i__ >= i__2; --i__) {
1118             im1 = i__ - 1;
1119             ic = iwork[inode + im1];
1120             nl = iwork[ndiml + im1];
1121             nr = iwork[ndimr + im1];
1122             nlf = ic - nl;
1123             nrf = ic + 1;
1124             if (i__ == ll) {
1125                 sqre = 0;
1126             } else {
1127                 sqre = 1;
1128             }
1129             ++j;
1130             clals0_(icompq, &nl, &nr, &sqre, nrhs, &b[nlf + b_dim1], ldb, &bx[
1131                     nlf + bx_dim1], ldbx, &perm[nlf + lvl * perm_dim1], &
1132                     givptr[j], &givcol[nlf + lvl2 * givcol_dim1], ldgcol, &
1133                     givnum[nlf + lvl2 * givnum_dim1], ldu, &poles[nlf + lvl2 *
1134                      poles_dim1], &difl[nlf + lvl * difl_dim1], &difr[nlf + 
1135                     lvl2 * difr_dim1], &z__[nlf + lvl * z_dim1], &k[j], &c__[
1136                     j], &s[j], &rwork[1], info);
1137 /* L180: */
1138         }
1139 /* L190: */
1140     }
1141
1142 /*     The nodes on the bottom level of the tree were solved */
1143 /*     by SLASDQ. The corresponding right singular vector */
1144 /*     matrices are in explicit form. Apply them back. */
1145
1146     ndb1 = (nd + 1) / 2;
1147     i__1 = nd;
1148     for (i__ = ndb1; i__ <= i__1; ++i__) {
1149         i1 = i__ - 1;
1150         ic = iwork[inode + i1];
1151         nl = iwork[ndiml + i1];
1152         nr = iwork[ndimr + i1];
1153         nlp1 = nl + 1;
1154         if (i__ == nd) {
1155             nrp1 = nr;
1156         } else {
1157             nrp1 = nr + 1;
1158         }
1159         nlf = ic - nl;
1160         nrf = ic + 1;
1161
1162 /*        Since B and BX are complex, the following call to SGEMM is */
1163 /*        performed in two steps (real and imaginary parts). */
1164
1165 /*        CALL SGEMM( 'T', 'N', NLP1, NRHS, NLP1, ONE, VT( NLF, 1 ), LDU, */
1166 /*    $               B( NLF, 1 ), LDB, ZERO, BX( NLF, 1 ), LDBX ) */
1167
1168         j = nlp1 * *nrhs << 1;
1169         i__2 = *nrhs;
1170         for (jcol = 1; jcol <= i__2; ++jcol) {
1171             i__3 = nlf + nlp1 - 1;
1172             for (jrow = nlf; jrow <= i__3; ++jrow) {
1173                 ++j;
1174                 i__4 = jrow + jcol * b_dim1;
1175                 rwork[j] = b[i__4].r;
1176 /* L200: */
1177             }
1178 /* L210: */
1179         }
1180         sgemm_("T", "N", &nlp1, nrhs, &nlp1, &c_b9, &vt[nlf + vt_dim1], ldu, &
1181                 rwork[(nlp1 * *nrhs << 1) + 1], &nlp1, &c_b10, &rwork[1], &
1182                 nlp1);
1183         j = nlp1 * *nrhs << 1;
1184         i__2 = *nrhs;
1185         for (jcol = 1; jcol <= i__2; ++jcol) {
1186             i__3 = nlf + nlp1 - 1;
1187             for (jrow = nlf; jrow <= i__3; ++jrow) {
1188                 ++j;
1189                 rwork[j] = r_imag(&b[jrow + jcol * b_dim1]);
1190 /* L220: */
1191             }
1192 /* L230: */
1193         }
1194         sgemm_("T", "N", &nlp1, nrhs, &nlp1, &c_b9, &vt[nlf + vt_dim1], ldu, &
1195                 rwork[(nlp1 * *nrhs << 1) + 1], &nlp1, &c_b10, &rwork[nlp1 * *
1196                 nrhs + 1], &nlp1);
1197         jreal = 0;
1198         jimag = nlp1 * *nrhs;
1199         i__2 = *nrhs;
1200         for (jcol = 1; jcol <= i__2; ++jcol) {
1201             i__3 = nlf + nlp1 - 1;
1202             for (jrow = nlf; jrow <= i__3; ++jrow) {
1203                 ++jreal;
1204                 ++jimag;
1205                 i__4 = jrow + jcol * bx_dim1;
1206                 i__5 = jreal;
1207                 i__6 = jimag;
1208                 q__1.r = rwork[i__5], q__1.i = rwork[i__6];
1209                 bx[i__4].r = q__1.r, bx[i__4].i = q__1.i;
1210 /* L240: */
1211             }
1212 /* L250: */
1213         }
1214
1215 /*        Since B and BX are complex, the following call to SGEMM is */
1216 /*        performed in two steps (real and imaginary parts). */
1217
1218 /*        CALL SGEMM( 'T', 'N', NRP1, NRHS, NRP1, ONE, VT( NRF, 1 ), LDU, */
1219 /*    $               B( NRF, 1 ), LDB, ZERO, BX( NRF, 1 ), LDBX ) */
1220
1221         j = nrp1 * *nrhs << 1;
1222         i__2 = *nrhs;
1223         for (jcol = 1; jcol <= i__2; ++jcol) {
1224             i__3 = nrf + nrp1 - 1;
1225             for (jrow = nrf; jrow <= i__3; ++jrow) {
1226                 ++j;
1227                 i__4 = jrow + jcol * b_dim1;
1228                 rwork[j] = b[i__4].r;
1229 /* L260: */
1230             }
1231 /* L270: */
1232         }
1233         sgemm_("T", "N", &nrp1, nrhs, &nrp1, &c_b9, &vt[nrf + vt_dim1], ldu, &
1234                 rwork[(nrp1 * *nrhs << 1) + 1], &nrp1, &c_b10, &rwork[1], &
1235                 nrp1);
1236         j = nrp1 * *nrhs << 1;
1237         i__2 = *nrhs;
1238         for (jcol = 1; jcol <= i__2; ++jcol) {
1239             i__3 = nrf + nrp1 - 1;
1240             for (jrow = nrf; jrow <= i__3; ++jrow) {
1241                 ++j;
1242                 rwork[j] = r_imag(&b[jrow + jcol * b_dim1]);
1243 /* L280: */
1244             }
1245 /* L290: */
1246         }
1247         sgemm_("T", "N", &nrp1, nrhs, &nrp1, &c_b9, &vt[nrf + vt_dim1], ldu, &
1248                 rwork[(nrp1 * *nrhs << 1) + 1], &nrp1, &c_b10, &rwork[nrp1 * *
1249                 nrhs + 1], &nrp1);
1250         jreal = 0;
1251         jimag = nrp1 * *nrhs;
1252         i__2 = *nrhs;
1253         for (jcol = 1; jcol <= i__2; ++jcol) {
1254             i__3 = nrf + nrp1 - 1;
1255             for (jrow = nrf; jrow <= i__3; ++jrow) {
1256                 ++jreal;
1257                 ++jimag;
1258                 i__4 = jrow + jcol * bx_dim1;
1259                 i__5 = jreal;
1260                 i__6 = jimag;
1261                 q__1.r = rwork[i__5], q__1.i = rwork[i__6];
1262                 bx[i__4].r = q__1.r, bx[i__4].i = q__1.i;
1263 /* L300: */
1264             }
1265 /* L310: */
1266         }
1267
1268 /* L320: */
1269     }
1270
1271 L330:
1272
1273     return 0;
1274
1275 /*     End of CLALSA */
1276
1277 } /* clalsa_ */
1278