C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dlasd3.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
514 /* Table of constant values */
515
516 static integer c__1 = 1;
517 static integer c__0 = 0;
518 static doublereal c_b13 = 1.;
519 static doublereal c_b26 = 0.;
520
521 /* > \brief \b DLASD3 finds all square roots of the roots of the secular equation, as defined by the values in
522  D and Z, and then updates the singular vectors by matrix multiplication. Used by sbdsdc. */
523
524 /*  =========== DOCUMENTATION =========== */
525
526 /* Online html documentation available at */
527 /*            http://www.netlib.org/lapack/explore-html/ */
528
529 /* > \htmlonly */
530 /* > Download DLASD3 + dependencies */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlasd3.
532 f"> */
533 /* > [TGZ]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlasd3.
535 f"> */
536 /* > [ZIP]</a> */
537 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlasd3.
538 f"> */
539 /* > [TXT]</a> */
540 /* > \endhtmlonly */
541
542 /*  Definition: */
543 /*  =========== */
544
545 /*       SUBROUTINE DLASD3( NL, NR, SQRE, K, D, Q, LDQ, DSIGMA, U, LDU, U2, */
546 /*                          LDU2, VT, LDVT, VT2, LDVT2, IDXC, CTOT, Z, */
547 /*                          INFO ) */
548
549 /*       INTEGER            INFO, K, LDQ, LDU, LDU2, LDVT, LDVT2, NL, NR, */
550 /*      $                   SQRE */
551 /*       INTEGER            CTOT( * ), IDXC( * ) */
552 /*       DOUBLE PRECISION   D( * ), DSIGMA( * ), Q( LDQ, * ), U( LDU, * ), */
553 /*      $                   U2( LDU2, * ), VT( LDVT, * ), VT2( LDVT2, * ), */
554 /*      $                   Z( * ) */
555
556
557 /* > \par Purpose: */
558 /*  ============= */
559 /* > */
560 /* > \verbatim */
561 /* > */
562 /* > DLASD3 finds all the square roots of the roots of the secular */
563 /* > equation, as defined by the values in D and Z.  It makes the */
564 /* > appropriate calls to DLASD4 and then updates the singular */
565 /* > vectors by matrix multiplication. */
566 /* > */
567 /* > This code makes very mild assumptions about floating point */
568 /* > arithmetic. It will work on machines with a guard digit in */
569 /* > add/subtract, or on those binary machines without guard digits */
570 /* > which subtract like the Cray XMP, Cray YMP, Cray C 90, or Cray 2. */
571 /* > It could conceivably fail on hexadecimal or decimal machines */
572 /* > without guard digits, but we know of none. */
573 /* > */
574 /* > DLASD3 is called from DLASD1. */
575 /* > \endverbatim */
576
577 /*  Arguments: */
578 /*  ========== */
579
580 /* > \param[in] NL */
581 /* > \verbatim */
582 /* >          NL is INTEGER */
583 /* >         The row dimension of the upper block.  NL >= 1. */
584 /* > \endverbatim */
585 /* > */
586 /* > \param[in] NR */
587 /* > \verbatim */
588 /* >          NR is INTEGER */
589 /* >         The row dimension of the lower block.  NR >= 1. */
590 /* > \endverbatim */
591 /* > */
592 /* > \param[in] SQRE */
593 /* > \verbatim */
594 /* >          SQRE is INTEGER */
595 /* >         = 0: the lower block is an NR-by-NR square matrix. */
596 /* >         = 1: the lower block is an NR-by-(NR+1) rectangular matrix. */
597 /* > */
598 /* >         The bidiagonal matrix has N = NL + NR + 1 rows and */
599 /* >         M = N + SQRE >= N columns. */
600 /* > \endverbatim */
601 /* > */
602 /* > \param[in] K */
603 /* > \verbatim */
604 /* >          K is INTEGER */
605 /* >         The size of the secular equation, 1 =< K = < N. */
606 /* > \endverbatim */
607 /* > */
608 /* > \param[out] D */
609 /* > \verbatim */
610 /* >          D is DOUBLE PRECISION array, dimension(K) */
611 /* >         On exit the square roots of the roots of the secular equation, */
612 /* >         in ascending order. */
613 /* > \endverbatim */
614 /* > */
615 /* > \param[out] Q */
616 /* > \verbatim */
617 /* >          Q is DOUBLE PRECISION array, dimension (LDQ,K) */
618 /* > \endverbatim */
619 /* > */
620 /* > \param[in] LDQ */
621 /* > \verbatim */
622 /* >          LDQ is INTEGER */
623 /* >         The leading dimension of the array Q.  LDQ >= K. */
624 /* > \endverbatim */
625 /* > */
626 /* > \param[in,out] DSIGMA */
627 /* > \verbatim */
628 /* >          DSIGMA is DOUBLE PRECISION array, dimension(K) */
629 /* >         The first K elements of this array contain the old roots */
630 /* >         of the deflated updating problem.  These are the poles */
631 /* >         of the secular equation. */
632 /* > \endverbatim */
633 /* > */
634 /* > \param[out] U */
635 /* > \verbatim */
636 /* >          U is DOUBLE PRECISION array, dimension (LDU, N) */
637 /* >         The last N - K columns of this matrix contain the deflated */
638 /* >         left singular vectors. */
639 /* > \endverbatim */
640 /* > */
641 /* > \param[in] LDU */
642 /* > \verbatim */
643 /* >          LDU is INTEGER */
644 /* >         The leading dimension of the array U.  LDU >= N. */
645 /* > \endverbatim */
646 /* > */
647 /* > \param[in] U2 */
648 /* > \verbatim */
649 /* >          U2 is DOUBLE PRECISION array, dimension (LDU2, N) */
650 /* >         The first K columns of this matrix contain the non-deflated */
651 /* >         left singular vectors for the split problem. */
652 /* > \endverbatim */
653 /* > */
654 /* > \param[in] LDU2 */
655 /* > \verbatim */
656 /* >          LDU2 is INTEGER */
657 /* >         The leading dimension of the array U2.  LDU2 >= N. */
658 /* > \endverbatim */
659 /* > */
660 /* > \param[out] VT */
661 /* > \verbatim */
662 /* >          VT is DOUBLE PRECISION array, dimension (LDVT, M) */
663 /* >         The last M - K columns of VT**T contain the deflated */
664 /* >         right singular vectors. */
665 /* > \endverbatim */
666 /* > */
667 /* > \param[in] LDVT */
668 /* > \verbatim */
669 /* >          LDVT is INTEGER */
670 /* >         The leading dimension of the array VT.  LDVT >= N. */
671 /* > \endverbatim */
672 /* > */
673 /* > \param[in,out] VT2 */
674 /* > \verbatim */
675 /* >          VT2 is DOUBLE PRECISION array, dimension (LDVT2, N) */
676 /* >         The first K columns of VT2**T contain the non-deflated */
677 /* >         right singular vectors for the split problem. */
678 /* > \endverbatim */
679 /* > */
680 /* > \param[in] LDVT2 */
681 /* > \verbatim */
682 /* >          LDVT2 is INTEGER */
683 /* >         The leading dimension of the array VT2.  LDVT2 >= N. */
684 /* > \endverbatim */
685 /* > */
686 /* > \param[in] IDXC */
687 /* > \verbatim */
688 /* >          IDXC is INTEGER array, dimension ( N ) */
689 /* >         The permutation used to arrange the columns of U (and rows of */
690 /* >         VT) into three groups:  the first group contains non-zero */
691 /* >         entries only at and above (or before) NL +1; the second */
692 /* >         contains non-zero entries only at and below (or after) NL+2; */
693 /* >         and the third is dense. The first column of U and the row of */
694 /* >         VT are treated separately, however. */
695 /* > */
696 /* >         The rows of the singular vectors found by DLASD4 */
697 /* >         must be likewise permuted before the matrix multiplies can */
698 /* >         take place. */
699 /* > \endverbatim */
700 /* > */
701 /* > \param[in] CTOT */
702 /* > \verbatim */
703 /* >          CTOT is INTEGER array, dimension ( 4 ) */
704 /* >         A count of the total number of the various types of columns */
705 /* >         in U (or rows in VT), as described in IDXC. The fourth column */
706 /* >         type is any column which has been deflated. */
707 /* > \endverbatim */
708 /* > */
709 /* > \param[in,out] Z */
710 /* > \verbatim */
711 /* >          Z is DOUBLE PRECISION array, dimension (K) */
712 /* >         The first K elements of this array contain the components */
713 /* >         of the deflation-adjusted updating row vector. */
714 /* > \endverbatim */
715 /* > */
716 /* > \param[out] INFO */
717 /* > \verbatim */
718 /* >          INFO is INTEGER */
719 /* >         = 0:  successful exit. */
720 /* >         < 0:  if INFO = -i, the i-th argument had an illegal value. */
721 /* >         > 0:  if INFO = 1, a singular value did not converge */
722 /* > \endverbatim */
723
724 /*  Authors: */
725 /*  ======== */
726
727 /* > \author Univ. of Tennessee */
728 /* > \author Univ. of California Berkeley */
729 /* > \author Univ. of Colorado Denver */
730 /* > \author NAG Ltd. */
731
732 /* > \date June 2017 */
733
734 /* > \ingroup OTHERauxiliary */
735
736 /* > \par Contributors: */
737 /*  ================== */
738 /* > */
739 /* >     Ming Gu and Huan Ren, Computer Science Division, University of */
740 /* >     California at Berkeley, USA */
741 /* > */
742 /*  ===================================================================== */
743 /* Subroutine */ int dlasd3_(integer *nl, integer *nr, integer *sqre, integer 
744         *k, doublereal *d__, doublereal *q, integer *ldq, doublereal *dsigma, 
745         doublereal *u, integer *ldu, doublereal *u2, integer *ldu2, 
746         doublereal *vt, integer *ldvt, doublereal *vt2, integer *ldvt2, 
747         integer *idxc, integer *ctot, doublereal *z__, integer *info)
748 {
749     /* System generated locals */
750     integer q_dim1, q_offset, u_dim1, u_offset, u2_dim1, u2_offset, vt_dim1, 
751             vt_offset, vt2_dim1, vt2_offset, i__1, i__2;
752     doublereal d__1, d__2;
753
754     /* Local variables */
755     doublereal temp;
756     extern doublereal dnrm2_(integer *, doublereal *, integer *);
757     integer i__, j, m, n;
758     extern /* Subroutine */ int dgemm_(char *, char *, integer *, integer *, 
759             integer *, doublereal *, doublereal *, integer *, doublereal *, 
760             integer *, doublereal *, doublereal *, integer *);
761     integer ctemp;
762     extern /* Subroutine */ int dcopy_(integer *, doublereal *, integer *, 
763             doublereal *, integer *);
764     integer ktemp;
765     extern doublereal dlamc3_(doublereal *, doublereal *);
766     extern /* Subroutine */ int dlasd4_(integer *, integer *, doublereal *, 
767             doublereal *, doublereal *, doublereal *, doublereal *, 
768             doublereal *, integer *);
769     integer jc;
770     extern /* Subroutine */ int dlascl_(char *, integer *, integer *, 
771             doublereal *, doublereal *, integer *, integer *, doublereal *, 
772             integer *, integer *), dlacpy_(char *, integer *, integer 
773             *, doublereal *, integer *, doublereal *, integer *), 
774             xerbla_(char *, integer *, ftnlen);
775     doublereal rho;
776     integer nlp1, nlp2, nrp1;
777
778
779 /*  -- LAPACK auxiliary routine (version 3.7.1) -- */
780 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
781 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
782 /*     June 2017 */
783
784
785 /*  ===================================================================== */
786
787
788 /*     Test the input parameters. */
789
790     /* Parameter adjustments */
791     --d__;
792     q_dim1 = *ldq;
793     q_offset = 1 + q_dim1 * 1;
794     q -= q_offset;
795     --dsigma;
796     u_dim1 = *ldu;
797     u_offset = 1 + u_dim1 * 1;
798     u -= u_offset;
799     u2_dim1 = *ldu2;
800     u2_offset = 1 + u2_dim1 * 1;
801     u2 -= u2_offset;
802     vt_dim1 = *ldvt;
803     vt_offset = 1 + vt_dim1 * 1;
804     vt -= vt_offset;
805     vt2_dim1 = *ldvt2;
806     vt2_offset = 1 + vt2_dim1 * 1;
807     vt2 -= vt2_offset;
808     --idxc;
809     --ctot;
810     --z__;
811
812     /* Function Body */
813     *info = 0;
814
815     if (*nl < 1) {
816         *info = -1;
817     } else if (*nr < 1) {
818         *info = -2;
819     } else if (*sqre != 1 && *sqre != 0) {
820         *info = -3;
821     }
822
823     n = *nl + *nr + 1;
824     m = n + *sqre;
825     nlp1 = *nl + 1;
826     nlp2 = *nl + 2;
827
828     if (*k < 1 || *k > n) {
829         *info = -4;
830     } else if (*ldq < *k) {
831         *info = -7;
832     } else if (*ldu < n) {
833         *info = -10;
834     } else if (*ldu2 < n) {
835         *info = -12;
836     } else if (*ldvt < m) {
837         *info = -14;
838     } else if (*ldvt2 < m) {
839         *info = -16;
840     }
841     if (*info != 0) {
842         i__1 = -(*info);
843         xerbla_("DLASD3", &i__1, (ftnlen)6);
844         return 0;
845     }
846
847 /*     Quick return if possible */
848
849     if (*k == 1) {
850         d__[1] = abs(z__[1]);
851         dcopy_(&m, &vt2[vt2_dim1 + 1], ldvt2, &vt[vt_dim1 + 1], ldvt);
852         if (z__[1] > 0.) {
853             dcopy_(&n, &u2[u2_dim1 + 1], &c__1, &u[u_dim1 + 1], &c__1);
854         } else {
855             i__1 = n;
856             for (i__ = 1; i__ <= i__1; ++i__) {
857                 u[i__ + u_dim1] = -u2[i__ + u2_dim1];
858 /* L10: */
859             }
860         }
861         return 0;
862     }
863
864 /*     Modify values DSIGMA(i) to make sure all DSIGMA(i)-DSIGMA(j) can */
865 /*     be computed with high relative accuracy (barring over/underflow). */
866 /*     This is a problem on machines without a guard digit in */
867 /*     add/subtract (Cray XMP, Cray YMP, Cray C 90 and Cray 2). */
868 /*     The following code replaces DSIGMA(I) by 2*DSIGMA(I)-DSIGMA(I), */
869 /*     which on any of these machines zeros out the bottommost */
870 /*     bit of DSIGMA(I) if it is 1; this makes the subsequent */
871 /*     subtractions DSIGMA(I)-DSIGMA(J) unproblematic when cancellation */
872 /*     occurs. On binary machines with a guard digit (almost all */
873 /*     machines) it does not change DSIGMA(I) at all. On hexadecimal */
874 /*     and decimal machines with a guard digit, it slightly */
875 /*     changes the bottommost bits of DSIGMA(I). It does not account */
876 /*     for hexadecimal or decimal machines without guard digits */
877 /*     (we know of none). We use a subroutine call to compute */
878 /*     2*DSIGMA(I) to prevent optimizing compilers from eliminating */
879 /*     this code. */
880
881     i__1 = *k;
882     for (i__ = 1; i__ <= i__1; ++i__) {
883         dsigma[i__] = dlamc3_(&dsigma[i__], &dsigma[i__]) - dsigma[i__];
884 /* L20: */
885     }
886
887 /*     Keep a copy of Z. */
888
889     dcopy_(k, &z__[1], &c__1, &q[q_offset], &c__1);
890
891 /*     Normalize Z. */
892
893     rho = dnrm2_(k, &z__[1], &c__1);
894     dlascl_("G", &c__0, &c__0, &rho, &c_b13, k, &c__1, &z__[1], k, info);
895     rho *= rho;
896
897 /*     Find the new singular values. */
898
899     i__1 = *k;
900     for (j = 1; j <= i__1; ++j) {
901         dlasd4_(k, &j, &dsigma[1], &z__[1], &u[j * u_dim1 + 1], &rho, &d__[j],
902                  &vt[j * vt_dim1 + 1], info);
903
904 /*        If the zero finder fails, report the convergence failure. */
905
906         if (*info != 0) {
907             return 0;
908         }
909 /* L30: */
910     }
911
912 /*     Compute updated Z. */
913
914     i__1 = *k;
915     for (i__ = 1; i__ <= i__1; ++i__) {
916         z__[i__] = u[i__ + *k * u_dim1] * vt[i__ + *k * vt_dim1];
917         i__2 = i__ - 1;
918         for (j = 1; j <= i__2; ++j) {
919             z__[i__] *= u[i__ + j * u_dim1] * vt[i__ + j * vt_dim1] / (dsigma[
920                     i__] - dsigma[j]) / (dsigma[i__] + dsigma[j]);
921 /* L40: */
922         }
923         i__2 = *k - 1;
924         for (j = i__; j <= i__2; ++j) {
925             z__[i__] *= u[i__ + j * u_dim1] * vt[i__ + j * vt_dim1] / (dsigma[
926                     i__] - dsigma[j + 1]) / (dsigma[i__] + dsigma[j + 1]);
927 /* L50: */
928         }
929         d__2 = sqrt((d__1 = z__[i__], abs(d__1)));
930         z__[i__] = d_sign(&d__2, &q[i__ + q_dim1]);
931 /* L60: */
932     }
933
934 /*     Compute left singular vectors of the modified diagonal matrix, */
935 /*     and store related information for the right singular vectors. */
936
937     i__1 = *k;
938     for (i__ = 1; i__ <= i__1; ++i__) {
939         vt[i__ * vt_dim1 + 1] = z__[1] / u[i__ * u_dim1 + 1] / vt[i__ * 
940                 vt_dim1 + 1];
941         u[i__ * u_dim1 + 1] = -1.;
942         i__2 = *k;
943         for (j = 2; j <= i__2; ++j) {
944             vt[j + i__ * vt_dim1] = z__[j] / u[j + i__ * u_dim1] / vt[j + i__ 
945                     * vt_dim1];
946             u[j + i__ * u_dim1] = dsigma[j] * vt[j + i__ * vt_dim1];
947 /* L70: */
948         }
949         temp = dnrm2_(k, &u[i__ * u_dim1 + 1], &c__1);
950         q[i__ * q_dim1 + 1] = u[i__ * u_dim1 + 1] / temp;
951         i__2 = *k;
952         for (j = 2; j <= i__2; ++j) {
953             jc = idxc[j];
954             q[j + i__ * q_dim1] = u[jc + i__ * u_dim1] / temp;
955 /* L80: */
956         }
957 /* L90: */
958     }
959
960 /*     Update the left singular vector matrix. */
961
962     if (*k == 2) {
963         dgemm_("N", "N", &n, k, k, &c_b13, &u2[u2_offset], ldu2, &q[q_offset],
964                  ldq, &c_b26, &u[u_offset], ldu);
965         goto L100;
966     }
967     if (ctot[1] > 0) {
968         dgemm_("N", "N", nl, k, &ctot[1], &c_b13, &u2[(u2_dim1 << 1) + 1], 
969                 ldu2, &q[q_dim1 + 2], ldq, &c_b26, &u[u_dim1 + 1], ldu);
970         if (ctot[3] > 0) {
971             ktemp = ctot[1] + 2 + ctot[2];
972             dgemm_("N", "N", nl, k, &ctot[3], &c_b13, &u2[ktemp * u2_dim1 + 1]
973                     , ldu2, &q[ktemp + q_dim1], ldq, &c_b13, &u[u_dim1 + 1], 
974                     ldu);
975         }
976     } else if (ctot[3] > 0) {
977         ktemp = ctot[1] + 2 + ctot[2];
978         dgemm_("N", "N", nl, k, &ctot[3], &c_b13, &u2[ktemp * u2_dim1 + 1], 
979                 ldu2, &q[ktemp + q_dim1], ldq, &c_b26, &u[u_dim1 + 1], ldu);
980     } else {
981         dlacpy_("F", nl, k, &u2[u2_offset], ldu2, &u[u_offset], ldu);
982     }
983     dcopy_(k, &q[q_dim1 + 1], ldq, &u[nlp1 + u_dim1], ldu);
984     ktemp = ctot[1] + 2;
985     ctemp = ctot[2] + ctot[3];
986     dgemm_("N", "N", nr, k, &ctemp, &c_b13, &u2[nlp2 + ktemp * u2_dim1], ldu2,
987              &q[ktemp + q_dim1], ldq, &c_b26, &u[nlp2 + u_dim1], ldu);
988
989 /*     Generate the right singular vectors. */
990
991 L100:
992     i__1 = *k;
993     for (i__ = 1; i__ <= i__1; ++i__) {
994         temp = dnrm2_(k, &vt[i__ * vt_dim1 + 1], &c__1);
995         q[i__ + q_dim1] = vt[i__ * vt_dim1 + 1] / temp;
996         i__2 = *k;
997         for (j = 2; j <= i__2; ++j) {
998             jc = idxc[j];
999             q[i__ + j * q_dim1] = vt[jc + i__ * vt_dim1] / temp;
1000 /* L110: */
1001         }
1002 /* L120: */
1003     }
1004
1005 /*     Update the right singular vector matrix. */
1006
1007     if (*k == 2) {
1008         dgemm_("N", "N", k, &m, k, &c_b13, &q[q_offset], ldq, &vt2[vt2_offset]
1009                 , ldvt2, &c_b26, &vt[vt_offset], ldvt);
1010         return 0;
1011     }
1012     ktemp = ctot[1] + 1;
1013     dgemm_("N", "N", k, &nlp1, &ktemp, &c_b13, &q[q_dim1 + 1], ldq, &vt2[
1014             vt2_dim1 + 1], ldvt2, &c_b26, &vt[vt_dim1 + 1], ldvt);
1015     ktemp = ctot[1] + 2 + ctot[2];
1016     if (ktemp <= *ldvt2) {
1017         dgemm_("N", "N", k, &nlp1, &ctot[3], &c_b13, &q[ktemp * q_dim1 + 1], 
1018                 ldq, &vt2[ktemp + vt2_dim1], ldvt2, &c_b13, &vt[vt_dim1 + 1], 
1019                 ldvt);
1020     }
1021
1022     ktemp = ctot[1] + 1;
1023     nrp1 = *nr + *sqre;
1024     if (ktemp > 1) {
1025         i__1 = *k;
1026         for (i__ = 1; i__ <= i__1; ++i__) {
1027             q[i__ + ktemp * q_dim1] = q[i__ + q_dim1];
1028 /* L130: */
1029         }
1030         i__1 = m;
1031         for (i__ = nlp2; i__ <= i__1; ++i__) {
1032             vt2[ktemp + i__ * vt2_dim1] = vt2[i__ * vt2_dim1 + 1];
1033 /* L140: */
1034         }
1035     }
1036     ctemp = ctot[2] + 1 + ctot[3];
1037     dgemm_("N", "N", k, &nrp1, &ctemp, &c_b13, &q[ktemp * q_dim1 + 1], ldq, &
1038             vt2[ktemp + nlp2 * vt2_dim1], ldvt2, &c_b26, &vt[nlp2 * vt_dim1 + 
1039             1], ldvt);
1040
1041     return 0;
1042
1043 /*     End of DLASD3 */
1044
1045 } /* dlasd3_ */
1046