C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dbbcsd.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 doublereal c_b10 = -.125;
516 static doublereal c_b35 = -1.;
517 static integer c__1 = 1;
518
519 /* > \brief \b DBBCSD */
520
521 /*  =========== DOCUMENTATION =========== */
522
523 /* Online html documentation available at */
524 /*            http://www.netlib.org/lapack/explore-html/ */
525
526 /* > \htmlonly */
527 /* > Download DBBCSD + dependencies */
528 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dbbcsd.
529 f"> */
530 /* > [TGZ]</a> */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dbbcsd.
532 f"> */
533 /* > [ZIP]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dbbcsd.
535 f"> */
536 /* > [TXT]</a> */
537 /* > \endhtmlonly */
538
539 /*  Definition: */
540 /*  =========== */
541
542 /*       SUBROUTINE DBBCSD( JOBU1, JOBU2, JOBV1T, JOBV2T, TRANS, M, P, Q, */
543 /*                          THETA, PHI, U1, LDU1, U2, LDU2, V1T, LDV1T, */
544 /*                          V2T, LDV2T, B11D, B11E, B12D, B12E, B21D, B21E, */
545 /*                          B22D, B22E, WORK, LWORK, INFO ) */
546
547 /*       CHARACTER          JOBU1, JOBU2, JOBV1T, JOBV2T, TRANS */
548 /*       INTEGER            INFO, LDU1, LDU2, LDV1T, LDV2T, LWORK, M, P, Q */
549 /*       DOUBLE PRECISION   B11D( * ), B11E( * ), B12D( * ), B12E( * ), */
550 /*      $                   B21D( * ), B21E( * ), B22D( * ), B22E( * ), */
551 /*      $                   PHI( * ), THETA( * ), WORK( * ) */
552 /*       DOUBLE PRECISION   U1( LDU1, * ), U2( LDU2, * ), V1T( LDV1T, * ), */
553 /*      $                   V2T( LDV2T, * ) */
554
555
556 /* > \par Purpose: */
557 /*  ============= */
558 /* > */
559 /* > \verbatim */
560 /* > */
561 /* > DBBCSD computes the CS decomposition of an orthogonal matrix in */
562 /* > bidiagonal-block form, */
563 /* > */
564 /* > */
565 /* >     [ B11 | B12 0  0 ] */
566 /* >     [  0  |  0 -I  0 ] */
567 /* > X = [----------------] */
568 /* >     [ B21 | B22 0  0 ] */
569 /* >     [  0  |  0  0  I ] */
570 /* > */
571 /* >                               [  C | -S  0  0 ] */
572 /* >                   [ U1 |    ] [  0 |  0 -I  0 ] [ V1 |    ]**T */
573 /* >                 = [---------] [---------------] [---------]   . */
574 /* >                   [    | U2 ] [  S |  C  0  0 ] [    | V2 ] */
575 /* >                               [  0 |  0  0  I ] */
576 /* > */
577 /* > X is M-by-M, its top-left block is P-by-Q, and Q must be no larger */
578 /* > than P, M-P, or M-Q. (If Q is not the smallest index, then X must be */
579 /* > transposed and/or permuted. This can be done in constant time using */
580 /* > the TRANS and SIGNS options. See DORCSD for details.) */
581 /* > */
582 /* > The bidiagonal matrices B11, B12, B21, and B22 are represented */
583 /* > implicitly by angles THETA(1:Q) and PHI(1:Q-1). */
584 /* > */
585 /* > The orthogonal matrices U1, U2, V1T, and V2T are input/output. */
586 /* > The input matrices are pre- or post-multiplied by the appropriate */
587 /* > singular vector matrices. */
588 /* > \endverbatim */
589
590 /*  Arguments: */
591 /*  ========== */
592
593 /* > \param[in] JOBU1 */
594 /* > \verbatim */
595 /* >          JOBU1 is CHARACTER */
596 /* >          = 'Y':      U1 is updated; */
597 /* >          otherwise:  U1 is not updated. */
598 /* > \endverbatim */
599 /* > */
600 /* > \param[in] JOBU2 */
601 /* > \verbatim */
602 /* >          JOBU2 is CHARACTER */
603 /* >          = 'Y':      U2 is updated; */
604 /* >          otherwise:  U2 is not updated. */
605 /* > \endverbatim */
606 /* > */
607 /* > \param[in] JOBV1T */
608 /* > \verbatim */
609 /* >          JOBV1T is CHARACTER */
610 /* >          = 'Y':      V1T is updated; */
611 /* >          otherwise:  V1T is not updated. */
612 /* > \endverbatim */
613 /* > */
614 /* > \param[in] JOBV2T */
615 /* > \verbatim */
616 /* >          JOBV2T is CHARACTER */
617 /* >          = 'Y':      V2T is updated; */
618 /* >          otherwise:  V2T is not updated. */
619 /* > \endverbatim */
620 /* > */
621 /* > \param[in] TRANS */
622 /* > \verbatim */
623 /* >          TRANS is CHARACTER */
624 /* >          = 'T':      X, U1, U2, V1T, and V2T are stored in row-major */
625 /* >                      order; */
626 /* >          otherwise:  X, U1, U2, V1T, and V2T are stored in column- */
627 /* >                      major order. */
628 /* > \endverbatim */
629 /* > */
630 /* > \param[in] M */
631 /* > \verbatim */
632 /* >          M is INTEGER */
633 /* >          The number of rows and columns in X, the orthogonal matrix in */
634 /* >          bidiagonal-block form. */
635 /* > \endverbatim */
636 /* > */
637 /* > \param[in] P */
638 /* > \verbatim */
639 /* >          P is INTEGER */
640 /* >          The number of rows in the top-left block of X. 0 <= P <= M. */
641 /* > \endverbatim */
642 /* > */
643 /* > \param[in] Q */
644 /* > \verbatim */
645 /* >          Q is INTEGER */
646 /* >          The number of columns in the top-left block of X. */
647 /* >          0 <= Q <= MIN(P,M-P,M-Q). */
648 /* > \endverbatim */
649 /* > */
650 /* > \param[in,out] THETA */
651 /* > \verbatim */
652 /* >          THETA is DOUBLE PRECISION array, dimension (Q) */
653 /* >          On entry, the angles THETA(1),...,THETA(Q) that, along with */
654 /* >          PHI(1), ...,PHI(Q-1), define the matrix in bidiagonal-block */
655 /* >          form. On exit, the angles whose cosines and sines define the */
656 /* >          diagonal blocks in the CS decomposition. */
657 /* > \endverbatim */
658 /* > */
659 /* > \param[in,out] PHI */
660 /* > \verbatim */
661 /* >          PHI is DOUBLE PRECISION array, dimension (Q-1) */
662 /* >          The angles PHI(1),...,PHI(Q-1) that, along with THETA(1),..., */
663 /* >          THETA(Q), define the matrix in bidiagonal-block form. */
664 /* > \endverbatim */
665 /* > */
666 /* > \param[in,out] U1 */
667 /* > \verbatim */
668 /* >          U1 is DOUBLE PRECISION array, dimension (LDU1,P) */
669 /* >          On entry, a P-by-P matrix. On exit, U1 is postmultiplied */
670 /* >          by the left singular vector matrix common to [ B11 ; 0 ] and */
671 /* >          [ B12 0 0 ; 0 -I 0 0 ]. */
672 /* > \endverbatim */
673 /* > */
674 /* > \param[in] LDU1 */
675 /* > \verbatim */
676 /* >          LDU1 is INTEGER */
677 /* >          The leading dimension of the array U1, LDU1 >= MAX(1,P). */
678 /* > \endverbatim */
679 /* > */
680 /* > \param[in,out] U2 */
681 /* > \verbatim */
682 /* >          U2 is DOUBLE PRECISION array, dimension (LDU2,M-P) */
683 /* >          On entry, an (M-P)-by-(M-P) matrix. On exit, U2 is */
684 /* >          postmultiplied by the left singular vector matrix common to */
685 /* >          [ B21 ; 0 ] and [ B22 0 0 ; 0 0 I ]. */
686 /* > \endverbatim */
687 /* > */
688 /* > \param[in] LDU2 */
689 /* > \verbatim */
690 /* >          LDU2 is INTEGER */
691 /* >          The leading dimension of the array U2, LDU2 >= MAX(1,M-P). */
692 /* > \endverbatim */
693 /* > */
694 /* > \param[in,out] V1T */
695 /* > \verbatim */
696 /* >          V1T is DOUBLE PRECISION array, dimension (LDV1T,Q) */
697 /* >          On entry, a Q-by-Q matrix. On exit, V1T is premultiplied */
698 /* >          by the transpose of the right singular vector */
699 /* >          matrix common to [ B11 ; 0 ] and [ B21 ; 0 ]. */
700 /* > \endverbatim */
701 /* > */
702 /* > \param[in] LDV1T */
703 /* > \verbatim */
704 /* >          LDV1T is INTEGER */
705 /* >          The leading dimension of the array V1T, LDV1T >= MAX(1,Q). */
706 /* > \endverbatim */
707 /* > */
708 /* > \param[in,out] V2T */
709 /* > \verbatim */
710 /* >          V2T is DOUBLE PRECISION array, dimension (LDV2T,M-Q) */
711 /* >          On entry, an (M-Q)-by-(M-Q) matrix. On exit, V2T is */
712 /* >          premultiplied by the transpose of the right */
713 /* >          singular vector matrix common to [ B12 0 0 ; 0 -I 0 ] and */
714 /* >          [ B22 0 0 ; 0 0 I ]. */
715 /* > \endverbatim */
716 /* > */
717 /* > \param[in] LDV2T */
718 /* > \verbatim */
719 /* >          LDV2T is INTEGER */
720 /* >          The leading dimension of the array V2T, LDV2T >= MAX(1,M-Q). */
721 /* > \endverbatim */
722 /* > */
723 /* > \param[out] B11D */
724 /* > \verbatim */
725 /* >          B11D is DOUBLE PRECISION array, dimension (Q) */
726 /* >          When DBBCSD converges, B11D contains the cosines of THETA(1), */
727 /* >          ..., THETA(Q). If DBBCSD fails to converge, then B11D */
728 /* >          contains the diagonal of the partially reduced top-left */
729 /* >          block. */
730 /* > \endverbatim */
731 /* > */
732 /* > \param[out] B11E */
733 /* > \verbatim */
734 /* >          B11E is DOUBLE PRECISION array, dimension (Q-1) */
735 /* >          When DBBCSD converges, B11E contains zeros. If DBBCSD fails */
736 /* >          to converge, then B11E contains the superdiagonal of the */
737 /* >          partially reduced top-left block. */
738 /* > \endverbatim */
739 /* > */
740 /* > \param[out] B12D */
741 /* > \verbatim */
742 /* >          B12D is DOUBLE PRECISION array, dimension (Q) */
743 /* >          When DBBCSD converges, B12D contains the negative sines of */
744 /* >          THETA(1), ..., THETA(Q). If DBBCSD fails to converge, then */
745 /* >          B12D contains the diagonal of the partially reduced top-right */
746 /* >          block. */
747 /* > \endverbatim */
748 /* > */
749 /* > \param[out] B12E */
750 /* > \verbatim */
751 /* >          B12E is DOUBLE PRECISION array, dimension (Q-1) */
752 /* >          When DBBCSD converges, B12E contains zeros. If DBBCSD fails */
753 /* >          to converge, then B12E contains the subdiagonal of the */
754 /* >          partially reduced top-right block. */
755 /* > \endverbatim */
756 /* > */
757 /* > \param[out] B21D */
758 /* > \verbatim */
759 /* >          B21D is DOUBLE PRECISION  array, dimension (Q) */
760 /* >          When DBBCSD converges, B21D contains the negative sines of */
761 /* >          THETA(1), ..., THETA(Q). If DBBCSD fails to converge, then */
762 /* >          B21D contains the diagonal of the partially reduced bottom-left */
763 /* >          block. */
764 /* > \endverbatim */
765 /* > */
766 /* > \param[out] B21E */
767 /* > \verbatim */
768 /* >          B21E is DOUBLE PRECISION  array, dimension (Q-1) */
769 /* >          When DBBCSD converges, B21E contains zeros. If DBBCSD fails */
770 /* >          to converge, then B21E contains the subdiagonal of the */
771 /* >          partially reduced bottom-left block. */
772 /* > \endverbatim */
773 /* > */
774 /* > \param[out] B22D */
775 /* > \verbatim */
776 /* >          B22D is DOUBLE PRECISION  array, dimension (Q) */
777 /* >          When DBBCSD converges, B22D contains the negative sines of */
778 /* >          THETA(1), ..., THETA(Q). If DBBCSD fails to converge, then */
779 /* >          B22D contains the diagonal of the partially reduced bottom-right */
780 /* >          block. */
781 /* > \endverbatim */
782 /* > */
783 /* > \param[out] B22E */
784 /* > \verbatim */
785 /* >          B22E is DOUBLE PRECISION  array, dimension (Q-1) */
786 /* >          When DBBCSD converges, B22E contains zeros. If DBBCSD fails */
787 /* >          to converge, then B22E contains the subdiagonal of the */
788 /* >          partially reduced bottom-right block. */
789 /* > \endverbatim */
790 /* > */
791 /* > \param[out] WORK */
792 /* > \verbatim */
793 /* >          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
794 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
795 /* > \endverbatim */
796 /* > */
797 /* > \param[in] LWORK */
798 /* > \verbatim */
799 /* >          LWORK is INTEGER */
800 /* >          The dimension of the array WORK. LWORK >= MAX(1,8*Q). */
801 /* > */
802 /* >          If LWORK = -1, then a workspace query is assumed; the */
803 /* >          routine only calculates the optimal size of the WORK array, */
804 /* >          returns this value as the first entry of the work array, and */
805 /* >          no error message related to LWORK is issued by XERBLA. */
806 /* > \endverbatim */
807 /* > */
808 /* > \param[out] INFO */
809 /* > \verbatim */
810 /* >          INFO is INTEGER */
811 /* >          = 0:  successful exit. */
812 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
813 /* >          > 0:  if DBBCSD did not converge, INFO specifies the number */
814 /* >                of nonzero entries in PHI, and B11D, B11E, etc., */
815 /* >                contain the partially reduced matrix. */
816 /* > \endverbatim */
817
818 /* > \par Internal Parameters: */
819 /*  ========================= */
820 /* > */
821 /* > \verbatim */
822 /* >  TOLMUL  DOUBLE PRECISION, default = MAX(10,MIN(100,EPS**(-1/8))) */
823 /* >          TOLMUL controls the convergence criterion of the QR loop. */
824 /* >          Angles THETA(i), PHI(i) are rounded to 0 or PI/2 when they */
825 /* >          are within TOLMUL*EPS of either bound. */
826 /* > \endverbatim */
827
828 /* > \par References: */
829 /*  ================ */
830 /* > */
831 /* >  [1] Brian D. Sutton. Computing the complete CS decomposition. Numer. */
832 /* >      Algorithms, 50(1):33-65, 2009. */
833
834 /*  Authors: */
835 /*  ======== */
836
837 /* > \author Univ. of Tennessee */
838 /* > \author Univ. of California Berkeley */
839 /* > \author Univ. of Colorado Denver */
840 /* > \author NAG Ltd. */
841
842 /* > \date June 2016 */
843
844 /* > \ingroup doubleOTHERcomputational */
845
846 /*  ===================================================================== */
847 /* Subroutine */ int dbbcsd_(char *jobu1, char *jobu2, char *jobv1t, char *
848         jobv2t, char *trans, integer *m, integer *p, integer *q, doublereal *
849         theta, doublereal *phi, doublereal *u1, integer *ldu1, doublereal *u2,
850          integer *ldu2, doublereal *v1t, integer *ldv1t, doublereal *v2t, 
851         integer *ldv2t, doublereal *b11d, doublereal *b11e, doublereal *b12d, 
852         doublereal *b12e, doublereal *b21d, doublereal *b21e, doublereal *
853         b22d, doublereal *b22e, doublereal *work, integer *lwork, integer *
854         info)
855 {
856     /* System generated locals */
857     integer u1_dim1, u1_offset, u2_dim1, u2_offset, v1t_dim1, v1t_offset, 
858             v2t_dim1, v2t_offset, i__1, i__2;
859     doublereal d__1, d__2, d__3, d__4;
860
861     /* Local variables */
862     integer imin, mini, imax, iter;
863     doublereal unfl, temp;
864     logical colmajor;
865     doublereal thetamin, thetamax;
866     logical restart11, restart12, restart21, restart22;
867     extern /* Subroutine */ int dlas2_(doublereal *, doublereal *, doublereal 
868             *, doublereal *, doublereal *);
869     integer lworkmin, iu1cs, iu2cs, iu1sn, iu2sn, lworkopt, i__, j;
870     doublereal r__;
871     extern /* Subroutine */ int dscal_(integer *, doublereal *, doublereal *, 
872             integer *);
873     extern logical lsame_(char *, char *);
874     extern /* Subroutine */ int dlasr_(char *, char *, char *, integer *, 
875             integer *, doublereal *, doublereal *, doublereal *, integer *), dswap_(integer *, doublereal *, integer *
876             , doublereal *, integer *);
877     integer maxit;
878     doublereal dummy, x1, x2, y1, y2;
879     integer iv1tcs, iv2tcs;
880     logical wantu1, wantu2;
881     integer iv1tsn, iv2tsn;
882     extern doublereal dlamch_(char *);
883     doublereal mu, nu, sigma11, sigma21;
884     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
885     doublereal thresh, tolmul;
886     extern /* Subroutine */ int mecago_();
887     logical lquery;
888     doublereal b11bulge;
889     logical wantv1t, wantv2t;
890     doublereal b12bulge, b21bulge, b22bulge, eps, tol;
891     extern /* Subroutine */ int dlartgp_(doublereal *, doublereal *, 
892             doublereal *, doublereal *, doublereal *), dlartgs_(doublereal *, 
893             doublereal *, doublereal *, doublereal *, doublereal *);
894
895
896 /*  -- LAPACK computational routine (version 3.7.1) -- */
897 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
898 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
899 /*     June 2016 */
900
901
902 /*  =================================================================== */
903
904
905
906 /*     Test input arguments */
907
908     /* Parameter adjustments */
909     --theta;
910     --phi;
911     u1_dim1 = *ldu1;
912     u1_offset = 1 + u1_dim1 * 1;
913     u1 -= u1_offset;
914     u2_dim1 = *ldu2;
915     u2_offset = 1 + u2_dim1 * 1;
916     u2 -= u2_offset;
917     v1t_dim1 = *ldv1t;
918     v1t_offset = 1 + v1t_dim1 * 1;
919     v1t -= v1t_offset;
920     v2t_dim1 = *ldv2t;
921     v2t_offset = 1 + v2t_dim1 * 1;
922     v2t -= v2t_offset;
923     --b11d;
924     --b11e;
925     --b12d;
926     --b12e;
927     --b21d;
928     --b21e;
929     --b22d;
930     --b22e;
931     --work;
932
933     /* Function Body */
934     *info = 0;
935     lquery = *lwork == -1;
936     wantu1 = lsame_(jobu1, "Y");
937     wantu2 = lsame_(jobu2, "Y");
938     wantv1t = lsame_(jobv1t, "Y");
939     wantv2t = lsame_(jobv2t, "Y");
940     colmajor = ! lsame_(trans, "T");
941
942     if (*m < 0) {
943         *info = -6;
944     } else if (*p < 0 || *p > *m) {
945         *info = -7;
946     } else if (*q < 0 || *q > *m) {
947         *info = -8;
948     } else if (*q > *p || *q > *m - *p || *q > *m - *q) {
949         *info = -8;
950     } else if (wantu1 && *ldu1 < *p) {
951         *info = -12;
952     } else if (wantu2 && *ldu2 < *m - *p) {
953         *info = -14;
954     } else if (wantv1t && *ldv1t < *q) {
955         *info = -16;
956     } else if (wantv2t && *ldv2t < *m - *q) {
957         *info = -18;
958     }
959
960 /*     Quick return if Q = 0 */
961
962     if (*info == 0 && *q == 0) {
963         lworkmin = 1;
964         work[1] = (doublereal) lworkmin;
965         return 0;
966     }
967
968 /*     Compute workspace */
969
970     if (*info == 0) {
971         iu1cs = 1;
972         iu1sn = iu1cs + *q;
973         iu2cs = iu1sn + *q;
974         iu2sn = iu2cs + *q;
975         iv1tcs = iu2sn + *q;
976         iv1tsn = iv1tcs + *q;
977         iv2tcs = iv1tsn + *q;
978         iv2tsn = iv2tcs + *q;
979         lworkopt = iv2tsn + *q - 1;
980         lworkmin = lworkopt;
981         work[1] = (doublereal) lworkopt;
982         if (*lwork < lworkmin && ! lquery) {
983             *info = -28;
984         }
985     }
986
987     if (*info != 0) {
988         i__1 = -(*info);
989         xerbla_("DBBCSD", &i__1, (ftnlen)6);
990         return 0;
991     } else if (lquery) {
992         return 0;
993     }
994
995 /*     Get machine constants */
996
997     eps = dlamch_("Epsilon");
998     unfl = dlamch_("Safe minimum");
999 /* Computing MAX */
1000 /* Computing MIN */
1001     d__3 = 100., d__4 = pow_dd(&eps, &c_b10);
1002     d__1 = 10., d__2 = f2cmin(d__3,d__4);
1003     tolmul = f2cmax(d__1,d__2);
1004     tol = tolmul * eps;
1005 /* Computing MAX */
1006     d__1 = tol, d__2 = *q * 6 * *q * unfl;
1007     thresh = f2cmax(d__1,d__2);
1008
1009 /*     Test for negligible sines or cosines */
1010
1011     i__1 = *q;
1012     for (i__ = 1; i__ <= i__1; ++i__) {
1013         if (theta[i__] < thresh) {
1014             theta[i__] = 0.;
1015         } else if (theta[i__] > 1.57079632679489662 - thresh) {
1016             theta[i__] = 1.57079632679489662;
1017         }
1018     }
1019     i__1 = *q - 1;
1020     for (i__ = 1; i__ <= i__1; ++i__) {
1021         if (phi[i__] < thresh) {
1022             phi[i__] = 0.;
1023         } else if (phi[i__] > 1.57079632679489662 - thresh) {
1024             phi[i__] = 1.57079632679489662;
1025         }
1026     }
1027
1028 /*     Initial deflation */
1029
1030     imax = *q;
1031     while(imax > 1) {
1032         if (phi[imax - 1] != 0.) {
1033             myexit_();
1034         }
1035         --imax;
1036     }
1037     imin = imax - 1;
1038     if (imin > 1) {
1039         while(phi[imin - 1] != 0.) {
1040             --imin;
1041             if (imin <= 1) {
1042                 myexit_();
1043             }
1044         }
1045     }
1046
1047 /*     Initialize iteration counter */
1048
1049     maxit = *q * 6 * *q;
1050     iter = 0;
1051
1052 /*     Begin main iteration loop */
1053
1054     while(imax > 1) {
1055
1056 /*        Compute the matrix entries */
1057
1058         b11d[imin] = cos(theta[imin]);
1059         b21d[imin] = -sin(theta[imin]);
1060         i__1 = imax - 1;
1061         for (i__ = imin; i__ <= i__1; ++i__) {
1062             b11e[i__] = -sin(theta[i__]) * sin(phi[i__]);
1063             b11d[i__ + 1] = cos(theta[i__ + 1]) * cos(phi[i__]);
1064             b12d[i__] = sin(theta[i__]) * cos(phi[i__]);
1065             b12e[i__] = cos(theta[i__ + 1]) * sin(phi[i__]);
1066             b21e[i__] = -cos(theta[i__]) * sin(phi[i__]);
1067             b21d[i__ + 1] = -sin(theta[i__ + 1]) * cos(phi[i__]);
1068             b22d[i__] = cos(theta[i__]) * cos(phi[i__]);
1069             b22e[i__] = -sin(theta[i__ + 1]) * sin(phi[i__]);
1070         }
1071         b12d[imax] = sin(theta[imax]);
1072         b22d[imax] = cos(theta[imax]);
1073
1074 /*        Abort if not converging; otherwise, increment ITER */
1075
1076         if (iter > maxit) {
1077             *info = 0;
1078             i__1 = *q;
1079             for (i__ = 1; i__ <= i__1; ++i__) {
1080                 if (phi[i__] != 0.) {
1081                     ++(*info);
1082                 }
1083             }
1084             return 0;
1085         }
1086
1087         iter = iter + imax - imin;
1088
1089 /*        Compute shifts */
1090
1091         thetamax = theta[imin];
1092         thetamin = theta[imin];
1093         i__1 = imax;
1094         for (i__ = imin + 1; i__ <= i__1; ++i__) {
1095             if (theta[i__] > thetamax) {
1096                 thetamax = theta[i__];
1097             }
1098             if (theta[i__] < thetamin) {
1099                 thetamin = theta[i__];
1100             }
1101         }
1102
1103         if (thetamax > 1.57079632679489662 - thresh) {
1104
1105 /*           Zero on diagonals of B11 and B22; induce deflation with a */
1106 /*           zero shift */
1107
1108             mu = 0.;
1109             nu = 1.;
1110
1111         } else if (thetamin < thresh) {
1112
1113 /*           Zero on diagonals of B12 and B22; induce deflation with a */
1114 /*           zero shift */
1115
1116             mu = 1.;
1117             nu = 0.;
1118
1119         } else {
1120
1121 /*           Compute shifts for B11 and B21 and use the lesser */
1122
1123             dlas2_(&b11d[imax - 1], &b11e[imax - 1], &b11d[imax], &sigma11, &
1124                     dummy);
1125             dlas2_(&b21d[imax - 1], &b21e[imax - 1], &b21d[imax], &sigma21, &
1126                     dummy);
1127
1128             if (sigma11 <= sigma21) {
1129                 mu = sigma11;
1130 /* Computing 2nd power */
1131                 d__1 = mu;
1132                 nu = sqrt(1. - d__1 * d__1);
1133                 if (mu < thresh) {
1134                     mu = 0.;
1135                     nu = 1.;
1136                 }
1137             } else {
1138                 nu = sigma21;
1139 /* Computing 2nd power */
1140                 d__1 = nu;
1141                 mu = sqrt(1.f - d__1 * d__1);
1142                 if (nu < thresh) {
1143                     mu = 1.;
1144                     nu = 0.;
1145                 }
1146             }
1147         }
1148
1149 /*        Rotate to produce bulges in B11 and B21 */
1150
1151         if (mu <= nu) {
1152             dlartgs_(&b11d[imin], &b11e[imin], &mu, &work[iv1tcs + imin - 1], 
1153                     &work[iv1tsn + imin - 1]);
1154         } else {
1155             dlartgs_(&b21d[imin], &b21e[imin], &nu, &work[iv1tcs + imin - 1], 
1156                     &work[iv1tsn + imin - 1]);
1157         }
1158
1159         temp = work[iv1tcs + imin - 1] * b11d[imin] + work[iv1tsn + imin - 1] 
1160                 * b11e[imin];
1161         b11e[imin] = work[iv1tcs + imin - 1] * b11e[imin] - work[iv1tsn + 
1162                 imin - 1] * b11d[imin];
1163         b11d[imin] = temp;
1164         b11bulge = work[iv1tsn + imin - 1] * b11d[imin + 1];
1165         b11d[imin + 1] = work[iv1tcs + imin - 1] * b11d[imin + 1];
1166         temp = work[iv1tcs + imin - 1] * b21d[imin] + work[iv1tsn + imin - 1] 
1167                 * b21e[imin];
1168         b21e[imin] = work[iv1tcs + imin - 1] * b21e[imin] - work[iv1tsn + 
1169                 imin - 1] * b21d[imin];
1170         b21d[imin] = temp;
1171         b21bulge = work[iv1tsn + imin - 1] * b21d[imin + 1];
1172         b21d[imin + 1] = work[iv1tcs + imin - 1] * b21d[imin + 1];
1173
1174 /*        Compute THETA(IMIN) */
1175
1176 /* Computing 2nd power */
1177         d__1 = b21d[imin];
1178 /* Computing 2nd power */
1179         d__2 = b21bulge;
1180 /* Computing 2nd power */
1181         d__3 = b11d[imin];
1182 /* Computing 2nd power */
1183         d__4 = b11bulge;
1184         theta[imin] = atan2(sqrt(d__1 * d__1 + d__2 * d__2), sqrt(d__3 * d__3 
1185                 + d__4 * d__4));
1186
1187 /*        Chase the bulges in B11(IMIN+1,IMIN) and B21(IMIN+1,IMIN) */
1188
1189 /* Computing 2nd power */
1190         d__1 = b11d[imin];
1191 /* Computing 2nd power */
1192         d__2 = b11bulge;
1193 /* Computing 2nd power */
1194         d__3 = thresh;
1195         if (d__1 * d__1 + d__2 * d__2 > d__3 * d__3) {
1196             dlartgp_(&b11bulge, &b11d[imin], &work[iu1sn + imin - 1], &work[
1197                     iu1cs + imin - 1], &r__);
1198         } else if (mu <= nu) {
1199             dlartgs_(&b11e[imin], &b11d[imin + 1], &mu, &work[iu1cs + imin - 
1200                     1], &work[iu1sn + imin - 1]);
1201         } else {
1202             dlartgs_(&b12d[imin], &b12e[imin], &nu, &work[iu1cs + imin - 1], &
1203                     work[iu1sn + imin - 1]);
1204         }
1205 /* Computing 2nd power */
1206         d__1 = b21d[imin];
1207 /* Computing 2nd power */
1208         d__2 = b21bulge;
1209 /* Computing 2nd power */
1210         d__3 = thresh;
1211         if (d__1 * d__1 + d__2 * d__2 > d__3 * d__3) {
1212             dlartgp_(&b21bulge, &b21d[imin], &work[iu2sn + imin - 1], &work[
1213                     iu2cs + imin - 1], &r__);
1214         } else if (nu < mu) {
1215             dlartgs_(&b21e[imin], &b21d[imin + 1], &nu, &work[iu2cs + imin - 
1216                     1], &work[iu2sn + imin - 1]);
1217         } else {
1218             dlartgs_(&b22d[imin], &b22e[imin], &mu, &work[iu2cs + imin - 1], &
1219                     work[iu2sn + imin - 1]);
1220         }
1221         work[iu2cs + imin - 1] = -work[iu2cs + imin - 1];
1222         work[iu2sn + imin - 1] = -work[iu2sn + imin - 1];
1223
1224         temp = work[iu1cs + imin - 1] * b11e[imin] + work[iu1sn + imin - 1] * 
1225                 b11d[imin + 1];
1226         b11d[imin + 1] = work[iu1cs + imin - 1] * b11d[imin + 1] - work[iu1sn 
1227                 + imin - 1] * b11e[imin];
1228         b11e[imin] = temp;
1229         if (imax > imin + 1) {
1230             b11bulge = work[iu1sn + imin - 1] * b11e[imin + 1];
1231             b11e[imin + 1] = work[iu1cs + imin - 1] * b11e[imin + 1];
1232         }
1233         temp = work[iu1cs + imin - 1] * b12d[imin] + work[iu1sn + imin - 1] * 
1234                 b12e[imin];
1235         b12e[imin] = work[iu1cs + imin - 1] * b12e[imin] - work[iu1sn + imin 
1236                 - 1] * b12d[imin];
1237         b12d[imin] = temp;
1238         b12bulge = work[iu1sn + imin - 1] * b12d[imin + 1];
1239         b12d[imin + 1] = work[iu1cs + imin - 1] * b12d[imin + 1];
1240         temp = work[iu2cs + imin - 1] * b21e[imin] + work[iu2sn + imin - 1] * 
1241                 b21d[imin + 1];
1242         b21d[imin + 1] = work[iu2cs + imin - 1] * b21d[imin + 1] - work[iu2sn 
1243                 + imin - 1] * b21e[imin];
1244         b21e[imin] = temp;
1245         if (imax > imin + 1) {
1246             b21bulge = work[iu2sn + imin - 1] * b21e[imin + 1];
1247             b21e[imin + 1] = work[iu2cs + imin - 1] * b21e[imin + 1];
1248         }
1249         temp = work[iu2cs + imin - 1] * b22d[imin] + work[iu2sn + imin - 1] * 
1250                 b22e[imin];
1251         b22e[imin] = work[iu2cs + imin - 1] * b22e[imin] - work[iu2sn + imin 
1252                 - 1] * b22d[imin];
1253         b22d[imin] = temp;
1254         b22bulge = work[iu2sn + imin - 1] * b22d[imin + 1];
1255         b22d[imin + 1] = work[iu2cs + imin - 1] * b22d[imin + 1];
1256
1257 /*        Inner loop: chase bulges from B11(IMIN,IMIN+2), */
1258 /*        B12(IMIN,IMIN+1), B21(IMIN,IMIN+2), and B22(IMIN,IMIN+1) to */
1259 /*        bottom-right */
1260
1261         i__1 = imax - 1;
1262         for (i__ = imin + 1; i__ <= i__1; ++i__) {
1263
1264 /*           Compute PHI(I-1) */
1265
1266             x1 = sin(theta[i__ - 1]) * b11e[i__ - 1] + cos(theta[i__ - 1]) * 
1267                     b21e[i__ - 1];
1268             x2 = sin(theta[i__ - 1]) * b11bulge + cos(theta[i__ - 1]) * 
1269                     b21bulge;
1270             y1 = sin(theta[i__ - 1]) * b12d[i__ - 1] + cos(theta[i__ - 1]) * 
1271                     b22d[i__ - 1];
1272             y2 = sin(theta[i__ - 1]) * b12bulge + cos(theta[i__ - 1]) * 
1273                     b22bulge;
1274
1275 /* Computing 2nd power */
1276             d__1 = x1;
1277 /* Computing 2nd power */
1278             d__2 = x2;
1279 /* Computing 2nd power */
1280             d__3 = y1;
1281 /* Computing 2nd power */
1282             d__4 = y2;
1283             phi[i__ - 1] = atan2(sqrt(d__1 * d__1 + d__2 * d__2), sqrt(d__3 * 
1284                     d__3 + d__4 * d__4));
1285
1286 /*           Determine if there are bulges to chase or if a new direct */
1287 /*           summand has been reached */
1288
1289 /* Computing 2nd power */
1290             d__1 = b11e[i__ - 1];
1291 /* Computing 2nd power */
1292             d__2 = b11bulge;
1293 /* Computing 2nd power */
1294             d__3 = thresh;
1295             restart11 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1296 /* Computing 2nd power */
1297             d__1 = b21e[i__ - 1];
1298 /* Computing 2nd power */
1299             d__2 = b21bulge;
1300 /* Computing 2nd power */
1301             d__3 = thresh;
1302             restart21 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1303 /* Computing 2nd power */
1304             d__1 = b12d[i__ - 1];
1305 /* Computing 2nd power */
1306             d__2 = b12bulge;
1307 /* Computing 2nd power */
1308             d__3 = thresh;
1309             restart12 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1310 /* Computing 2nd power */
1311             d__1 = b22d[i__ - 1];
1312 /* Computing 2nd power */
1313             d__2 = b22bulge;
1314 /* Computing 2nd power */
1315             d__3 = thresh;
1316             restart22 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1317
1318 /*           If possible, chase bulges from B11(I-1,I+1), B12(I-1,I), */
1319 /*           B21(I-1,I+1), and B22(I-1,I). If necessary, restart bulge- */
1320 /*           chasing by applying the original shift again. */
1321
1322             if (! restart11 && ! restart21) {
1323                 dlartgp_(&x2, &x1, &work[iv1tsn + i__ - 1], &work[iv1tcs + 
1324                         i__ - 1], &r__);
1325             } else if (! restart11 && restart21) {
1326                 dlartgp_(&b11bulge, &b11e[i__ - 1], &work[iv1tsn + i__ - 1], &
1327                         work[iv1tcs + i__ - 1], &r__);
1328             } else if (restart11 && ! restart21) {
1329                 dlartgp_(&b21bulge, &b21e[i__ - 1], &work[iv1tsn + i__ - 1], &
1330                         work[iv1tcs + i__ - 1], &r__);
1331             } else if (mu <= nu) {
1332                 dlartgs_(&b11d[i__], &b11e[i__], &mu, &work[iv1tcs + i__ - 1],
1333                          &work[iv1tsn + i__ - 1]);
1334             } else {
1335                 dlartgs_(&b21d[i__], &b21e[i__], &nu, &work[iv1tcs + i__ - 1],
1336                          &work[iv1tsn + i__ - 1]);
1337             }
1338             work[iv1tcs + i__ - 1] = -work[iv1tcs + i__ - 1];
1339             work[iv1tsn + i__ - 1] = -work[iv1tsn + i__ - 1];
1340             if (! restart12 && ! restart22) {
1341                 dlartgp_(&y2, &y1, &work[iv2tsn + i__ - 2], &work[iv2tcs + 
1342                         i__ - 2], &r__);
1343             } else if (! restart12 && restart22) {
1344                 dlartgp_(&b12bulge, &b12d[i__ - 1], &work[iv2tsn + i__ - 2], &
1345                         work[iv2tcs + i__ - 2], &r__);
1346             } else if (restart12 && ! restart22) {
1347                 dlartgp_(&b22bulge, &b22d[i__ - 1], &work[iv2tsn + i__ - 2], &
1348                         work[iv2tcs + i__ - 2], &r__);
1349             } else if (nu < mu) {
1350                 dlartgs_(&b12e[i__ - 1], &b12d[i__], &nu, &work[iv2tcs + i__ 
1351                         - 2], &work[iv2tsn + i__ - 2]);
1352             } else {
1353                 dlartgs_(&b22e[i__ - 1], &b22d[i__], &mu, &work[iv2tcs + i__ 
1354                         - 2], &work[iv2tsn + i__ - 2]);
1355             }
1356
1357             temp = work[iv1tcs + i__ - 1] * b11d[i__] + work[iv1tsn + i__ - 1]
1358                      * b11e[i__];
1359             b11e[i__] = work[iv1tcs + i__ - 1] * b11e[i__] - work[iv1tsn + 
1360                     i__ - 1] * b11d[i__];
1361             b11d[i__] = temp;
1362             b11bulge = work[iv1tsn + i__ - 1] * b11d[i__ + 1];
1363             b11d[i__ + 1] = work[iv1tcs + i__ - 1] * b11d[i__ + 1];
1364             temp = work[iv1tcs + i__ - 1] * b21d[i__] + work[iv1tsn + i__ - 1]
1365                      * b21e[i__];
1366             b21e[i__] = work[iv1tcs + i__ - 1] * b21e[i__] - work[iv1tsn + 
1367                     i__ - 1] * b21d[i__];
1368             b21d[i__] = temp;
1369             b21bulge = work[iv1tsn + i__ - 1] * b21d[i__ + 1];
1370             b21d[i__ + 1] = work[iv1tcs + i__ - 1] * b21d[i__ + 1];
1371             temp = work[iv2tcs + i__ - 2] * b12e[i__ - 1] + work[iv2tsn + i__ 
1372                     - 2] * b12d[i__];
1373             b12d[i__] = work[iv2tcs + i__ - 2] * b12d[i__] - work[iv2tsn + 
1374                     i__ - 2] * b12e[i__ - 1];
1375             b12e[i__ - 1] = temp;
1376             b12bulge = work[iv2tsn + i__ - 2] * b12e[i__];
1377             b12e[i__] = work[iv2tcs + i__ - 2] * b12e[i__];
1378             temp = work[iv2tcs + i__ - 2] * b22e[i__ - 1] + work[iv2tsn + i__ 
1379                     - 2] * b22d[i__];
1380             b22d[i__] = work[iv2tcs + i__ - 2] * b22d[i__] - work[iv2tsn + 
1381                     i__ - 2] * b22e[i__ - 1];
1382             b22e[i__ - 1] = temp;
1383             b22bulge = work[iv2tsn + i__ - 2] * b22e[i__];
1384             b22e[i__] = work[iv2tcs + i__ - 2] * b22e[i__];
1385
1386 /*           Compute THETA(I) */
1387
1388             x1 = cos(phi[i__ - 1]) * b11d[i__] + sin(phi[i__ - 1]) * b12e[i__ 
1389                     - 1];
1390             x2 = cos(phi[i__ - 1]) * b11bulge + sin(phi[i__ - 1]) * b12bulge;
1391             y1 = cos(phi[i__ - 1]) * b21d[i__] + sin(phi[i__ - 1]) * b22e[i__ 
1392                     - 1];
1393             y2 = cos(phi[i__ - 1]) * b21bulge + sin(phi[i__ - 1]) * b22bulge;
1394
1395 /* Computing 2nd power */
1396             d__1 = y1;
1397 /* Computing 2nd power */
1398             d__2 = y2;
1399 /* Computing 2nd power */
1400             d__3 = x1;
1401 /* Computing 2nd power */
1402             d__4 = x2;
1403             theta[i__] = atan2(sqrt(d__1 * d__1 + d__2 * d__2), sqrt(d__3 * 
1404                     d__3 + d__4 * d__4));
1405
1406 /*           Determine if there are bulges to chase or if a new direct */
1407 /*           summand has been reached */
1408
1409 /* Computing 2nd power */
1410             d__1 = b11d[i__];
1411 /* Computing 2nd power */
1412             d__2 = b11bulge;
1413 /* Computing 2nd power */
1414             d__3 = thresh;
1415             restart11 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1416 /* Computing 2nd power */
1417             d__1 = b12e[i__ - 1];
1418 /* Computing 2nd power */
1419             d__2 = b12bulge;
1420 /* Computing 2nd power */
1421             d__3 = thresh;
1422             restart12 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1423 /* Computing 2nd power */
1424             d__1 = b21d[i__];
1425 /* Computing 2nd power */
1426             d__2 = b21bulge;
1427 /* Computing 2nd power */
1428             d__3 = thresh;
1429             restart21 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1430 /* Computing 2nd power */
1431             d__1 = b22e[i__ - 1];
1432 /* Computing 2nd power */
1433             d__2 = b22bulge;
1434 /* Computing 2nd power */
1435             d__3 = thresh;
1436             restart22 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1437
1438 /*           If possible, chase bulges from B11(I+1,I), B12(I+1,I-1), */
1439 /*           B21(I+1,I), and B22(I+1,I-1). If necessary, restart bulge- */
1440 /*           chasing by applying the original shift again. */
1441
1442             if (! restart11 && ! restart12) {
1443                 dlartgp_(&x2, &x1, &work[iu1sn + i__ - 1], &work[iu1cs + i__ 
1444                         - 1], &r__);
1445             } else if (! restart11 && restart12) {
1446                 dlartgp_(&b11bulge, &b11d[i__], &work[iu1sn + i__ - 1], &work[
1447                         iu1cs + i__ - 1], &r__);
1448             } else if (restart11 && ! restart12) {
1449                 dlartgp_(&b12bulge, &b12e[i__ - 1], &work[iu1sn + i__ - 1], &
1450                         work[iu1cs + i__ - 1], &r__);
1451             } else if (mu <= nu) {
1452                 dlartgs_(&b11e[i__], &b11d[i__ + 1], &mu, &work[iu1cs + i__ - 
1453                         1], &work[iu1sn + i__ - 1]);
1454             } else {
1455                 dlartgs_(&b12d[i__], &b12e[i__], &nu, &work[iu1cs + i__ - 1], 
1456                         &work[iu1sn + i__ - 1]);
1457             }
1458             if (! restart21 && ! restart22) {
1459                 dlartgp_(&y2, &y1, &work[iu2sn + i__ - 1], &work[iu2cs + i__ 
1460                         - 1], &r__);
1461             } else if (! restart21 && restart22) {
1462                 dlartgp_(&b21bulge, &b21d[i__], &work[iu2sn + i__ - 1], &work[
1463                         iu2cs + i__ - 1], &r__);
1464             } else if (restart21 && ! restart22) {
1465                 dlartgp_(&b22bulge, &b22e[i__ - 1], &work[iu2sn + i__ - 1], &
1466                         work[iu2cs + i__ - 1], &r__);
1467             } else if (nu < mu) {
1468                 dlartgs_(&b21e[i__], &b21e[i__ + 1], &nu, &work[iu2cs + i__ - 
1469                         1], &work[iu2sn + i__ - 1]);
1470             } else {
1471                 dlartgs_(&b22d[i__], &b22e[i__], &mu, &work[iu2cs + i__ - 1], 
1472                         &work[iu2sn + i__ - 1]);
1473             }
1474             work[iu2cs + i__ - 1] = -work[iu2cs + i__ - 1];
1475             work[iu2sn + i__ - 1] = -work[iu2sn + i__ - 1];
1476
1477             temp = work[iu1cs + i__ - 1] * b11e[i__] + work[iu1sn + i__ - 1] *
1478                      b11d[i__ + 1];
1479             b11d[i__ + 1] = work[iu1cs + i__ - 1] * b11d[i__ + 1] - work[
1480                     iu1sn + i__ - 1] * b11e[i__];
1481             b11e[i__] = temp;
1482             if (i__ < imax - 1) {
1483                 b11bulge = work[iu1sn + i__ - 1] * b11e[i__ + 1];
1484                 b11e[i__ + 1] = work[iu1cs + i__ - 1] * b11e[i__ + 1];
1485             }
1486             temp = work[iu2cs + i__ - 1] * b21e[i__] + work[iu2sn + i__ - 1] *
1487                      b21d[i__ + 1];
1488             b21d[i__ + 1] = work[iu2cs + i__ - 1] * b21d[i__ + 1] - work[
1489                     iu2sn + i__ - 1] * b21e[i__];
1490             b21e[i__] = temp;
1491             if (i__ < imax - 1) {
1492                 b21bulge = work[iu2sn + i__ - 1] * b21e[i__ + 1];
1493                 b21e[i__ + 1] = work[iu2cs + i__ - 1] * b21e[i__ + 1];
1494             }
1495             temp = work[iu1cs + i__ - 1] * b12d[i__] + work[iu1sn + i__ - 1] *
1496                      b12e[i__];
1497             b12e[i__] = work[iu1cs + i__ - 1] * b12e[i__] - work[iu1sn + i__ 
1498                     - 1] * b12d[i__];
1499             b12d[i__] = temp;
1500             b12bulge = work[iu1sn + i__ - 1] * b12d[i__ + 1];
1501             b12d[i__ + 1] = work[iu1cs + i__ - 1] * b12d[i__ + 1];
1502             temp = work[iu2cs + i__ - 1] * b22d[i__] + work[iu2sn + i__ - 1] *
1503                      b22e[i__];
1504             b22e[i__] = work[iu2cs + i__ - 1] * b22e[i__] - work[iu2sn + i__ 
1505                     - 1] * b22d[i__];
1506             b22d[i__] = temp;
1507             b22bulge = work[iu2sn + i__ - 1] * b22d[i__ + 1];
1508             b22d[i__ + 1] = work[iu2cs + i__ - 1] * b22d[i__ + 1];
1509
1510         }
1511
1512 /*        Compute PHI(IMAX-1) */
1513
1514         x1 = sin(theta[imax - 1]) * b11e[imax - 1] + cos(theta[imax - 1]) * 
1515                 b21e[imax - 1];
1516         y1 = sin(theta[imax - 1]) * b12d[imax - 1] + cos(theta[imax - 1]) * 
1517                 b22d[imax - 1];
1518         y2 = sin(theta[imax - 1]) * b12bulge + cos(theta[imax - 1]) * 
1519                 b22bulge;
1520
1521 /* Computing 2nd power */
1522         d__1 = y1;
1523 /* Computing 2nd power */
1524         d__2 = y2;
1525         phi[imax - 1] = atan2((abs(x1)), sqrt(d__1 * d__1 + d__2 * d__2));
1526
1527 /*        Chase bulges from B12(IMAX-1,IMAX) and B22(IMAX-1,IMAX) */
1528
1529 /* Computing 2nd power */
1530         d__1 = b12d[imax - 1];
1531 /* Computing 2nd power */
1532         d__2 = b12bulge;
1533 /* Computing 2nd power */
1534         d__3 = thresh;
1535         restart12 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1536 /* Computing 2nd power */
1537         d__1 = b22d[imax - 1];
1538 /* Computing 2nd power */
1539         d__2 = b22bulge;
1540 /* Computing 2nd power */
1541         d__3 = thresh;
1542         restart22 = d__1 * d__1 + d__2 * d__2 <= d__3 * d__3;
1543
1544         if (! restart12 && ! restart22) {
1545             dlartgp_(&y2, &y1, &work[iv2tsn + imax - 2], &work[iv2tcs + imax 
1546                     - 2], &r__);
1547         } else if (! restart12 && restart22) {
1548             dlartgp_(&b12bulge, &b12d[imax - 1], &work[iv2tsn + imax - 2], &
1549                     work[iv2tcs + imax - 2], &r__);
1550         } else if (restart12 && ! restart22) {
1551             dlartgp_(&b22bulge, &b22d[imax - 1], &work[iv2tsn + imax - 2], &
1552                     work[iv2tcs + imax - 2], &r__);
1553         } else if (nu < mu) {
1554             dlartgs_(&b12e[imax - 1], &b12d[imax], &nu, &work[iv2tcs + imax - 
1555                     2], &work[iv2tsn + imax - 2]);
1556         } else {
1557             dlartgs_(&b22e[imax - 1], &b22d[imax], &mu, &work[iv2tcs + imax - 
1558                     2], &work[iv2tsn + imax - 2]);
1559         }
1560
1561         temp = work[iv2tcs + imax - 2] * b12e[imax - 1] + work[iv2tsn + imax 
1562                 - 2] * b12d[imax];
1563         b12d[imax] = work[iv2tcs + imax - 2] * b12d[imax] - work[iv2tsn + 
1564                 imax - 2] * b12e[imax - 1];
1565         b12e[imax - 1] = temp;
1566         temp = work[iv2tcs + imax - 2] * b22e[imax - 1] + work[iv2tsn + imax 
1567                 - 2] * b22d[imax];
1568         b22d[imax] = work[iv2tcs + imax - 2] * b22d[imax] - work[iv2tsn + 
1569                 imax - 2] * b22e[imax - 1];
1570         b22e[imax - 1] = temp;
1571
1572 /*        Update singular vectors */
1573
1574         if (wantu1) {
1575             if (colmajor) {
1576                 i__1 = imax - imin + 1;
1577                 dlasr_("R", "V", "F", p, &i__1, &work[iu1cs + imin - 1], &
1578                         work[iu1sn + imin - 1], &u1[imin * u1_dim1 + 1], ldu1);
1579             } else {
1580                 i__1 = imax - imin + 1;
1581                 dlasr_("L", "V", "F", &i__1, p, &work[iu1cs + imin - 1], &
1582                         work[iu1sn + imin - 1], &u1[imin + u1_dim1], ldu1);
1583             }
1584         }
1585         if (wantu2) {
1586             if (colmajor) {
1587                 i__1 = *m - *p;
1588                 i__2 = imax - imin + 1;
1589                 dlasr_("R", "V", "F", &i__1, &i__2, &work[iu2cs + imin - 1], &
1590                         work[iu2sn + imin - 1], &u2[imin * u2_dim1 + 1], ldu2);
1591             } else {
1592                 i__1 = imax - imin + 1;
1593                 i__2 = *m - *p;
1594                 dlasr_("L", "V", "F", &i__1, &i__2, &work[iu2cs + imin - 1], &
1595                         work[iu2sn + imin - 1], &u2[imin + u2_dim1], ldu2);
1596             }
1597         }
1598         if (wantv1t) {
1599             if (colmajor) {
1600                 i__1 = imax - imin + 1;
1601                 dlasr_("L", "V", "F", &i__1, q, &work[iv1tcs + imin - 1], &
1602                         work[iv1tsn + imin - 1], &v1t[imin + v1t_dim1], ldv1t);
1603             } else {
1604                 i__1 = imax - imin + 1;
1605                 dlasr_("R", "V", "F", q, &i__1, &work[iv1tcs + imin - 1], &
1606                         work[iv1tsn + imin - 1], &v1t[imin * v1t_dim1 + 1], 
1607                         ldv1t);
1608             }
1609         }
1610         if (wantv2t) {
1611             if (colmajor) {
1612                 i__1 = imax - imin + 1;
1613                 i__2 = *m - *q;
1614                 dlasr_("L", "V", "F", &i__1, &i__2, &work[iv2tcs + imin - 1], 
1615                         &work[iv2tsn + imin - 1], &v2t[imin + v2t_dim1], 
1616                         ldv2t);
1617             } else {
1618                 i__1 = *m - *q;
1619                 i__2 = imax - imin + 1;
1620                 dlasr_("R", "V", "F", &i__1, &i__2, &work[iv2tcs + imin - 1], 
1621                         &work[iv2tsn + imin - 1], &v2t[imin * v2t_dim1 + 1], 
1622                         ldv2t);
1623             }
1624         }
1625
1626 /*        Fix signs on B11(IMAX-1,IMAX) and B21(IMAX-1,IMAX) */
1627
1628         if (b11e[imax - 1] + b21e[imax - 1] > 0.) {
1629             b11d[imax] = -b11d[imax];
1630             b21d[imax] = -b21d[imax];
1631             if (wantv1t) {
1632                 if (colmajor) {
1633                     dscal_(q, &c_b35, &v1t[imax + v1t_dim1], ldv1t);
1634                 } else {
1635                     dscal_(q, &c_b35, &v1t[imax * v1t_dim1 + 1], &c__1);
1636                 }
1637             }
1638         }
1639
1640 /*        Compute THETA(IMAX) */
1641
1642         x1 = cos(phi[imax - 1]) * b11d[imax] + sin(phi[imax - 1]) * b12e[imax 
1643                 - 1];
1644         y1 = cos(phi[imax - 1]) * b21d[imax] + sin(phi[imax - 1]) * b22e[imax 
1645                 - 1];
1646
1647         theta[imax] = atan2((abs(y1)), (abs(x1)));
1648
1649 /*        Fix signs on B11(IMAX,IMAX), B12(IMAX,IMAX-1), B21(IMAX,IMAX), */
1650 /*        and B22(IMAX,IMAX-1) */
1651
1652         if (b11d[imax] + b12e[imax - 1] < 0.) {
1653             b12d[imax] = -b12d[imax];
1654             if (wantu1) {
1655                 if (colmajor) {
1656                     dscal_(p, &c_b35, &u1[imax * u1_dim1 + 1], &c__1);
1657                 } else {
1658                     dscal_(p, &c_b35, &u1[imax + u1_dim1], ldu1);
1659                 }
1660             }
1661         }
1662         if (b21d[imax] + b22e[imax - 1] > 0.) {
1663             b22d[imax] = -b22d[imax];
1664             if (wantu2) {
1665                 if (colmajor) {
1666                     i__1 = *m - *p;
1667                     dscal_(&i__1, &c_b35, &u2[imax * u2_dim1 + 1], &c__1);
1668                 } else {
1669                     i__1 = *m - *p;
1670                     dscal_(&i__1, &c_b35, &u2[imax + u2_dim1], ldu2);
1671                 }
1672             }
1673         }
1674
1675 /*        Fix signs on B12(IMAX,IMAX) and B22(IMAX,IMAX) */
1676
1677         if (b12d[imax] + b22d[imax] < 0.) {
1678             if (wantv2t) {
1679                 if (colmajor) {
1680                     i__1 = *m - *q;
1681                     dscal_(&i__1, &c_b35, &v2t[imax + v2t_dim1], ldv2t);
1682                 } else {
1683                     i__1 = *m - *q;
1684                     dscal_(&i__1, &c_b35, &v2t[imax * v2t_dim1 + 1], &c__1);
1685                 }
1686             }
1687         }
1688
1689 /*        Test for negligible sines or cosines */
1690
1691         i__1 = imax;
1692         for (i__ = imin; i__ <= i__1; ++i__) {
1693             if (theta[i__] < thresh) {
1694                 theta[i__] = 0.;
1695             } else if (theta[i__] > 1.57079632679489662 - thresh) {
1696                 theta[i__] = 1.57079632679489662;
1697             }
1698         }
1699         i__1 = imax - 1;
1700         for (i__ = imin; i__ <= i__1; ++i__) {
1701             if (phi[i__] < thresh) {
1702                 phi[i__] = 0.;
1703             } else if (phi[i__] > 1.57079632679489662 - thresh) {
1704                 phi[i__] = 1.57079632679489662;
1705             }
1706         }
1707
1708 /*        Deflate */
1709
1710         if (imax > 1) {
1711             while(phi[imax - 1] == 0.) {
1712                 --imax;
1713                 if (imax <= 1) {
1714                     myexit_();
1715                 }
1716             }
1717         }
1718         if (imin > imax - 1) {
1719             imin = imax - 1;
1720         }
1721         if (imin > 1) {
1722             while(phi[imin - 1] != 0.) {
1723                 --imin;
1724                 if (imin <= 1) {
1725                     myexit_();
1726                 }
1727             }
1728         }
1729
1730 /*        Repeat main iteration loop */
1731
1732     }
1733
1734 /*     Postprocessing: order THETA from least to greatest */
1735
1736     i__1 = *q;
1737     for (i__ = 1; i__ <= i__1; ++i__) {
1738
1739         mini = i__;
1740         thetamin = theta[i__];
1741         i__2 = *q;
1742         for (j = i__ + 1; j <= i__2; ++j) {
1743             if (theta[j] < thetamin) {
1744                 mini = j;
1745                 thetamin = theta[j];
1746             }
1747         }
1748
1749         if (mini != i__) {
1750             theta[mini] = theta[i__];
1751             theta[i__] = thetamin;
1752             if (colmajor) {
1753                 if (wantu1) {
1754                     dswap_(p, &u1[i__ * u1_dim1 + 1], &c__1, &u1[mini * 
1755                             u1_dim1 + 1], &c__1);
1756                 }
1757                 if (wantu2) {
1758                     i__2 = *m - *p;
1759                     dswap_(&i__2, &u2[i__ * u2_dim1 + 1], &c__1, &u2[mini * 
1760                             u2_dim1 + 1], &c__1);
1761                 }
1762                 if (wantv1t) {
1763                     dswap_(q, &v1t[i__ + v1t_dim1], ldv1t, &v1t[mini + 
1764                             v1t_dim1], ldv1t);
1765                 }
1766                 if (wantv2t) {
1767                     i__2 = *m - *q;
1768                     dswap_(&i__2, &v2t[i__ + v2t_dim1], ldv2t, &v2t[mini + 
1769                             v2t_dim1], ldv2t);
1770                 }
1771             } else {
1772                 if (wantu1) {
1773                     dswap_(p, &u1[i__ + u1_dim1], ldu1, &u1[mini + u1_dim1], 
1774                             ldu1);
1775                 }
1776                 if (wantu2) {
1777                     i__2 = *m - *p;
1778                     dswap_(&i__2, &u2[i__ + u2_dim1], ldu2, &u2[mini + 
1779                             u2_dim1], ldu2);
1780                 }
1781                 if (wantv1t) {
1782                     dswap_(q, &v1t[i__ * v1t_dim1 + 1], &c__1, &v1t[mini * 
1783                             v1t_dim1 + 1], &c__1);
1784                 }
1785                 if (wantv2t) {
1786                     i__2 = *m - *q;
1787                     dswap_(&i__2, &v2t[i__ * v2t_dim1 + 1], &c__1, &v2t[mini *
1788                              v2t_dim1 + 1], &c__1);
1789                 }
1790             }
1791         }
1792
1793     }
1794
1795     return 0;
1796
1797 /*     End of DBBCSD */
1798
1799 } /* dbbcsd_ */
1800