1047 lines
30 KiB
C
1047 lines
30 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 blasint logical;
|
|
|
|
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++ */
|
|
|
|
|
|
#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 integer c__1 = 1;
|
|
static integer c__0 = 0;
|
|
static doublereal c_b13 = 1.;
|
|
static doublereal c_b26 = 0.;
|
|
|
|
/* > \brief \b DLASD3 finds all square roots of the roots of the secular equation, as defined by the values in
|
|
D and Z, and then updates the singular vectors by matrix multiplication. Used by sbdsdc. */
|
|
|
|
/* =========== DOCUMENTATION =========== */
|
|
|
|
/* Online html documentation available at */
|
|
/* http://www.netlib.org/lapack/explore-html/ */
|
|
|
|
/* > \htmlonly */
|
|
/* > Download DLASD3 + dependencies */
|
|
/* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlasd3.
|
|
f"> */
|
|
/* > [TGZ]</a> */
|
|
/* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlasd3.
|
|
f"> */
|
|
/* > [ZIP]</a> */
|
|
/* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlasd3.
|
|
f"> */
|
|
/* > [TXT]</a> */
|
|
/* > \endhtmlonly */
|
|
|
|
/* Definition: */
|
|
/* =========== */
|
|
|
|
/* SUBROUTINE DLASD3( NL, NR, SQRE, K, D, Q, LDQ, DSIGMA, U, LDU, U2, */
|
|
/* LDU2, VT, LDVT, VT2, LDVT2, IDXC, CTOT, Z, */
|
|
/* INFO ) */
|
|
|
|
/* INTEGER INFO, K, LDQ, LDU, LDU2, LDVT, LDVT2, NL, NR, */
|
|
/* $ SQRE */
|
|
/* INTEGER CTOT( * ), IDXC( * ) */
|
|
/* DOUBLE PRECISION D( * ), DSIGMA( * ), Q( LDQ, * ), U( LDU, * ), */
|
|
/* $ U2( LDU2, * ), VT( LDVT, * ), VT2( LDVT2, * ), */
|
|
/* $ Z( * ) */
|
|
|
|
|
|
/* > \par Purpose: */
|
|
/* ============= */
|
|
/* > */
|
|
/* > \verbatim */
|
|
/* > */
|
|
/* > DLASD3 finds all the square roots of the roots of the secular */
|
|
/* > equation, as defined by the values in D and Z. It makes the */
|
|
/* > appropriate calls to DLASD4 and then updates the singular */
|
|
/* > vectors by matrix multiplication. */
|
|
/* > */
|
|
/* > This code makes very mild assumptions about floating point */
|
|
/* > arithmetic. It will work on machines with a guard digit in */
|
|
/* > add/subtract, or on those binary machines without guard digits */
|
|
/* > which subtract like the Cray XMP, Cray YMP, Cray C 90, or Cray 2. */
|
|
/* > It could conceivably fail on hexadecimal or decimal machines */
|
|
/* > without guard digits, but we know of none. */
|
|
/* > */
|
|
/* > DLASD3 is called from DLASD1. */
|
|
/* > \endverbatim */
|
|
|
|
/* Arguments: */
|
|
/* ========== */
|
|
|
|
/* > \param[in] NL */
|
|
/* > \verbatim */
|
|
/* > NL is INTEGER */
|
|
/* > The row dimension of the upper block. NL >= 1. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] NR */
|
|
/* > \verbatim */
|
|
/* > NR is INTEGER */
|
|
/* > The row dimension of the lower block. NR >= 1. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] SQRE */
|
|
/* > \verbatim */
|
|
/* > SQRE is INTEGER */
|
|
/* > = 0: the lower block is an NR-by-NR square matrix. */
|
|
/* > = 1: the lower block is an NR-by-(NR+1) rectangular matrix. */
|
|
/* > */
|
|
/* > The bidiagonal matrix has N = NL + NR + 1 rows and */
|
|
/* > M = N + SQRE >= N columns. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] K */
|
|
/* > \verbatim */
|
|
/* > K is INTEGER */
|
|
/* > The size of the secular equation, 1 =< K = < N. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[out] D */
|
|
/* > \verbatim */
|
|
/* > D is DOUBLE PRECISION array, dimension(K) */
|
|
/* > On exit the square roots of the roots of the secular equation, */
|
|
/* > in ascending order. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[out] Q */
|
|
/* > \verbatim */
|
|
/* > Q is DOUBLE PRECISION array, dimension (LDQ,K) */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] LDQ */
|
|
/* > \verbatim */
|
|
/* > LDQ is INTEGER */
|
|
/* > The leading dimension of the array Q. LDQ >= K. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in,out] DSIGMA */
|
|
/* > \verbatim */
|
|
/* > DSIGMA is DOUBLE PRECISION array, dimension(K) */
|
|
/* > The first K elements of this array contain the old roots */
|
|
/* > of the deflated updating problem. These are the poles */
|
|
/* > of the secular equation. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[out] U */
|
|
/* > \verbatim */
|
|
/* > U is DOUBLE PRECISION array, dimension (LDU, N) */
|
|
/* > The last N - K columns of this matrix contain the deflated */
|
|
/* > left singular vectors. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] LDU */
|
|
/* > \verbatim */
|
|
/* > LDU is INTEGER */
|
|
/* > The leading dimension of the array U. LDU >= N. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] U2 */
|
|
/* > \verbatim */
|
|
/* > U2 is DOUBLE PRECISION array, dimension (LDU2, N) */
|
|
/* > The first K columns of this matrix contain the non-deflated */
|
|
/* > left singular vectors for the split problem. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] LDU2 */
|
|
/* > \verbatim */
|
|
/* > LDU2 is INTEGER */
|
|
/* > The leading dimension of the array U2. LDU2 >= N. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[out] VT */
|
|
/* > \verbatim */
|
|
/* > VT is DOUBLE PRECISION array, dimension (LDVT, M) */
|
|
/* > The last M - K columns of VT**T contain the deflated */
|
|
/* > right singular vectors. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] LDVT */
|
|
/* > \verbatim */
|
|
/* > LDVT is INTEGER */
|
|
/* > The leading dimension of the array VT. LDVT >= N. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in,out] VT2 */
|
|
/* > \verbatim */
|
|
/* > VT2 is DOUBLE PRECISION array, dimension (LDVT2, N) */
|
|
/* > The first K columns of VT2**T contain the non-deflated */
|
|
/* > right singular vectors for the split problem. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] LDVT2 */
|
|
/* > \verbatim */
|
|
/* > LDVT2 is INTEGER */
|
|
/* > The leading dimension of the array VT2. LDVT2 >= N. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] IDXC */
|
|
/* > \verbatim */
|
|
/* > IDXC is INTEGER array, dimension ( N ) */
|
|
/* > The permutation used to arrange the columns of U (and rows of */
|
|
/* > VT) into three groups: the first group contains non-zero */
|
|
/* > entries only at and above (or before) NL +1; the second */
|
|
/* > contains non-zero entries only at and below (or after) NL+2; */
|
|
/* > and the third is dense. The first column of U and the row of */
|
|
/* > VT are treated separately, however. */
|
|
/* > */
|
|
/* > The rows of the singular vectors found by DLASD4 */
|
|
/* > must be likewise permuted before the matrix multiplies can */
|
|
/* > take place. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in] CTOT */
|
|
/* > \verbatim */
|
|
/* > CTOT is INTEGER array, dimension ( 4 ) */
|
|
/* > A count of the total number of the various types of columns */
|
|
/* > in U (or rows in VT), as described in IDXC. The fourth column */
|
|
/* > type is any column which has been deflated. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[in,out] Z */
|
|
/* > \verbatim */
|
|
/* > Z is DOUBLE PRECISION array, dimension (K) */
|
|
/* > The first K elements of this array contain the components */
|
|
/* > of the deflation-adjusted updating row vector. */
|
|
/* > \endverbatim */
|
|
/* > */
|
|
/* > \param[out] INFO */
|
|
/* > \verbatim */
|
|
/* > INFO is INTEGER */
|
|
/* > = 0: successful exit. */
|
|
/* > < 0: if INFO = -i, the i-th argument had an illegal value. */
|
|
/* > > 0: if INFO = 1, a singular value did not converge */
|
|
/* > \endverbatim */
|
|
|
|
/* Authors: */
|
|
/* ======== */
|
|
|
|
/* > \author Univ. of Tennessee */
|
|
/* > \author Univ. of California Berkeley */
|
|
/* > \author Univ. of Colorado Denver */
|
|
/* > \author NAG Ltd. */
|
|
|
|
/* > \date June 2017 */
|
|
|
|
/* > \ingroup OTHERauxiliary */
|
|
|
|
/* > \par Contributors: */
|
|
/* ================== */
|
|
/* > */
|
|
/* > Ming Gu and Huan Ren, Computer Science Division, University of */
|
|
/* > California at Berkeley, USA */
|
|
/* > */
|
|
/* ===================================================================== */
|
|
/* Subroutine */ void dlasd3_(integer *nl, integer *nr, integer *sqre, integer
|
|
*k, doublereal *d__, doublereal *q, integer *ldq, doublereal *dsigma,
|
|
doublereal *u, integer *ldu, doublereal *u2, integer *ldu2,
|
|
doublereal *vt, integer *ldvt, doublereal *vt2, integer *ldvt2,
|
|
integer *idxc, integer *ctot, doublereal *z__, integer *info)
|
|
{
|
|
/* System generated locals */
|
|
integer q_dim1, q_offset, u_dim1, u_offset, u2_dim1, u2_offset, vt_dim1,
|
|
vt_offset, vt2_dim1, vt2_offset, i__1, i__2;
|
|
doublereal d__1, d__2;
|
|
|
|
/* Local variables */
|
|
doublereal temp;
|
|
extern doublereal dnrm2_(integer *, doublereal *, integer *);
|
|
integer i__, j, m, n;
|
|
extern /* Subroutine */ void dgemm_(char *, char *, integer *, integer *,
|
|
integer *, doublereal *, doublereal *, integer *, doublereal *,
|
|
integer *, doublereal *, doublereal *, integer *);
|
|
integer ctemp;
|
|
extern /* Subroutine */ void dcopy_(integer *, doublereal *, integer *,
|
|
doublereal *, integer *);
|
|
integer ktemp;
|
|
extern doublereal dlamc3_(doublereal *, doublereal *);
|
|
extern /* Subroutine */ void dlasd4_(integer *, integer *, doublereal *,
|
|
doublereal *, doublereal *, doublereal *, doublereal *,
|
|
doublereal *, integer *);
|
|
integer jc;
|
|
extern /* Subroutine */ void dlascl_(char *, integer *, integer *,
|
|
doublereal *, doublereal *, integer *, integer *, doublereal *,
|
|
integer *, integer *), dlacpy_(char *, integer *, integer
|
|
*, doublereal *, integer *, doublereal *, integer *);
|
|
extern int xerbla_(char *, integer *, ftnlen);
|
|
doublereal rho;
|
|
integer nlp1, nlp2, nrp1;
|
|
|
|
|
|
/* -- LAPACK auxiliary routine (version 3.7.1) -- */
|
|
/* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
|
|
/* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
|
|
/* June 2017 */
|
|
|
|
|
|
/* ===================================================================== */
|
|
|
|
|
|
/* Test the input parameters. */
|
|
|
|
/* Parameter adjustments */
|
|
--d__;
|
|
q_dim1 = *ldq;
|
|
q_offset = 1 + q_dim1 * 1;
|
|
q -= q_offset;
|
|
--dsigma;
|
|
u_dim1 = *ldu;
|
|
u_offset = 1 + u_dim1 * 1;
|
|
u -= u_offset;
|
|
u2_dim1 = *ldu2;
|
|
u2_offset = 1 + u2_dim1 * 1;
|
|
u2 -= u2_offset;
|
|
vt_dim1 = *ldvt;
|
|
vt_offset = 1 + vt_dim1 * 1;
|
|
vt -= vt_offset;
|
|
vt2_dim1 = *ldvt2;
|
|
vt2_offset = 1 + vt2_dim1 * 1;
|
|
vt2 -= vt2_offset;
|
|
--idxc;
|
|
--ctot;
|
|
--z__;
|
|
|
|
/* Function Body */
|
|
*info = 0;
|
|
|
|
if (*nl < 1) {
|
|
*info = -1;
|
|
} else if (*nr < 1) {
|
|
*info = -2;
|
|
} else if (*sqre != 1 && *sqre != 0) {
|
|
*info = -3;
|
|
}
|
|
|
|
n = *nl + *nr + 1;
|
|
m = n + *sqre;
|
|
nlp1 = *nl + 1;
|
|
nlp2 = *nl + 2;
|
|
|
|
if (*k < 1 || *k > n) {
|
|
*info = -4;
|
|
} else if (*ldq < *k) {
|
|
*info = -7;
|
|
} else if (*ldu < n) {
|
|
*info = -10;
|
|
} else if (*ldu2 < n) {
|
|
*info = -12;
|
|
} else if (*ldvt < m) {
|
|
*info = -14;
|
|
} else if (*ldvt2 < m) {
|
|
*info = -16;
|
|
}
|
|
if (*info != 0) {
|
|
i__1 = -(*info);
|
|
xerbla_("DLASD3", &i__1, (ftnlen)6);
|
|
return;
|
|
}
|
|
|
|
/* Quick return if possible */
|
|
|
|
if (*k == 1) {
|
|
d__[1] = abs(z__[1]);
|
|
dcopy_(&m, &vt2[vt2_dim1 + 1], ldvt2, &vt[vt_dim1 + 1], ldvt);
|
|
if (z__[1] > 0.) {
|
|
dcopy_(&n, &u2[u2_dim1 + 1], &c__1, &u[u_dim1 + 1], &c__1);
|
|
} else {
|
|
i__1 = n;
|
|
for (i__ = 1; i__ <= i__1; ++i__) {
|
|
u[i__ + u_dim1] = -u2[i__ + u2_dim1];
|
|
/* L10: */
|
|
}
|
|
}
|
|
return;
|
|
}
|
|
|
|
/* Modify values DSIGMA(i) to make sure all DSIGMA(i)-DSIGMA(j) can */
|
|
/* be computed with high relative accuracy (barring over/underflow). */
|
|
/* This is a problem on machines without a guard digit in */
|
|
/* add/subtract (Cray XMP, Cray YMP, Cray C 90 and Cray 2). */
|
|
/* The following code replaces DSIGMA(I) by 2*DSIGMA(I)-DSIGMA(I), */
|
|
/* which on any of these machines zeros out the bottommost */
|
|
/* bit of DSIGMA(I) if it is 1; this makes the subsequent */
|
|
/* subtractions DSIGMA(I)-DSIGMA(J) unproblematic when cancellation */
|
|
/* occurs. On binary machines with a guard digit (almost all */
|
|
/* machines) it does not change DSIGMA(I) at all. On hexadecimal */
|
|
/* and decimal machines with a guard digit, it slightly */
|
|
/* changes the bottommost bits of DSIGMA(I). It does not account */
|
|
/* for hexadecimal or decimal machines without guard digits */
|
|
/* (we know of none). We use a subroutine call to compute */
|
|
/* 2*DSIGMA(I) to prevent optimizing compilers from eliminating */
|
|
/* this code. */
|
|
|
|
i__1 = *k;
|
|
for (i__ = 1; i__ <= i__1; ++i__) {
|
|
dsigma[i__] = dlamc3_(&dsigma[i__], &dsigma[i__]) - dsigma[i__];
|
|
/* L20: */
|
|
}
|
|
|
|
/* Keep a copy of Z. */
|
|
|
|
dcopy_(k, &z__[1], &c__1, &q[q_offset], &c__1);
|
|
|
|
/* Normalize Z. */
|
|
|
|
rho = dnrm2_(k, &z__[1], &c__1);
|
|
dlascl_("G", &c__0, &c__0, &rho, &c_b13, k, &c__1, &z__[1], k, info);
|
|
rho *= rho;
|
|
|
|
/* Find the new singular values. */
|
|
|
|
i__1 = *k;
|
|
for (j = 1; j <= i__1; ++j) {
|
|
dlasd4_(k, &j, &dsigma[1], &z__[1], &u[j * u_dim1 + 1], &rho, &d__[j],
|
|
&vt[j * vt_dim1 + 1], info);
|
|
|
|
/* If the zero finder fails, report the convergence failure. */
|
|
|
|
if (*info != 0) {
|
|
return;
|
|
}
|
|
/* L30: */
|
|
}
|
|
|
|
/* Compute updated Z. */
|
|
|
|
i__1 = *k;
|
|
for (i__ = 1; i__ <= i__1; ++i__) {
|
|
z__[i__] = u[i__ + *k * u_dim1] * vt[i__ + *k * vt_dim1];
|
|
i__2 = i__ - 1;
|
|
for (j = 1; j <= i__2; ++j) {
|
|
z__[i__] *= u[i__ + j * u_dim1] * vt[i__ + j * vt_dim1] / (dsigma[
|
|
i__] - dsigma[j]) / (dsigma[i__] + dsigma[j]);
|
|
/* L40: */
|
|
}
|
|
i__2 = *k - 1;
|
|
for (j = i__; j <= i__2; ++j) {
|
|
z__[i__] *= u[i__ + j * u_dim1] * vt[i__ + j * vt_dim1] / (dsigma[
|
|
i__] - dsigma[j + 1]) / (dsigma[i__] + dsigma[j + 1]);
|
|
/* L50: */
|
|
}
|
|
d__2 = sqrt((d__1 = z__[i__], abs(d__1)));
|
|
z__[i__] = d_sign(&d__2, &q[i__ + q_dim1]);
|
|
/* L60: */
|
|
}
|
|
|
|
/* Compute left singular vectors of the modified diagonal matrix, */
|
|
/* and store related information for the right singular vectors. */
|
|
|
|
i__1 = *k;
|
|
for (i__ = 1; i__ <= i__1; ++i__) {
|
|
vt[i__ * vt_dim1 + 1] = z__[1] / u[i__ * u_dim1 + 1] / vt[i__ *
|
|
vt_dim1 + 1];
|
|
u[i__ * u_dim1 + 1] = -1.;
|
|
i__2 = *k;
|
|
for (j = 2; j <= i__2; ++j) {
|
|
vt[j + i__ * vt_dim1] = z__[j] / u[j + i__ * u_dim1] / vt[j + i__
|
|
* vt_dim1];
|
|
u[j + i__ * u_dim1] = dsigma[j] * vt[j + i__ * vt_dim1];
|
|
/* L70: */
|
|
}
|
|
temp = dnrm2_(k, &u[i__ * u_dim1 + 1], &c__1);
|
|
q[i__ * q_dim1 + 1] = u[i__ * u_dim1 + 1] / temp;
|
|
i__2 = *k;
|
|
for (j = 2; j <= i__2; ++j) {
|
|
jc = idxc[j];
|
|
q[j + i__ * q_dim1] = u[jc + i__ * u_dim1] / temp;
|
|
/* L80: */
|
|
}
|
|
/* L90: */
|
|
}
|
|
|
|
/* Update the left singular vector matrix. */
|
|
|
|
if (*k == 2) {
|
|
dgemm_("N", "N", &n, k, k, &c_b13, &u2[u2_offset], ldu2, &q[q_offset],
|
|
ldq, &c_b26, &u[u_offset], ldu);
|
|
goto L100;
|
|
}
|
|
if (ctot[1] > 0) {
|
|
dgemm_("N", "N", nl, k, &ctot[1], &c_b13, &u2[(u2_dim1 << 1) + 1],
|
|
ldu2, &q[q_dim1 + 2], ldq, &c_b26, &u[u_dim1 + 1], ldu);
|
|
if (ctot[3] > 0) {
|
|
ktemp = ctot[1] + 2 + ctot[2];
|
|
dgemm_("N", "N", nl, k, &ctot[3], &c_b13, &u2[ktemp * u2_dim1 + 1]
|
|
, ldu2, &q[ktemp + q_dim1], ldq, &c_b13, &u[u_dim1 + 1],
|
|
ldu);
|
|
}
|
|
} else if (ctot[3] > 0) {
|
|
ktemp = ctot[1] + 2 + ctot[2];
|
|
dgemm_("N", "N", nl, k, &ctot[3], &c_b13, &u2[ktemp * u2_dim1 + 1],
|
|
ldu2, &q[ktemp + q_dim1], ldq, &c_b26, &u[u_dim1 + 1], ldu);
|
|
} else {
|
|
dlacpy_("F", nl, k, &u2[u2_offset], ldu2, &u[u_offset], ldu);
|
|
}
|
|
dcopy_(k, &q[q_dim1 + 1], ldq, &u[nlp1 + u_dim1], ldu);
|
|
ktemp = ctot[1] + 2;
|
|
ctemp = ctot[2] + ctot[3];
|
|
dgemm_("N", "N", nr, k, &ctemp, &c_b13, &u2[nlp2 + ktemp * u2_dim1], ldu2,
|
|
&q[ktemp + q_dim1], ldq, &c_b26, &u[nlp2 + u_dim1], ldu);
|
|
|
|
/* Generate the right singular vectors. */
|
|
|
|
L100:
|
|
i__1 = *k;
|
|
for (i__ = 1; i__ <= i__1; ++i__) {
|
|
temp = dnrm2_(k, &vt[i__ * vt_dim1 + 1], &c__1);
|
|
q[i__ + q_dim1] = vt[i__ * vt_dim1 + 1] / temp;
|
|
i__2 = *k;
|
|
for (j = 2; j <= i__2; ++j) {
|
|
jc = idxc[j];
|
|
q[i__ + j * q_dim1] = vt[jc + i__ * vt_dim1] / temp;
|
|
/* L110: */
|
|
}
|
|
/* L120: */
|
|
}
|
|
|
|
/* Update the right singular vector matrix. */
|
|
|
|
if (*k == 2) {
|
|
dgemm_("N", "N", k, &m, k, &c_b13, &q[q_offset], ldq, &vt2[vt2_offset]
|
|
, ldvt2, &c_b26, &vt[vt_offset], ldvt);
|
|
return;
|
|
}
|
|
ktemp = ctot[1] + 1;
|
|
dgemm_("N", "N", k, &nlp1, &ktemp, &c_b13, &q[q_dim1 + 1], ldq, &vt2[
|
|
vt2_dim1 + 1], ldvt2, &c_b26, &vt[vt_dim1 + 1], ldvt);
|
|
ktemp = ctot[1] + 2 + ctot[2];
|
|
if (ktemp <= *ldvt2) {
|
|
dgemm_("N", "N", k, &nlp1, &ctot[3], &c_b13, &q[ktemp * q_dim1 + 1],
|
|
ldq, &vt2[ktemp + vt2_dim1], ldvt2, &c_b13, &vt[vt_dim1 + 1],
|
|
ldvt);
|
|
}
|
|
|
|
ktemp = ctot[1] + 1;
|
|
nrp1 = *nr + *sqre;
|
|
if (ktemp > 1) {
|
|
i__1 = *k;
|
|
for (i__ = 1; i__ <= i__1; ++i__) {
|
|
q[i__ + ktemp * q_dim1] = q[i__ + q_dim1];
|
|
/* L130: */
|
|
}
|
|
i__1 = m;
|
|
for (i__ = nlp2; i__ <= i__1; ++i__) {
|
|
vt2[ktemp + i__ * vt2_dim1] = vt2[i__ * vt2_dim1 + 1];
|
|
/* L140: */
|
|
}
|
|
}
|
|
ctemp = ctot[2] + 1 + ctot[3];
|
|
dgemm_("N", "N", k, &nrp1, &ctemp, &c_b13, &q[ktemp * q_dim1 + 1], ldq, &
|
|
vt2[ktemp + nlp2 * vt2_dim1], ldvt2, &c_b26, &vt[nlp2 * vt_dim1 +
|
|
1], ldvt);
|
|
|
|
return;
|
|
|
|
/* End of DLASD3 */
|
|
|
|
} /* dlasd3_ */
|
|
|