1090 lines
		
	
	
		
			32 KiB
		
	
	
	
		
			C
		
	
	
	
			
		
		
	
	
			1090 lines
		
	
	
		
			32 KiB
		
	
	
	
		
			C
		
	
	
	
| #include <math.h>
 | |
| #include <stdlib.h>
 | |
| #include <string.h>
 | |
| #include <stdio.h>
 | |
| #include <complex.h>
 | |
| #ifdef complex
 | |
| #undef complex
 | |
| #endif
 | |
| #ifdef I
 | |
| #undef I
 | |
| #endif
 | |
| 
 | |
| #if defined(_WIN64)
 | |
| typedef long long BLASLONG;
 | |
| typedef unsigned long long BLASULONG;
 | |
| #else
 | |
| typedef long BLASLONG;
 | |
| typedef unsigned long BLASULONG;
 | |
| #endif
 | |
| 
 | |
| #ifdef LAPACK_ILP64
 | |
| typedef BLASLONG blasint;
 | |
| #if defined(_WIN64)
 | |
| #define blasabs(x) llabs(x)
 | |
| #else
 | |
| #define blasabs(x) labs(x)
 | |
| #endif
 | |
| #else
 | |
| typedef int blasint;
 | |
| #define blasabs(x) abs(x)
 | |
| #endif
 | |
| 
 | |
| typedef blasint integer;
 | |
| 
 | |
| typedef unsigned int uinteger;
 | |
| typedef char *address;
 | |
| typedef short int shortint;
 | |
| typedef float real;
 | |
| typedef double doublereal;
 | |
| typedef struct { real r, i; } complex;
 | |
| typedef struct { doublereal r, i; } doublecomplex;
 | |
| #ifdef _MSC_VER
 | |
| static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
 | |
| static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
 | |
| static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
 | |
| static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
 | |
| #else
 | |
| static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
 | |
| static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
 | |
| static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
 | |
| static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
 | |
| #endif
 | |
| #define pCf(z) (*_pCf(z))
 | |
| #define pCd(z) (*_pCd(z))
 | |
| typedef int logical;
 | |
| typedef short int shortlogical;
 | |
| typedef char logical1;
 | |
| typedef char integer1;
 | |
| 
 | |
| #define TRUE_ (1)
 | |
| #define FALSE_ (0)
 | |
| 
 | |
| /* Extern is for use with -E */
 | |
| #ifndef Extern
 | |
| #define Extern extern
 | |
| #endif
 | |
| 
 | |
| /* I/O stuff */
 | |
| 
 | |
| typedef int flag;
 | |
| typedef int ftnlen;
 | |
| typedef int ftnint;
 | |
| 
 | |
| /*external read, write*/
 | |
| typedef struct
 | |
| {	flag cierr;
 | |
| 	ftnint ciunit;
 | |
| 	flag ciend;
 | |
| 	char *cifmt;
 | |
| 	ftnint cirec;
 | |
| } cilist;
 | |
| 
 | |
| /*internal read, write*/
 | |
| typedef struct
 | |
| {	flag icierr;
 | |
| 	char *iciunit;
 | |
| 	flag iciend;
 | |
| 	char *icifmt;
 | |
| 	ftnint icirlen;
 | |
| 	ftnint icirnum;
 | |
| } icilist;
 | |
| 
 | |
| /*open*/
 | |
| typedef struct
 | |
| {	flag oerr;
 | |
| 	ftnint ounit;
 | |
| 	char *ofnm;
 | |
| 	ftnlen ofnmlen;
 | |
| 	char *osta;
 | |
| 	char *oacc;
 | |
| 	char *ofm;
 | |
| 	ftnint orl;
 | |
| 	char *oblnk;
 | |
| } olist;
 | |
| 
 | |
| /*close*/
 | |
| typedef struct
 | |
| {	flag cerr;
 | |
| 	ftnint cunit;
 | |
| 	char *csta;
 | |
| } cllist;
 | |
| 
 | |
| /*rewind, backspace, endfile*/
 | |
| typedef struct
 | |
| {	flag aerr;
 | |
| 	ftnint aunit;
 | |
| } alist;
 | |
| 
 | |
| /* inquire */
 | |
| typedef struct
 | |
| {	flag inerr;
 | |
| 	ftnint inunit;
 | |
| 	char *infile;
 | |
| 	ftnlen infilen;
 | |
| 	ftnint	*inex;	/*parameters in standard's order*/
 | |
| 	ftnint	*inopen;
 | |
| 	ftnint	*innum;
 | |
| 	ftnint	*innamed;
 | |
| 	char	*inname;
 | |
| 	ftnlen	innamlen;
 | |
| 	char	*inacc;
 | |
| 	ftnlen	inacclen;
 | |
| 	char	*inseq;
 | |
| 	ftnlen	inseqlen;
 | |
| 	char 	*indir;
 | |
| 	ftnlen	indirlen;
 | |
| 	char	*infmt;
 | |
| 	ftnlen	infmtlen;
 | |
| 	char	*inform;
 | |
| 	ftnint	informlen;
 | |
| 	char	*inunf;
 | |
| 	ftnlen	inunflen;
 | |
| 	ftnint	*inrecl;
 | |
| 	ftnint	*innrec;
 | |
| 	char	*inblank;
 | |
| 	ftnlen	inblanklen;
 | |
| } inlist;
 | |
| 
 | |
| #define VOID void
 | |
| 
 | |
| union Multitype {	/* for multiple entry points */
 | |
| 	integer1 g;
 | |
| 	shortint h;
 | |
| 	integer i;
 | |
| 	/* longint j; */
 | |
| 	real r;
 | |
| 	doublereal d;
 | |
| 	complex c;
 | |
| 	doublecomplex z;
 | |
| 	};
 | |
| 
 | |
| typedef union Multitype Multitype;
 | |
| 
 | |
| struct Vardesc {	/* for Namelist */
 | |
| 	char *name;
 | |
| 	char *addr;
 | |
| 	ftnlen *dims;
 | |
| 	int  type;
 | |
| 	};
 | |
| typedef struct Vardesc Vardesc;
 | |
| 
 | |
| struct Namelist {
 | |
| 	char *name;
 | |
| 	Vardesc **vars;
 | |
| 	int nvars;
 | |
| 	};
 | |
| typedef struct Namelist Namelist;
 | |
| 
 | |
| #define abs(x) ((x) >= 0 ? (x) : -(x))
 | |
| #define dabs(x) (fabs(x))
 | |
| #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
 | |
| #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
 | |
| #define dmin(a,b) (f2cmin(a,b))
 | |
| #define dmax(a,b) (f2cmax(a,b))
 | |
| #define bit_test(a,b)	((a) >> (b) & 1)
 | |
| #define bit_clear(a,b)	((a) & ~((uinteger)1 << (b)))
 | |
| #define bit_set(a,b)	((a) |  ((uinteger)1 << (b)))
 | |
| 
 | |
| #define abort_() { sig_die("Fortran abort routine called", 1); }
 | |
| #define c_abs(z) (cabsf(Cf(z)))
 | |
| #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
 | |
| #ifdef _MSC_VER
 | |
| #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]);}
 | |
| #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]);}
 | |
| #else
 | |
| #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
 | |
| #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
 | |
| #endif
 | |
| #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
 | |
| #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
 | |
| #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
 | |
| //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
 | |
| #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
 | |
| #define d_abs(x) (fabs(*(x)))
 | |
| #define d_acos(x) (acos(*(x)))
 | |
| #define d_asin(x) (asin(*(x)))
 | |
| #define d_atan(x) (atan(*(x)))
 | |
| #define d_atn2(x, y) (atan2(*(x),*(y)))
 | |
| #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
 | |
| #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
 | |
| #define d_cos(x) (cos(*(x)))
 | |
| #define d_cosh(x) (cosh(*(x)))
 | |
| #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
 | |
| #define d_exp(x) (exp(*(x)))
 | |
| #define d_imag(z) (cimag(Cd(z)))
 | |
| #define r_imag(z) (cimagf(Cf(z)))
 | |
| #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
 | |
| #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
 | |
| #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
 | |
| #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
 | |
| #define d_log(x) (log(*(x)))
 | |
| #define d_mod(x, y) (fmod(*(x), *(y)))
 | |
| #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
 | |
| #define d_nint(x) u_nint(*(x))
 | |
| #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
 | |
| #define d_sign(a,b) u_sign(*(a),*(b))
 | |
| #define r_sign(a,b) u_sign(*(a),*(b))
 | |
| #define d_sin(x) (sin(*(x)))
 | |
| #define d_sinh(x) (sinh(*(x)))
 | |
| #define d_sqrt(x) (sqrt(*(x)))
 | |
| #define d_tan(x) (tan(*(x)))
 | |
| #define d_tanh(x) (tanh(*(x)))
 | |
| #define i_abs(x) abs(*(x))
 | |
| #define i_dnnt(x) ((integer)u_nint(*(x)))
 | |
| #define i_len(s, n) (n)
 | |
| #define i_nint(x) ((integer)u_nint(*(x)))
 | |
| #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
 | |
| #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
 | |
| #define pow_si(B,E) spow_ui(*(B),*(E))
 | |
| #define pow_ri(B,E) spow_ui(*(B),*(E))
 | |
| #define pow_di(B,E) dpow_ui(*(B),*(E))
 | |
| #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
 | |
| #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
 | |
| #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
 | |
| #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++ = ' '; }
 | |
| #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
 | |
| #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]; }
 | |
| #define sig_die(s, kill) { exit(1); }
 | |
| #define s_stop(s, n) {exit(0);}
 | |
| static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
 | |
| #define z_abs(z) (cabs(Cd(z)))
 | |
| #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
 | |
| #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
 | |
| #define myexit_() break;
 | |
| #define mycycle() continue;
 | |
| #define myceiling(w) {ceil(w)}
 | |
| #define myhuge(w) {HUGE_VAL}
 | |
| //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
 | |
| #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
 | |
| 
 | |
| /* procedure parameter types for -A and -C++ */
 | |
| 
 | |
| #define F2C_proc_par_types 1
 | |
| #ifdef __cplusplus
 | |
| typedef logical (*L_fp)(...);
 | |
| #else
 | |
| typedef logical (*L_fp)();
 | |
| #endif
 | |
| 
 | |
| static float spow_ui(float x, integer n) {
 | |
| 	float pow=1.0; unsigned long int u;
 | |
| 	if(n != 0) {
 | |
| 		if(n < 0) n = -n, x = 1/x;
 | |
| 		for(u = n; ; ) {
 | |
| 			if(u & 01) pow *= x;
 | |
| 			if(u >>= 1) x *= x;
 | |
| 			else break;
 | |
| 		}
 | |
| 	}
 | |
| 	return pow;
 | |
| }
 | |
| static double dpow_ui(double x, integer n) {
 | |
| 	double pow=1.0; unsigned long int u;
 | |
| 	if(n != 0) {
 | |
| 		if(n < 0) n = -n, x = 1/x;
 | |
| 		for(u = n; ; ) {
 | |
| 			if(u & 01) pow *= x;
 | |
| 			if(u >>= 1) x *= x;
 | |
| 			else break;
 | |
| 		}
 | |
| 	}
 | |
| 	return pow;
 | |
| }
 | |
| #ifdef _MSC_VER
 | |
| static _Fcomplex cpow_ui(complex x, integer n) {
 | |
| 	complex pow={1.0,0.0}; unsigned long int u;
 | |
| 		if(n != 0) {
 | |
| 		if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
 | |
| 		for(u = n; ; ) {
 | |
| 			if(u & 01) pow.r *= x.r, pow.i *= x.i;
 | |
| 			if(u >>= 1) x.r *= x.r, x.i *= x.i;
 | |
| 			else break;
 | |
| 		}
 | |
| 	}
 | |
| 	_Fcomplex p={pow.r, pow.i};
 | |
| 	return p;
 | |
| }
 | |
| #else
 | |
| static _Complex float cpow_ui(_Complex float x, integer n) {
 | |
| 	_Complex float pow=1.0; unsigned long int u;
 | |
| 	if(n != 0) {
 | |
| 		if(n < 0) n = -n, x = 1/x;
 | |
| 		for(u = n; ; ) {
 | |
| 			if(u & 01) pow *= x;
 | |
| 			if(u >>= 1) x *= x;
 | |
| 			else break;
 | |
| 		}
 | |
| 	}
 | |
| 	return pow;
 | |
| }
 | |
| #endif
 | |
| #ifdef _MSC_VER
 | |
| static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
 | |
| 	_Dcomplex pow={1.0,0.0}; unsigned long int u;
 | |
| 	if(n != 0) {
 | |
| 		if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
 | |
| 		for(u = n; ; ) {
 | |
| 			if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
 | |
| 			if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
 | |
| 			else break;
 | |
| 		}
 | |
| 	}
 | |
| 	_Dcomplex p = {pow._Val[0], pow._Val[1]};
 | |
| 	return p;
 | |
| }
 | |
| #else
 | |
| static _Complex double zpow_ui(_Complex double x, integer n) {
 | |
| 	_Complex double pow=1.0; unsigned long int u;
 | |
| 	if(n != 0) {
 | |
| 		if(n < 0) n = -n, x = 1/x;
 | |
| 		for(u = n; ; ) {
 | |
| 			if(u & 01) pow *= x;
 | |
| 			if(u >>= 1) x *= x;
 | |
| 			else break;
 | |
| 		}
 | |
| 	}
 | |
| 	return pow;
 | |
| }
 | |
| #endif
 | |
| static integer pow_ii(integer x, integer n) {
 | |
| 	integer pow; unsigned long int u;
 | |
| 	if (n <= 0) {
 | |
| 		if (n == 0 || x == 1) pow = 1;
 | |
| 		else if (x != -1) pow = x == 0 ? 1/x : 0;
 | |
| 		else n = -n;
 | |
| 	}
 | |
| 	if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
 | |
| 		u = n;
 | |
| 		for(pow = 1; ; ) {
 | |
| 			if(u & 01) pow *= x;
 | |
| 			if(u >>= 1) x *= x;
 | |
| 			else break;
 | |
| 		}
 | |
| 	}
 | |
| 	return pow;
 | |
| }
 | |
| static integer dmaxloc_(double *w, integer s, integer e, integer *n)
 | |
| {
 | |
| 	double m; integer i, mi;
 | |
| 	for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
 | |
| 		if (w[i-1]>m) mi=i ,m=w[i-1];
 | |
| 	return mi-s+1;
 | |
| }
 | |
| static integer smaxloc_(float *w, integer s, integer e, integer *n)
 | |
| {
 | |
| 	float m; integer i, mi;
 | |
| 	for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
 | |
| 		if (w[i-1]>m) mi=i ,m=w[i-1];
 | |
| 	return mi-s+1;
 | |
| }
 | |
| static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
 | |
| 	integer n = *n_, incx = *incx_, incy = *incy_, i;
 | |
| #ifdef _MSC_VER
 | |
| 	_Fcomplex zdotc = {0.0, 0.0};
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
 | |
| 			zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
 | |
| 			zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
 | |
| 		}
 | |
| 	}
 | |
| 	pCf(z) = zdotc;
 | |
| }
 | |
| #else
 | |
| 	_Complex float zdotc = 0.0;
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
 | |
| 		}
 | |
| 	}
 | |
| 	pCf(z) = zdotc;
 | |
| }
 | |
| #endif
 | |
| static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
 | |
| 	integer n = *n_, incx = *incx_, incy = *incy_, i;
 | |
| #ifdef _MSC_VER
 | |
| 	_Dcomplex zdotc = {0.0, 0.0};
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
 | |
| 			zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
 | |
| 			zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
 | |
| 		}
 | |
| 	}
 | |
| 	pCd(z) = zdotc;
 | |
| }
 | |
| #else
 | |
| 	_Complex double zdotc = 0.0;
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
 | |
| 		}
 | |
| 	}
 | |
| 	pCd(z) = zdotc;
 | |
| }
 | |
| #endif	
 | |
| static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
 | |
| 	integer n = *n_, incx = *incx_, incy = *incy_, i;
 | |
| #ifdef _MSC_VER
 | |
| 	_Fcomplex zdotc = {0.0, 0.0};
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
 | |
| 			zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
 | |
| 			zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
 | |
| 		}
 | |
| 	}
 | |
| 	pCf(z) = zdotc;
 | |
| }
 | |
| #else
 | |
| 	_Complex float zdotc = 0.0;
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += Cf(&x[i]) * Cf(&y[i]);
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
 | |
| 		}
 | |
| 	}
 | |
| 	pCf(z) = zdotc;
 | |
| }
 | |
| #endif
 | |
| static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
 | |
| 	integer n = *n_, incx = *incx_, incy = *incy_, i;
 | |
| #ifdef _MSC_VER
 | |
| 	_Dcomplex zdotc = {0.0, 0.0};
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
 | |
| 			zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
 | |
| 			zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
 | |
| 		}
 | |
| 	}
 | |
| 	pCd(z) = zdotc;
 | |
| }
 | |
| #else
 | |
| 	_Complex double zdotc = 0.0;
 | |
| 	if (incx == 1 && incy == 1) {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += Cd(&x[i]) * Cd(&y[i]);
 | |
| 		}
 | |
| 	} else {
 | |
| 		for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
 | |
| 			zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
 | |
| 		}
 | |
| 	}
 | |
| 	pCd(z) = zdotc;
 | |
| }
 | |
| #endif
 | |
| /*  -- translated by f2c (version 20000121).
 | |
|    You must link the resulting object file with the libraries:
 | |
| 	-lf2c -lm   (in that order)
 | |
| */
 | |
| 
 | |
| 
 | |
| 
 | |
| 
 | |
| /* Table of constant values */
 | |
| 
 | |
| static doublecomplex c_b1 = {0.,0.};
 | |
| static doublecomplex c_b2 = {1.,0.};
 | |
| static integer c__1 = 1;
 | |
| static integer c_n1 = -1;
 | |
| static integer c__0 = 0;
 | |
| static integer c__2 = 2;
 | |
| 
 | |
| /* > \brief <b> ZGELSY solves overdetermined or underdetermined systems for GE matrices</b> */
 | |
| 
 | |
| /*  =========== DOCUMENTATION =========== */
 | |
| 
 | |
| /* Online html documentation available at */
 | |
| /*            http://www.netlib.org/lapack/explore-html/ */
 | |
| 
 | |
| /* > \htmlonly */
 | |
| /* > Download ZGELSY + dependencies */
 | |
| /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgelsy.
 | |
| f"> */
 | |
| /* > [TGZ]</a> */
 | |
| /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgelsy.
 | |
| f"> */
 | |
| /* > [ZIP]</a> */
 | |
| /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgelsy.
 | |
| f"> */
 | |
| /* > [TXT]</a> */
 | |
| /* > \endhtmlonly */
 | |
| 
 | |
| /*  Definition: */
 | |
| /*  =========== */
 | |
| 
 | |
| /*       SUBROUTINE ZGELSY( M, N, NRHS, A, LDA, B, LDB, JPVT, RCOND, RANK, */
 | |
| /*                          WORK, LWORK, RWORK, INFO ) */
 | |
| 
 | |
| /*       INTEGER            INFO, LDA, LDB, LWORK, M, N, NRHS, RANK */
 | |
| /*       DOUBLE PRECISION   RCOND */
 | |
| /*       INTEGER            JPVT( * ) */
 | |
| /*       DOUBLE PRECISION   RWORK( * ) */
 | |
| /*       COMPLEX*16         A( LDA, * ), B( LDB, * ), WORK( * ) */
 | |
| 
 | |
| 
 | |
| /* > \par Purpose: */
 | |
| /*  ============= */
 | |
| /* > */
 | |
| /* > \verbatim */
 | |
| /* > */
 | |
| /* > ZGELSY computes the minimum-norm solution to a complex linear least */
 | |
| /* > squares problem: */
 | |
| /* >     minimize || A * X - B || */
 | |
| /* > using a complete orthogonal factorization of A.  A is an M-by-N */
 | |
| /* > matrix which may be rank-deficient. */
 | |
| /* > */
 | |
| /* > Several right hand side vectors b and solution vectors x can be */
 | |
| /* > handled in a single call; they are stored as the columns of the */
 | |
| /* > M-by-NRHS right hand side matrix B and the N-by-NRHS solution */
 | |
| /* > matrix X. */
 | |
| /* > */
 | |
| /* > The routine first computes a QR factorization with column pivoting: */
 | |
| /* >     A * P = Q * [ R11 R12 ] */
 | |
| /* >                 [  0  R22 ] */
 | |
| /* > with R11 defined as the largest leading submatrix whose estimated */
 | |
| /* > condition number is less than 1/RCOND.  The order of R11, RANK, */
 | |
| /* > is the effective rank of A. */
 | |
| /* > */
 | |
| /* > Then, R22 is considered to be negligible, and R12 is annihilated */
 | |
| /* > by unitary transformations from the right, arriving at the */
 | |
| /* > complete orthogonal factorization: */
 | |
| /* >    A * P = Q * [ T11 0 ] * Z */
 | |
| /* >                [  0  0 ] */
 | |
| /* > The minimum-norm solution is then */
 | |
| /* >    X = P * Z**H [ inv(T11)*Q1**H*B ] */
 | |
| /* >                 [        0         ] */
 | |
| /* > where Q1 consists of the first RANK columns of Q. */
 | |
| /* > */
 | |
| /* > This routine is basically identical to the original xGELSX except */
 | |
| /* > three differences: */
 | |
| /* >   o The permutation of matrix B (the right hand side) is faster and */
 | |
| /* >     more simple. */
 | |
| /* >   o The call to the subroutine xGEQPF has been substituted by the */
 | |
| /* >     the call to the subroutine xGEQP3. This subroutine is a Blas-3 */
 | |
| /* >     version of the QR factorization with column pivoting. */
 | |
| /* >   o Matrix B (the right hand side) is updated with Blas-3. */
 | |
| /* > \endverbatim */
 | |
| 
 | |
| /*  Arguments: */
 | |
| /*  ========== */
 | |
| 
 | |
| /* > \param[in] M */
 | |
| /* > \verbatim */
 | |
| /* >          M is INTEGER */
 | |
| /* >          The number of rows of the matrix A.  M >= 0. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in] N */
 | |
| /* > \verbatim */
 | |
| /* >          N is INTEGER */
 | |
| /* >          The number of columns of the matrix A.  N >= 0. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in] NRHS */
 | |
| /* > \verbatim */
 | |
| /* >          NRHS is INTEGER */
 | |
| /* >          The number of right hand sides, i.e., the number of */
 | |
| /* >          columns of matrices B and X. NRHS >= 0. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in,out] A */
 | |
| /* > \verbatim */
 | |
| /* >          A is COMPLEX*16 array, dimension (LDA,N) */
 | |
| /* >          On entry, the M-by-N matrix A. */
 | |
| /* >          On exit, A has been overwritten by details of its */
 | |
| /* >          complete orthogonal factorization. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in] LDA */
 | |
| /* > \verbatim */
 | |
| /* >          LDA is INTEGER */
 | |
| /* >          The leading dimension of the array A.  LDA >= f2cmax(1,M). */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in,out] B */
 | |
| /* > \verbatim */
 | |
| /* >          B is COMPLEX*16 array, dimension (LDB,NRHS) */
 | |
| /* >          On entry, the M-by-NRHS right hand side matrix B. */
 | |
| /* >          On exit, the N-by-NRHS solution matrix X. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in] LDB */
 | |
| /* > \verbatim */
 | |
| /* >          LDB is INTEGER */
 | |
| /* >          The leading dimension of the array B. LDB >= f2cmax(1,M,N). */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in,out] JPVT */
 | |
| /* > \verbatim */
 | |
| /* >          JPVT is INTEGER array, dimension (N) */
 | |
| /* >          On entry, if JPVT(i) .ne. 0, the i-th column of A is permuted */
 | |
| /* >          to the front of AP, otherwise column i is a free column. */
 | |
| /* >          On exit, if JPVT(i) = k, then the i-th column of A*P */
 | |
| /* >          was the k-th column of A. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in] RCOND */
 | |
| /* > \verbatim */
 | |
| /* >          RCOND is DOUBLE PRECISION */
 | |
| /* >          RCOND is used to determine the effective rank of A, which */
 | |
| /* >          is defined as the order of the largest leading triangular */
 | |
| /* >          submatrix R11 in the QR factorization with pivoting of A, */
 | |
| /* >          whose estimated condition number < 1/RCOND. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[out] RANK */
 | |
| /* > \verbatim */
 | |
| /* >          RANK is INTEGER */
 | |
| /* >          The effective rank of A, i.e., the order of the submatrix */
 | |
| /* >          R11.  This is the same as the order of the submatrix T11 */
 | |
| /* >          in the complete orthogonal factorization of A. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[out] WORK */
 | |
| /* > \verbatim */
 | |
| /* >          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) */
 | |
| /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[in] LWORK */
 | |
| /* > \verbatim */
 | |
| /* >          LWORK is INTEGER */
 | |
| /* >          The dimension of the array WORK. */
 | |
| /* >          The unblocked strategy requires that: */
 | |
| /* >            LWORK >= MN + MAX( 2*MN, N+1, MN+NRHS ) */
 | |
| /* >          where MN = f2cmin(M,N). */
 | |
| /* >          The block algorithm requires that: */
 | |
| /* >            LWORK >= MN + MAX( 2*MN, NB*(N+1), MN+MN*NB, MN+NB*NRHS ) */
 | |
| /* >          where NB is an upper bound on the blocksize returned */
 | |
| /* >          by ILAENV for the routines ZGEQP3, ZTZRZF, CTZRQF, ZUNMQR, */
 | |
| /* >          and ZUNMRZ. */
 | |
| /* > */
 | |
| /* >          If LWORK = -1, then a workspace query is assumed; the routine */
 | |
| /* >          only calculates the optimal size of the WORK array, returns */
 | |
| /* >          this value as the first entry of the WORK array, and no error */
 | |
| /* >          message related to LWORK is issued by XERBLA. */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[out] RWORK */
 | |
| /* > \verbatim */
 | |
| /* >          RWORK is DOUBLE PRECISION array, dimension (2*N) */
 | |
| /* > \endverbatim */
 | |
| /* > */
 | |
| /* > \param[out] INFO */
 | |
| /* > \verbatim */
 | |
| /* >          INFO is INTEGER */
 | |
| /* >          = 0: successful exit */
 | |
| /* >          < 0: if INFO = -i, the i-th argument had an illegal value */
 | |
| /* > \endverbatim */
 | |
| 
 | |
| /*  Authors: */
 | |
| /*  ======== */
 | |
| 
 | |
| /* > \author Univ. of Tennessee */
 | |
| /* > \author Univ. of California Berkeley */
 | |
| /* > \author Univ. of Colorado Denver */
 | |
| /* > \author NAG Ltd. */
 | |
| 
 | |
| /* > \date December 2016 */
 | |
| 
 | |
| /* > \ingroup complex16GEsolve */
 | |
| 
 | |
| /* > \par Contributors: */
 | |
| /*  ================== */
 | |
| /* > */
 | |
| /* >    A. Petitet, Computer Science Dept., Univ. of Tenn., Knoxville, USA \n */
 | |
| /* >    E. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain \n */
 | |
| /* >    G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain \n */
 | |
| /* > */
 | |
| /*  ===================================================================== */
 | |
| /* Subroutine */ void zgelsy_(integer *m, integer *n, integer *nrhs, 
 | |
| 	doublecomplex *a, integer *lda, doublecomplex *b, integer *ldb, 
 | |
| 	integer *jpvt, doublereal *rcond, integer *rank, doublecomplex *work, 
 | |
| 	integer *lwork, doublereal *rwork, integer *info)
 | |
| {
 | |
|     /* System generated locals */
 | |
|     integer a_dim1, a_offset, b_dim1, b_offset, i__1, i__2, i__3, i__4;
 | |
|     doublereal d__1, d__2;
 | |
|     doublecomplex z__1;
 | |
| 
 | |
|     /* Local variables */
 | |
|     doublereal anrm, bnrm, smin, smax;
 | |
|     integer i__, j, iascl, ibscl, ismin, ismax;
 | |
|     doublecomplex c1, c2;
 | |
|     doublereal wsize;
 | |
|     doublecomplex s1, s2;
 | |
|     extern /* Subroutine */ void zcopy_(integer *, doublecomplex *, integer *, 
 | |
| 	    doublecomplex *, integer *), ztrsm_(char *, char *, char *, char *
 | |
| 	    , integer *, integer *, doublecomplex *, doublecomplex *, integer 
 | |
| 	    *, doublecomplex *, integer *), 
 | |
| 	    zlaic1_(integer *, integer *, doublecomplex *, doublereal *, 
 | |
| 	    doublecomplex *, doublecomplex *, doublereal *, doublecomplex *, 
 | |
| 	    doublecomplex *), dlabad_(doublereal *, doublereal *), zgeqp3_(
 | |
| 	    integer *, integer *, doublecomplex *, integer *, integer *, 
 | |
| 	    doublecomplex *, doublecomplex *, integer *, doublereal *, 
 | |
| 	    integer *);
 | |
|     integer nb;
 | |
|     extern doublereal dlamch_(char *);
 | |
|     integer mn;
 | |
|     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
 | |
|     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
 | |
| 	    integer *, integer *, ftnlen, ftnlen);
 | |
|     extern doublereal zlange_(char *, integer *, integer *, doublecomplex *, 
 | |
| 	    integer *, doublereal *);
 | |
|     doublereal bignum;
 | |
|     extern /* Subroutine */ void zlascl_(char *, integer *, integer *, 
 | |
| 	    doublereal *, doublereal *, integer *, integer *, doublecomplex *,
 | |
| 	     integer *, integer *);
 | |
|     integer nb1, nb2, nb3, nb4;
 | |
|     extern /* Subroutine */ void zlaset_(char *, integer *, integer *, 
 | |
| 	    doublecomplex *, doublecomplex *, doublecomplex *, integer *);
 | |
|     doublereal sminpr, smaxpr, smlnum;
 | |
|     integer lwkopt;
 | |
|     logical lquery;
 | |
|     extern /* Subroutine */ void zunmqr_(char *, char *, integer *, integer *, 
 | |
| 	    integer *, doublecomplex *, integer *, doublecomplex *, 
 | |
| 	    doublecomplex *, integer *, doublecomplex *, integer *, integer *), zunmrz_(char *, char *, integer *, integer *, 
 | |
| 	    integer *, integer *, doublecomplex *, integer *, doublecomplex *,
 | |
| 	     doublecomplex *, integer *, doublecomplex *, integer *, integer *
 | |
| 	    ), ztzrzf_(integer *, integer *, doublecomplex *, 
 | |
| 	    integer *, doublecomplex *, doublecomplex *, integer *, integer *)
 | |
| 	    ;
 | |
| 
 | |
| 
 | |
| /*  -- LAPACK driver routine (version 3.7.0) -- */
 | |
| /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
 | |
| /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
 | |
| /*     December 2016 */
 | |
| 
 | |
| 
 | |
| /*  ===================================================================== */
 | |
| 
 | |
| 
 | |
|     /* Parameter adjustments */
 | |
|     a_dim1 = *lda;
 | |
|     a_offset = 1 + a_dim1 * 1;
 | |
|     a -= a_offset;
 | |
|     b_dim1 = *ldb;
 | |
|     b_offset = 1 + b_dim1 * 1;
 | |
|     b -= b_offset;
 | |
|     --jpvt;
 | |
|     --work;
 | |
|     --rwork;
 | |
| 
 | |
|     /* Function Body */
 | |
|     mn = f2cmin(*m,*n);
 | |
|     ismin = mn + 1;
 | |
|     ismax = (mn << 1) + 1;
 | |
| 
 | |
| /*     Test the input arguments. */
 | |
| 
 | |
|     *info = 0;
 | |
|     nb1 = ilaenv_(&c__1, "ZGEQRF", " ", m, n, &c_n1, &c_n1, (ftnlen)6, (
 | |
| 	    ftnlen)1);
 | |
|     nb2 = ilaenv_(&c__1, "ZGERQF", " ", m, n, &c_n1, &c_n1, (ftnlen)6, (
 | |
| 	    ftnlen)1);
 | |
|     nb3 = ilaenv_(&c__1, "ZUNMQR", " ", m, n, nrhs, &c_n1, (ftnlen)6, (ftnlen)
 | |
| 	    1);
 | |
|     nb4 = ilaenv_(&c__1, "ZUNMRQ", " ", m, n, nrhs, &c_n1, (ftnlen)6, (ftnlen)
 | |
| 	    1);
 | |
| /* Computing MAX */
 | |
|     i__1 = f2cmax(nb1,nb2), i__1 = f2cmax(i__1,nb3);
 | |
|     nb = f2cmax(i__1,nb4);
 | |
| /* Computing MAX */
 | |
|     i__1 = 1, i__2 = mn + (*n << 1) + nb * (*n + 1), i__1 = f2cmax(i__1,i__2), 
 | |
| 	    i__2 = (mn << 1) + nb * *nrhs;
 | |
|     lwkopt = f2cmax(i__1,i__2);
 | |
|     z__1.r = (doublereal) lwkopt, z__1.i = 0.;
 | |
|     work[1].r = z__1.r, work[1].i = z__1.i;
 | |
|     lquery = *lwork == -1;
 | |
|     if (*m < 0) {
 | |
| 	*info = -1;
 | |
|     } else if (*n < 0) {
 | |
| 	*info = -2;
 | |
|     } else if (*nrhs < 0) {
 | |
| 	*info = -3;
 | |
|     } else if (*lda < f2cmax(1,*m)) {
 | |
| 	*info = -5;
 | |
|     } else /* if(complicated condition) */ {
 | |
| /* Computing MAX */
 | |
| 	i__1 = f2cmax(1,*m);
 | |
| 	if (*ldb < f2cmax(i__1,*n)) {
 | |
| 	    *info = -7;
 | |
| 	} else /* if(complicated condition) */ {
 | |
| /* Computing MAX */
 | |
| 	    i__1 = mn << 1, i__2 = *n + 1, i__1 = f2cmax(i__1,i__2), i__2 = mn + 
 | |
| 		    *nrhs;
 | |
| 	    if (*lwork < mn + f2cmax(i__1,i__2) && ! lquery) {
 | |
| 		*info = -12;
 | |
| 	    }
 | |
| 	}
 | |
|     }
 | |
| 
 | |
|     if (*info != 0) {
 | |
| 	i__1 = -(*info);
 | |
| 	xerbla_("ZGELSY", &i__1, (ftnlen)6);
 | |
| 	return;
 | |
|     } else if (lquery) {
 | |
| 	return;
 | |
|     }
 | |
| 
 | |
| /*     Quick return if possible */
 | |
| 
 | |
| /* Computing MIN */
 | |
|     i__1 = f2cmin(*m,*n);
 | |
|     if (f2cmin(i__1,*nrhs) == 0) {
 | |
| 	*rank = 0;
 | |
| 	return;
 | |
|     }
 | |
| 
 | |
| /*     Get machine parameters */
 | |
| 
 | |
|     smlnum = dlamch_("S") / dlamch_("P");
 | |
|     bignum = 1. / smlnum;
 | |
|     dlabad_(&smlnum, &bignum);
 | |
| 
 | |
| /*     Scale A, B if f2cmax entries outside range [SMLNUM,BIGNUM] */
 | |
| 
 | |
|     anrm = zlange_("M", m, n, &a[a_offset], lda, &rwork[1]);
 | |
|     iascl = 0;
 | |
|     if (anrm > 0. && anrm < smlnum) {
 | |
| 
 | |
| /*        Scale matrix norm up to SMLNUM */
 | |
| 
 | |
| 	zlascl_("G", &c__0, &c__0, &anrm, &smlnum, m, n, &a[a_offset], lda, 
 | |
| 		info);
 | |
| 	iascl = 1;
 | |
|     } else if (anrm > bignum) {
 | |
| 
 | |
| /*        Scale matrix norm down to BIGNUM */
 | |
| 
 | |
| 	zlascl_("G", &c__0, &c__0, &anrm, &bignum, m, n, &a[a_offset], lda, 
 | |
| 		info);
 | |
| 	iascl = 2;
 | |
|     } else if (anrm == 0.) {
 | |
| 
 | |
| /*        Matrix all zero. Return zero solution. */
 | |
| 
 | |
| 	i__1 = f2cmax(*m,*n);
 | |
| 	zlaset_("F", &i__1, nrhs, &c_b1, &c_b1, &b[b_offset], ldb);
 | |
| 	*rank = 0;
 | |
| 	goto L70;
 | |
|     }
 | |
| 
 | |
|     bnrm = zlange_("M", m, nrhs, &b[b_offset], ldb, &rwork[1]);
 | |
|     ibscl = 0;
 | |
|     if (bnrm > 0. && bnrm < smlnum) {
 | |
| 
 | |
| /*        Scale matrix norm up to SMLNUM */
 | |
| 
 | |
| 	zlascl_("G", &c__0, &c__0, &bnrm, &smlnum, m, nrhs, &b[b_offset], ldb,
 | |
| 		 info);
 | |
| 	ibscl = 1;
 | |
|     } else if (bnrm > bignum) {
 | |
| 
 | |
| /*        Scale matrix norm down to BIGNUM */
 | |
| 
 | |
| 	zlascl_("G", &c__0, &c__0, &bnrm, &bignum, m, nrhs, &b[b_offset], ldb,
 | |
| 		 info);
 | |
| 	ibscl = 2;
 | |
|     }
 | |
| 
 | |
| /*     Compute QR factorization with column pivoting of A: */
 | |
| /*        A * P = Q * R */
 | |
| 
 | |
|     i__1 = *lwork - mn;
 | |
|     zgeqp3_(m, n, &a[a_offset], lda, &jpvt[1], &work[1], &work[mn + 1], &i__1,
 | |
| 	     &rwork[1], info);
 | |
|     i__1 = mn + 1;
 | |
|     wsize = mn + work[i__1].r;
 | |
| 
 | |
| /*     complex workspace: MN+NB*(N+1). real workspace 2*N. */
 | |
| /*     Details of Householder rotations stored in WORK(1:MN). */
 | |
| 
 | |
| /*     Determine RANK using incremental condition estimation */
 | |
| 
 | |
|     i__1 = ismin;
 | |
|     work[i__1].r = 1., work[i__1].i = 0.;
 | |
|     i__1 = ismax;
 | |
|     work[i__1].r = 1., work[i__1].i = 0.;
 | |
|     smax = z_abs(&a[a_dim1 + 1]);
 | |
|     smin = smax;
 | |
|     if (z_abs(&a[a_dim1 + 1]) == 0.) {
 | |
| 	*rank = 0;
 | |
| 	i__1 = f2cmax(*m,*n);
 | |
| 	zlaset_("F", &i__1, nrhs, &c_b1, &c_b1, &b[b_offset], ldb);
 | |
| 	goto L70;
 | |
|     } else {
 | |
| 	*rank = 1;
 | |
|     }
 | |
| 
 | |
| L10:
 | |
|     if (*rank < mn) {
 | |
| 	i__ = *rank + 1;
 | |
| 	zlaic1_(&c__2, rank, &work[ismin], &smin, &a[i__ * a_dim1 + 1], &a[
 | |
| 		i__ + i__ * a_dim1], &sminpr, &s1, &c1);
 | |
| 	zlaic1_(&c__1, rank, &work[ismax], &smax, &a[i__ * a_dim1 + 1], &a[
 | |
| 		i__ + i__ * a_dim1], &smaxpr, &s2, &c2);
 | |
| 
 | |
| 	if (smaxpr * *rcond <= sminpr) {
 | |
| 	    i__1 = *rank;
 | |
| 	    for (i__ = 1; i__ <= i__1; ++i__) {
 | |
| 		i__2 = ismin + i__ - 1;
 | |
| 		i__3 = ismin + i__ - 1;
 | |
| 		z__1.r = s1.r * work[i__3].r - s1.i * work[i__3].i, z__1.i = 
 | |
| 			s1.r * work[i__3].i + s1.i * work[i__3].r;
 | |
| 		work[i__2].r = z__1.r, work[i__2].i = z__1.i;
 | |
| 		i__2 = ismax + i__ - 1;
 | |
| 		i__3 = ismax + i__ - 1;
 | |
| 		z__1.r = s2.r * work[i__3].r - s2.i * work[i__3].i, z__1.i = 
 | |
| 			s2.r * work[i__3].i + s2.i * work[i__3].r;
 | |
| 		work[i__2].r = z__1.r, work[i__2].i = z__1.i;
 | |
| /* L20: */
 | |
| 	    }
 | |
| 	    i__1 = ismin + *rank;
 | |
| 	    work[i__1].r = c1.r, work[i__1].i = c1.i;
 | |
| 	    i__1 = ismax + *rank;
 | |
| 	    work[i__1].r = c2.r, work[i__1].i = c2.i;
 | |
| 	    smin = sminpr;
 | |
| 	    smax = smaxpr;
 | |
| 	    ++(*rank);
 | |
| 	    goto L10;
 | |
| 	}
 | |
|     }
 | |
| 
 | |
| /*     complex workspace: 3*MN. */
 | |
| 
 | |
| /*     Logically partition R = [ R11 R12 ] */
 | |
| /*                             [  0  R22 ] */
 | |
| /*     where R11 = R(1:RANK,1:RANK) */
 | |
| 
 | |
| /*     [R11,R12] = [ T11, 0 ] * Y */
 | |
| 
 | |
|     if (*rank < *n) {
 | |
| 	i__1 = *lwork - (mn << 1);
 | |
| 	ztzrzf_(rank, n, &a[a_offset], lda, &work[mn + 1], &work[(mn << 1) + 
 | |
| 		1], &i__1, info);
 | |
|     }
 | |
| 
 | |
| /*     complex workspace: 2*MN. */
 | |
| /*     Details of Householder rotations stored in WORK(MN+1:2*MN) */
 | |
| 
 | |
| /*     B(1:M,1:NRHS) := Q**H * B(1:M,1:NRHS) */
 | |
| 
 | |
|     i__1 = *lwork - (mn << 1);
 | |
|     zunmqr_("Left", "Conjugate transpose", m, nrhs, &mn, &a[a_offset], lda, &
 | |
| 	    work[1], &b[b_offset], ldb, &work[(mn << 1) + 1], &i__1, info);
 | |
| /* Computing MAX */
 | |
|     i__1 = (mn << 1) + 1;
 | |
|     d__1 = wsize, d__2 = (mn << 1) + work[i__1].r;
 | |
|     wsize = f2cmax(d__1,d__2);
 | |
| 
 | |
| /*     complex workspace: 2*MN+NB*NRHS. */
 | |
| 
 | |
| /*     B(1:RANK,1:NRHS) := inv(T11) * B(1:RANK,1:NRHS) */
 | |
| 
 | |
|     ztrsm_("Left", "Upper", "No transpose", "Non-unit", rank, nrhs, &c_b2, &a[
 | |
| 	    a_offset], lda, &b[b_offset], ldb);
 | |
| 
 | |
|     i__1 = *nrhs;
 | |
|     for (j = 1; j <= i__1; ++j) {
 | |
| 	i__2 = *n;
 | |
| 	for (i__ = *rank + 1; i__ <= i__2; ++i__) {
 | |
| 	    i__3 = i__ + j * b_dim1;
 | |
| 	    b[i__3].r = 0., b[i__3].i = 0.;
 | |
| /* L30: */
 | |
| 	}
 | |
| /* L40: */
 | |
|     }
 | |
| 
 | |
| /*     B(1:N,1:NRHS) := Y**H * B(1:N,1:NRHS) */
 | |
| 
 | |
|     if (*rank < *n) {
 | |
| 	i__1 = *n - *rank;
 | |
| 	i__2 = *lwork - (mn << 1);
 | |
| 	zunmrz_("Left", "Conjugate transpose", n, nrhs, rank, &i__1, &a[
 | |
| 		a_offset], lda, &work[mn + 1], &b[b_offset], ldb, &work[(mn <<
 | |
| 		 1) + 1], &i__2, info);
 | |
|     }
 | |
| 
 | |
| /*     complex workspace: 2*MN+NRHS. */
 | |
| 
 | |
| /*     B(1:N,1:NRHS) := P * B(1:N,1:NRHS) */
 | |
| 
 | |
|     i__1 = *nrhs;
 | |
|     for (j = 1; j <= i__1; ++j) {
 | |
| 	i__2 = *n;
 | |
| 	for (i__ = 1; i__ <= i__2; ++i__) {
 | |
| 	    i__3 = jpvt[i__];
 | |
| 	    i__4 = i__ + j * b_dim1;
 | |
| 	    work[i__3].r = b[i__4].r, work[i__3].i = b[i__4].i;
 | |
| /* L50: */
 | |
| 	}
 | |
| 	zcopy_(n, &work[1], &c__1, &b[j * b_dim1 + 1], &c__1);
 | |
| /* L60: */
 | |
|     }
 | |
| 
 | |
| /*     complex workspace: N. */
 | |
| 
 | |
| /*     Undo scaling */
 | |
| 
 | |
|     if (iascl == 1) {
 | |
| 	zlascl_("G", &c__0, &c__0, &anrm, &smlnum, n, nrhs, &b[b_offset], ldb,
 | |
| 		 info);
 | |
| 	zlascl_("U", &c__0, &c__0, &smlnum, &anrm, rank, rank, &a[a_offset], 
 | |
| 		lda, info);
 | |
|     } else if (iascl == 2) {
 | |
| 	zlascl_("G", &c__0, &c__0, &anrm, &bignum, n, nrhs, &b[b_offset], ldb,
 | |
| 		 info);
 | |
| 	zlascl_("U", &c__0, &c__0, &bignum, &anrm, rank, rank, &a[a_offset], 
 | |
| 		lda, info);
 | |
|     }
 | |
|     if (ibscl == 1) {
 | |
| 	zlascl_("G", &c__0, &c__0, &smlnum, &bnrm, n, nrhs, &b[b_offset], ldb,
 | |
| 		 info);
 | |
|     } else if (ibscl == 2) {
 | |
| 	zlascl_("G", &c__0, &c__0, &bignum, &bnrm, n, nrhs, &b[b_offset], ldb,
 | |
| 		 info);
 | |
|     }
 | |
| 
 | |
| L70:
 | |
|     z__1.r = (doublereal) lwkopt, z__1.i = 0.;
 | |
|     work[1].r = z__1.r, work[1].i = z__1.i;
 | |
| 
 | |
|     return;
 | |
| 
 | |
| /*     End of ZGELSY */
 | |
| 
 | |
| } /* zgelsy_ */
 | |
| 
 |