C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dgglse.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static integer c__1 = 1;
516 static integer c_n1 = -1;
517 static doublereal c_b31 = -1.;
518 static doublereal c_b33 = 1.;
519
520 /* > \brief <b> DGGLSE solves overdetermined or underdetermined systems for OTHER matrices</b> */
521
522 /*  =========== DOCUMENTATION =========== */
523
524 /* Online html documentation available at */
525 /*            http://www.netlib.org/lapack/explore-html/ */
526
527 /* > \htmlonly */
528 /* > Download DGGLSE + dependencies */
529 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgglse.
530 f"> */
531 /* > [TGZ]</a> */
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgglse.
533 f"> */
534 /* > [ZIP]</a> */
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgglse.
536 f"> */
537 /* > [TXT]</a> */
538 /* > \endhtmlonly */
539
540 /*  Definition: */
541 /*  =========== */
542
543 /*       SUBROUTINE DGGLSE( M, N, P, A, LDA, B, LDB, C, D, X, WORK, LWORK, */
544 /*                          INFO ) */
545
546 /*       INTEGER            INFO, LDA, LDB, LWORK, M, N, P */
547 /*       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), C( * ), D( * ), */
548 /*      $                   WORK( * ), X( * ) */
549
550
551 /* > \par Purpose: */
552 /*  ============= */
553 /* > */
554 /* > \verbatim */
555 /* > */
556 /* > DGGLSE solves the linear equality-constrained least squares (LSE) */
557 /* > problem: */
558 /* > */
559 /* >         minimize || c - A*x ||_2   subject to   B*x = d */
560 /* > */
561 /* > where A is an M-by-N matrix, B is a P-by-N matrix, c is a given */
562 /* > M-vector, and d is a given P-vector. It is assumed that */
563 /* > P <= N <= M+P, and */
564 /* > */
565 /* >          rank(B) = P and  rank( (A) ) = N. */
566 /* >                               ( (B) ) */
567 /* > */
568 /* > These conditions ensure that the LSE problem has a unique solution, */
569 /* > which is obtained using a generalized RQ factorization of the */
570 /* > matrices (B, A) given by */
571 /* > */
572 /* >    B = (0 R)*Q,   A = Z*T*Q. */
573 /* > \endverbatim */
574
575 /*  Arguments: */
576 /*  ========== */
577
578 /* > \param[in] M */
579 /* > \verbatim */
580 /* >          M is INTEGER */
581 /* >          The number of rows of the matrix A.  M >= 0. */
582 /* > \endverbatim */
583 /* > */
584 /* > \param[in] N */
585 /* > \verbatim */
586 /* >          N is INTEGER */
587 /* >          The number of columns of the matrices A and B. N >= 0. */
588 /* > \endverbatim */
589 /* > */
590 /* > \param[in] P */
591 /* > \verbatim */
592 /* >          P is INTEGER */
593 /* >          The number of rows of the matrix B. 0 <= P <= N <= M+P. */
594 /* > \endverbatim */
595 /* > */
596 /* > \param[in,out] A */
597 /* > \verbatim */
598 /* >          A is DOUBLE PRECISION array, dimension (LDA,N) */
599 /* >          On entry, the M-by-N matrix A. */
600 /* >          On exit, the elements on and above the diagonal of the array */
601 /* >          contain the f2cmin(M,N)-by-N upper trapezoidal matrix T. */
602 /* > \endverbatim */
603 /* > */
604 /* > \param[in] LDA */
605 /* > \verbatim */
606 /* >          LDA is INTEGER */
607 /* >          The leading dimension of the array A. LDA >= f2cmax(1,M). */
608 /* > \endverbatim */
609 /* > */
610 /* > \param[in,out] B */
611 /* > \verbatim */
612 /* >          B is DOUBLE PRECISION array, dimension (LDB,N) */
613 /* >          On entry, the P-by-N matrix B. */
614 /* >          On exit, the upper triangle of the subarray B(1:P,N-P+1:N) */
615 /* >          contains the P-by-P upper triangular matrix R. */
616 /* > \endverbatim */
617 /* > */
618 /* > \param[in] LDB */
619 /* > \verbatim */
620 /* >          LDB is INTEGER */
621 /* >          The leading dimension of the array B. LDB >= f2cmax(1,P). */
622 /* > \endverbatim */
623 /* > */
624 /* > \param[in,out] C */
625 /* > \verbatim */
626 /* >          C is DOUBLE PRECISION array, dimension (M) */
627 /* >          On entry, C contains the right hand side vector for the */
628 /* >          least squares part of the LSE problem. */
629 /* >          On exit, the residual sum of squares for the solution */
630 /* >          is given by the sum of squares of elements N-P+1 to M of */
631 /* >          vector C. */
632 /* > \endverbatim */
633 /* > */
634 /* > \param[in,out] D */
635 /* > \verbatim */
636 /* >          D is DOUBLE PRECISION array, dimension (P) */
637 /* >          On entry, D contains the right hand side vector for the */
638 /* >          constrained equation. */
639 /* >          On exit, D is destroyed. */
640 /* > \endverbatim */
641 /* > */
642 /* > \param[out] X */
643 /* > \verbatim */
644 /* >          X is DOUBLE PRECISION array, dimension (N) */
645 /* >          On exit, X is the solution of the LSE problem. */
646 /* > \endverbatim */
647 /* > */
648 /* > \param[out] WORK */
649 /* > \verbatim */
650 /* >          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
651 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
652 /* > \endverbatim */
653 /* > */
654 /* > \param[in] LWORK */
655 /* > \verbatim */
656 /* >          LWORK is INTEGER */
657 /* >          The dimension of the array WORK. LWORK >= f2cmax(1,M+N+P). */
658 /* >          For optimum performance LWORK >= P+f2cmin(M,N)+f2cmax(M,N)*NB, */
659 /* >          where NB is an upper bound for the optimal blocksizes for */
660 /* >          DGEQRF, SGERQF, DORMQR and SORMRQ. */
661 /* > */
662 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
663 /* >          only calculates the optimal size of the WORK array, returns */
664 /* >          this value as the first entry of the WORK array, and no error */
665 /* >          message related to LWORK is issued by XERBLA. */
666 /* > \endverbatim */
667 /* > */
668 /* > \param[out] INFO */
669 /* > \verbatim */
670 /* >          INFO is INTEGER */
671 /* >          = 0:  successful exit. */
672 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
673 /* >          = 1:  the upper triangular factor R associated with B in the */
674 /* >                generalized RQ factorization of the pair (B, A) is */
675 /* >                singular, so that rank(B) < P; the least squares */
676 /* >                solution could not be computed. */
677 /* >          = 2:  the (N-P) by (N-P) part of the upper trapezoidal factor */
678 /* >                T associated with A in the generalized RQ factorization */
679 /* >                of the pair (B, A) is singular, so that */
680 /* >                rank( (A) ) < N; the least squares solution could not */
681 /* >                    ( (B) ) */
682 /* >                be computed. */
683 /* > \endverbatim */
684
685 /*  Authors: */
686 /*  ======== */
687
688 /* > \author Univ. of Tennessee */
689 /* > \author Univ. of California Berkeley */
690 /* > \author Univ. of Colorado Denver */
691 /* > \author NAG Ltd. */
692
693 /* > \date December 2016 */
694
695 /* > \ingroup doubleOTHERsolve */
696
697 /*  ===================================================================== */
698 /* Subroutine */ int dgglse_(integer *m, integer *n, integer *p, doublereal *
699         a, integer *lda, doublereal *b, integer *ldb, doublereal *c__, 
700         doublereal *d__, doublereal *x, doublereal *work, integer *lwork, 
701         integer *info)
702 {
703     /* System generated locals */
704     integer a_dim1, a_offset, b_dim1, b_offset, i__1, i__2;
705
706     /* Local variables */
707     integer lopt;
708     extern /* Subroutine */ int dgemv_(char *, integer *, integer *, 
709             doublereal *, doublereal *, integer *, doublereal *, integer *, 
710             doublereal *, doublereal *, integer *), dcopy_(integer *, 
711             doublereal *, integer *, doublereal *, integer *), daxpy_(integer 
712             *, doublereal *, doublereal *, integer *, doublereal *, integer *)
713             , dtrmv_(char *, char *, char *, integer *, doublereal *, integer 
714             *, doublereal *, integer *);
715     integer nb, mn, nr;
716     extern /* Subroutine */ int dggrqf_(integer *, integer *, integer *, 
717             doublereal *, integer *, doublereal *, doublereal *, integer *, 
718             doublereal *, doublereal *, integer *, integer *), xerbla_(char *,
719              integer *, ftnlen);
720     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
721             integer *, integer *, ftnlen, ftnlen);
722     integer lwkmin, nb1, nb2, nb3, nb4;
723     extern /* Subroutine */ int dormqr_(char *, char *, integer *, integer *, 
724             integer *, doublereal *, integer *, doublereal *, doublereal *, 
725             integer *, doublereal *, integer *, integer *), 
726             dormrq_(char *, char *, integer *, integer *, integer *, 
727             doublereal *, integer *, doublereal *, doublereal *, integer *, 
728             doublereal *, integer *, integer *);
729     integer lwkopt;
730     logical lquery;
731     extern /* Subroutine */ int dtrtrs_(char *, char *, char *, integer *, 
732             integer *, doublereal *, integer *, doublereal *, integer *, 
733             integer *);
734
735
736 /*  -- LAPACK driver routine (version 3.7.0) -- */
737 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
738 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
739 /*     December 2016 */
740
741
742 /*  ===================================================================== */
743
744
745 /*     Test the input parameters */
746
747     /* Parameter adjustments */
748     a_dim1 = *lda;
749     a_offset = 1 + a_dim1 * 1;
750     a -= a_offset;
751     b_dim1 = *ldb;
752     b_offset = 1 + b_dim1 * 1;
753     b -= b_offset;
754     --c__;
755     --d__;
756     --x;
757     --work;
758
759     /* Function Body */
760     *info = 0;
761     mn = f2cmin(*m,*n);
762     lquery = *lwork == -1;
763     if (*m < 0) {
764         *info = -1;
765     } else if (*n < 0) {
766         *info = -2;
767     } else if (*p < 0 || *p > *n || *p < *n - *m) {
768         *info = -3;
769     } else if (*lda < f2cmax(1,*m)) {
770         *info = -5;
771     } else if (*ldb < f2cmax(1,*p)) {
772         *info = -7;
773     }
774
775 /*     Calculate workspace */
776
777     if (*info == 0) {
778         if (*n == 0) {
779             lwkmin = 1;
780             lwkopt = 1;
781         } else {
782             nb1 = ilaenv_(&c__1, "DGEQRF", " ", m, n, &c_n1, &c_n1, (ftnlen)6,
783                      (ftnlen)1);
784             nb2 = ilaenv_(&c__1, "DGERQF", " ", m, n, &c_n1, &c_n1, (ftnlen)6,
785                      (ftnlen)1);
786             nb3 = ilaenv_(&c__1, "DORMQR", " ", m, n, p, &c_n1, (ftnlen)6, (
787                     ftnlen)1);
788             nb4 = ilaenv_(&c__1, "DORMRQ", " ", m, n, p, &c_n1, (ftnlen)6, (
789                     ftnlen)1);
790 /* Computing MAX */
791             i__1 = f2cmax(nb1,nb2), i__1 = f2cmax(i__1,nb3);
792             nb = f2cmax(i__1,nb4);
793             lwkmin = *m + *n + *p;
794             lwkopt = *p + mn + f2cmax(*m,*n) * nb;
795         }
796         work[1] = (doublereal) lwkopt;
797
798         if (*lwork < lwkmin && ! lquery) {
799             *info = -12;
800         }
801     }
802
803     if (*info != 0) {
804         i__1 = -(*info);
805         xerbla_("DGGLSE", &i__1, (ftnlen)6);
806         return 0;
807     } else if (lquery) {
808         return 0;
809     }
810
811 /*     Quick return if possible */
812
813     if (*n == 0) {
814         return 0;
815     }
816
817 /*     Compute the GRQ factorization of matrices B and A: */
818
819 /*            B*Q**T = (  0  T12 ) P   Z**T*A*Q**T = ( R11 R12 ) N-P */
820 /*                        N-P  P                     (  0  R22 ) M+P-N */
821 /*                                                      N-P  P */
822
823 /*     where T12 and R11 are upper triangular, and Q and Z are */
824 /*     orthogonal. */
825
826     i__1 = *lwork - *p - mn;
827     dggrqf_(p, m, n, &b[b_offset], ldb, &work[1], &a[a_offset], lda, &work[*p 
828             + 1], &work[*p + mn + 1], &i__1, info);
829     lopt = (integer) work[*p + mn + 1];
830
831 /*     Update c = Z**T *c = ( c1 ) N-P */
832 /*                          ( c2 ) M+P-N */
833
834     i__1 = f2cmax(1,*m);
835     i__2 = *lwork - *p - mn;
836     dormqr_("Left", "Transpose", m, &c__1, &mn, &a[a_offset], lda, &work[*p + 
837             1], &c__[1], &i__1, &work[*p + mn + 1], &i__2, info);
838 /* Computing MAX */
839     i__1 = lopt, i__2 = (integer) work[*p + mn + 1];
840     lopt = f2cmax(i__1,i__2);
841
842 /*     Solve T12*x2 = d for x2 */
843
844     if (*p > 0) {
845         dtrtrs_("Upper", "No transpose", "Non-unit", p, &c__1, &b[(*n - *p + 
846                 1) * b_dim1 + 1], ldb, &d__[1], p, info);
847
848         if (*info > 0) {
849             *info = 1;
850             return 0;
851         }
852
853 /*        Put the solution in X */
854
855         dcopy_(p, &d__[1], &c__1, &x[*n - *p + 1], &c__1);
856
857 /*        Update c1 */
858
859         i__1 = *n - *p;
860         dgemv_("No transpose", &i__1, p, &c_b31, &a[(*n - *p + 1) * a_dim1 + 
861                 1], lda, &d__[1], &c__1, &c_b33, &c__[1], &c__1);
862     }
863
864 /*     Solve R11*x1 = c1 for x1 */
865
866     if (*n > *p) {
867         i__1 = *n - *p;
868         i__2 = *n - *p;
869         dtrtrs_("Upper", "No transpose", "Non-unit", &i__1, &c__1, &a[
870                 a_offset], lda, &c__[1], &i__2, info);
871
872         if (*info > 0) {
873             *info = 2;
874             return 0;
875         }
876
877 /*        Put the solutions in X */
878
879         i__1 = *n - *p;
880         dcopy_(&i__1, &c__[1], &c__1, &x[1], &c__1);
881     }
882
883 /*     Compute the residual vector: */
884
885     if (*m < *n) {
886         nr = *m + *p - *n;
887         if (nr > 0) {
888             i__1 = *n - *m;
889             dgemv_("No transpose", &nr, &i__1, &c_b31, &a[*n - *p + 1 + (*m + 
890                     1) * a_dim1], lda, &d__[nr + 1], &c__1, &c_b33, &c__[*n - 
891                     *p + 1], &c__1);
892         }
893     } else {
894         nr = *p;
895     }
896     if (nr > 0) {
897         dtrmv_("Upper", "No transpose", "Non unit", &nr, &a[*n - *p + 1 + (*n 
898                 - *p + 1) * a_dim1], lda, &d__[1], &c__1);
899         daxpy_(&nr, &c_b31, &d__[1], &c__1, &c__[*n - *p + 1], &c__1);
900     }
901
902 /*     Backward transformation x = Q**T*x */
903
904     i__1 = *lwork - *p - mn;
905     dormrq_("Left", "Transpose", n, &c__1, p, &b[b_offset], ldb, &work[1], &x[
906             1], n, &work[*p + mn + 1], &i__1, info);
907 /* Computing MAX */
908     i__1 = lopt, i__2 = (integer) work[*p + mn + 1];
909     work[1] = (doublereal) (*p + mn + f2cmax(i__1,i__2));
910
911     return 0;
912
913 /*     End of DGGLSE */
914
915 } /* dgglse_ */
916