#include #include #include #include #include #include #include typedef double complex TCD; typedef float complex TCF; #undef complex #include "lapack-aux.h" #define MACRO(B) do {B} while (0) #define ERROR(CODE) MACRO(return CODE;) #define REQUIRES(COND, CODE) MACRO(if(!(COND)) {ERROR(CODE);}) #define MIN(A,B) ((A)<(B)?(A):(B)) #define MAX(A,B) ((A)>(B)?(A):(B)) // #define DBGL #ifdef DBGL #define DEBUGMSG(M) printf("\nLAPACK "M"\n"); #else #define DEBUGMSG(M) #endif #define OK return 0; // #ifdef DBGL // #define DEBUGMSG(M) printf("LAPACK Wrapper "M"\n: "); size_t t0 = time(NULL); // #define OK MACRO(printf("%ld s\n",time(0)-t0); return 0;); // #else // #define DEBUGMSG(M) // #define OK return 0; // #endif #define INFOMAT(M) printf("%dx%d %d:%d\n",M##r,M##c,M##Xr,M##Xc); #define TRACEMAT(M) {int q; printf(" %d x %d: ",M##r,M##c); \ for(q=0;q=1 && ar==ac && ar==br,BAD_SIZE); DEBUGMSG("linearSolveR_l"); integer * ipiv = (integer*)malloc(n*sizeof(integer)); integer res; dgesv_ (&n,&nhrs, ap, &n, ipiv, bp, &n, &res); if(res>0) { return SINGULAR; } CHECK(res,res); free(ipiv); OK } //////////////////// general complex linear system //////////// int zgesv_(integer *n, integer *nrhs, doublecomplex *a, integer *lda, integer *ipiv, doublecomplex *b, integer *ldb, integer * info); int linearSolveC_l(OCMAT(a),OCMAT(b)) { integer n = ar; integer nhrs = bc; REQUIRES(n>=1 && ar==ac && ar==br,BAD_SIZE); DEBUGMSG("linearSolveC_l"); integer * ipiv = (integer*)malloc(n*sizeof(integer)); integer res; zgesv_ (&n,&nhrs, ap, &n, ipiv, bp, &n, &res); if(res>0) { return SINGULAR; } CHECK(res,res); free(ipiv); OK } //////// symmetric positive definite real linear system using Cholesky //////////// int dpotrs_(char *uplo, integer *n, integer *nrhs, doublereal *a, integer *lda, doublereal *b, integer *ldb, integer * info); int cholSolveR_l(KODMAT(a),ODMAT(b)) { integer n = ar; integer lda = aXc; integer nhrs = bc; REQUIRES(n>=1 && ar==ac && ar==br,BAD_SIZE); DEBUGMSG("cholSolveR_l"); integer res; dpotrs_ ("U", &n,&nhrs, (double*)ap, &lda, bp, &n, &res); CHECK(res,res); OK } //////// Hermitian positive definite real linear system using Cholesky //////////// int zpotrs_(char *uplo, integer *n, integer *nrhs, doublecomplex *a, integer *lda, doublecomplex *b, integer *ldb, integer *info); int cholSolveC_l(KOCMAT(a),OCMAT(b)) { integer n = ar; integer lda = aXc; integer nhrs = bc; REQUIRES(n>=1 && ar==ac && ar==br,BAD_SIZE); DEBUGMSG("cholSolveC_l"); integer res; zpotrs_ ("U", &n,&nhrs, (doublecomplex*)ap, &lda, bp, &n, &res); CHECK(res,res); OK } //////////////////// least squares real linear system //////////// int dgels_(char *trans, integer *m, integer *n, integer * nrhs, doublereal *a, integer *lda, doublereal *b, integer *ldb, doublereal *work, integer *lwork, integer *info); int linearSolveLSR_l(ODMAT(a),ODMAT(b)) { integer m = ar; integer n = ac; integer nrhs = bc; integer ldb = bXc; REQUIRES(m>=1 && n>=1 && br==MAX(m,n), BAD_SIZE); DEBUGMSG("linearSolveLSR_l"); integer res; integer lwork = -1; double ans; dgels_ ("N",&m,&n,&nrhs, ap,&m, bp,&ldb, &ans,&lwork, &res); lwork = ceil(ans); double * work = (double*)malloc(lwork*sizeof(double)); dgels_ ("N",&m,&n,&nrhs, ap,&m, bp,&ldb, work,&lwork, &res); if(res>0) { return SINGULAR; } CHECK(res,res); free(work); OK } //////////////////// least squares complex linear system //////////// int zgels_(char *trans, integer *m, integer *n, integer * nrhs, doublecomplex *a, integer *lda, doublecomplex *b, integer *ldb, doublecomplex *work, integer *lwork, integer *info); int linearSolveLSC_l(OCMAT(a),OCMAT(b)) { integer m = ar; integer n = ac; integer nrhs = bc; integer ldb = bXc; REQUIRES(m>=1 && n>=1 && br==MAX(m,n), BAD_SIZE); DEBUGMSG("linearSolveLSC_l"); integer res; integer lwork = -1; doublecomplex ans; zgels_ ("N",&m,&n,&nrhs, ap,&m, bp,&ldb, &ans,&lwork, &res); lwork = ceil(ans.r); doublecomplex * work = (doublecomplex*)malloc(lwork*sizeof(doublecomplex)); zgels_ ("N",&m,&n,&nrhs, ap,&m, bp,&ldb, work,&lwork, &res); if(res>0) { return SINGULAR; } CHECK(res,res); free(work); OK } //////////////////// least squares real linear system using SVD //////////// int dgelss_(integer *m, integer *n, integer *nrhs, doublereal *a, integer *lda, doublereal *b, integer *ldb, doublereal * s, doublereal *rcond, integer *rank, doublereal *work, integer *lwork, integer *info); int linearSolveSVDR_l(double rcond,ODMAT(a),ODMAT(b)) { integer m = ar; integer n = ac; integer nrhs = bc; integer ldb = bXc; REQUIRES(m>=1 && n>=1 && br==MAX(m,n), BAD_SIZE); DEBUGMSG("linearSolveSVDR_l"); double*S = (double*)malloc(MIN(m,n)*sizeof(double)); integer res; integer lwork = -1; integer rank; double ans; dgelss_ (&m,&n,&nrhs, ap,&m, bp,&ldb, S, &rcond,&rank, &ans,&lwork, &res); lwork = ceil(ans); double * work = (double*)malloc(lwork*sizeof(double)); dgelss_ (&m,&n,&nrhs, ap,&m, bp,&ldb, S, &rcond,&rank, work,&lwork, &res); if(res>0) { return NOCONVER; } CHECK(res,res); free(work); free(S); OK } //////////////////// least squares complex linear system using SVD //////////// int zgelss_(integer *m, integer *n, integer *nhrs, doublecomplex *a, integer *lda, doublecomplex *b, integer *ldb, doublereal *s, doublereal *rcond, integer* rank, doublecomplex *work, integer* lwork, doublereal* rwork, integer *info); int linearSolveSVDC_l(double rcond, OCMAT(a),OCMAT(b)) { integer m = ar; integer n = ac; integer nrhs = bc; integer ldb = bXc; REQUIRES(m>=1 && n>=1 && br==MAX(m,n), BAD_SIZE); DEBUGMSG("linearSolveSVDC_l"); double*S = (double*)malloc(MIN(m,n)*sizeof(double)); double*RWORK = (double*)malloc(5*MIN(m,n)*sizeof(double)); integer res; integer lwork = -1; integer rank; doublecomplex ans; zgelss_ (&m,&n,&nrhs, ap,&m, bp,&ldb, S, &rcond,&rank, &ans,&lwork, RWORK, &res); lwork = ceil(ans.r); doublecomplex * work = (doublecomplex*)malloc(lwork*sizeof(doublecomplex)); zgelss_ (&m,&n,&nrhs, ap,&m, bp,&ldb, S, &rcond,&rank, work,&lwork, RWORK, &res); if(res>0) { return NOCONVER; } CHECK(res,res); free(work); free(RWORK); free(S); OK } //////////////////// Cholesky factorization ///////////////////////// int zpotrf_(char *uplo, integer *n, doublecomplex *a, integer *lda, integer *info); int chol_l_H(OCMAT(l)) { integer n = lr; REQUIRES(n>=1 && lc == n,BAD_SIZE); DEBUGMSG("chol_l_H"); char uplo = 'U'; integer res; zpotrf_ (&uplo,&n,lp,&n,&res); CHECK(res>0,NODEFPOS); CHECK(res,res); doublecomplex zero = {0.,0.}; int r,c; for (r=0; r=1 && lc == n,BAD_SIZE); DEBUGMSG("chol_l_S"); char uplo = 'U'; integer res; dpotrf_ (&uplo,&n,lp,&n,&res); CHECK(res>0,NODEFPOS); CHECK(res,res); int r,c; for (r=0; r=1 && n >=1 && taun == mn, BAD_SIZE); DEBUGMSG("qr_l_R"); double *WORK = (double*)malloc(n*sizeof(double)); CHECK(!WORK,MEM); integer res; dgeqr2_ (&m,&n,rp,&m,taup,WORK,&res); CHECK(res,res); free(WORK); OK } int zgeqr2_(integer *m, integer *n, doublecomplex *a, integer *lda, doublecomplex *tau, doublecomplex *work, integer *info); int qr_l_C(CVEC(tau), OCMAT(r)) { integer m = rr; integer n = rc; integer mn = MIN(m,n); REQUIRES(m>=1 && n >=1 && taun == mn, BAD_SIZE); DEBUGMSG("qr_l_C"); doublecomplex *WORK = (doublecomplex*)malloc(n*sizeof(doublecomplex)); CHECK(!WORK,MEM); integer res; zgeqr2_ (&m,&n,rp,&m,taup,WORK,&res); CHECK(res,res); free(WORK); OK } int dorgqr_(integer *m, integer *n, integer *k, doublereal * a, integer *lda, doublereal *tau, doublereal *work, integer *lwork, integer *info); int c_dorgqr(KDVEC(tau), ODMAT(r)) { integer m = rr; integer n = MIN(rc,rr); integer k = taun; DEBUGMSG("c_dorgqr"); integer lwork = 8*n; // FIXME double *WORK = (double*)malloc(lwork*sizeof(double)); CHECK(!WORK,MEM); integer res; dorgqr_ (&m,&n,&k,rp,&m,(double*)taup,WORK,&lwork,&res); CHECK(res,res); free(WORK); OK } int zungqr_(integer *m, integer *n, integer *k, doublecomplex *a, integer *lda, doublecomplex *tau, doublecomplex * work, integer *lwork, integer *info); int c_zungqr(KCVEC(tau), OCMAT(r)) { integer m = rr; integer n = MIN(rc,rr); integer k = taun; DEBUGMSG("z_ungqr"); integer lwork = 8*n; // FIXME doublecomplex *WORK = (doublecomplex*)malloc(lwork*sizeof(doublecomplex)); CHECK(!WORK,MEM); integer res; zungqr_ (&m,&n,&k,rp,&m,(doublecomplex*)taup,WORK,&lwork,&res); CHECK(res,res); free(WORK); OK } //////////////////// Hessenberg factorization ///////////////////////// int dgehrd_(integer *n, integer *ilo, integer *ihi, doublereal *a, integer *lda, doublereal *tau, doublereal *work, integer *lwork, integer *info); int hess_l_R(DVEC(tau), ODMAT(r)) { integer m = rr; integer n = rc; integer mn = MIN(m,n); REQUIRES(m>=1 && n == m && taun == mn-1, BAD_SIZE); DEBUGMSG("hess_l_R"); integer lwork = 5*n; // FIXME double *WORK = (double*)malloc(lwork*sizeof(double)); CHECK(!WORK,MEM); integer res; integer one = 1; dgehrd_ (&n,&one,&n,rp,&n,taup,WORK,&lwork,&res); CHECK(res,res); free(WORK); OK } int zgehrd_(integer *n, integer *ilo, integer *ihi, doublecomplex *a, integer *lda, doublecomplex *tau, doublecomplex * work, integer *lwork, integer *info); int hess_l_C(CVEC(tau), OCMAT(r)) { integer m = rr; integer n = rc; integer mn = MIN(m,n); REQUIRES(m>=1 && n == m && taun == mn-1, BAD_SIZE); DEBUGMSG("hess_l_C"); integer lwork = 5*n; // FIXME doublecomplex *WORK = (doublecomplex*)malloc(lwork*sizeof(doublecomplex)); CHECK(!WORK,MEM); integer res; integer one = 1; zgehrd_ (&n,&one,&n,rp,&n,taup,WORK,&lwork,&res); CHECK(res,res); free(WORK); OK } //////////////////// Schur factorization ///////////////////////// int dgees_(char *jobvs, char *sort, L_fp select, integer *n, doublereal *a, integer *lda, integer *sdim, doublereal *wr, doublereal *wi, doublereal *vs, integer *ldvs, doublereal *work, integer *lwork, logical *bwork, integer *info); int schur_l_R(ODMAT(u), ODMAT(s)) { integer m = sr; integer n = sc; REQUIRES(m>=1 && n==m && ur==n && uc==n, BAD_SIZE); DEBUGMSG("schur_l_R"); integer lwork = 6*n; // FIXME double *WORK = (double*)malloc(lwork*sizeof(double)); double *WR = (double*)malloc(n*sizeof(double)); double *WI = (double*)malloc(n*sizeof(double)); // WR and WI not really required in this call logical *BWORK = (logical*)malloc(n*sizeof(logical)); integer res; integer sdim; dgees_ ("V","N",NULL,&n,sp,&n,&sdim,WR,WI,up,&n,WORK,&lwork,BWORK,&res); if(res>0) { return NOCONVER; } CHECK(res,res); free(WR); free(WI); free(BWORK); free(WORK); OK } int zgees_(char *jobvs, char *sort, L_fp select, integer *n, doublecomplex *a, integer *lda, integer *sdim, doublecomplex *w, doublecomplex *vs, integer *ldvs, doublecomplex *work, integer *lwork, doublereal *rwork, logical *bwork, integer *info); int schur_l_C(OCMAT(u), OCMAT(s)) { integer m = sr; integer n = sc; REQUIRES(m>=1 && n==m && ur==n && uc==n, BAD_SIZE); DEBUGMSG("schur_l_C"); integer lwork = 6*n; // FIXME doublecomplex *WORK = (doublecomplex*)malloc(lwork*sizeof(doublecomplex)); doublecomplex *W = (doublecomplex*)malloc(n*sizeof(doublecomplex)); // W not really required in this call logical *BWORK = (logical*)malloc(n*sizeof(logical)); double *RWORK = (double*)malloc(n*sizeof(double)); integer res; integer sdim; zgees_ ("V","N",NULL,&n,sp,&n,&sdim,W, up,&n, WORK,&lwork,RWORK,BWORK,&res); if(res>0) { return NOCONVER; } CHECK(res,res); free(W); free(BWORK); free(WORK); OK } //////////////////// LU factorization ///////////////////////// int dgetrf_(integer *m, integer *n, doublereal *a, integer * lda, integer *ipiv, integer *info); int lu_l_R(DVEC(ipiv), ODMAT(r)) { integer m = rr; integer n = rc; integer mn = MIN(m,n); REQUIRES(m>=1 && n >=1 && ipivn == mn, BAD_SIZE); DEBUGMSG("lu_l_R"); integer* auxipiv = (integer*)malloc(mn*sizeof(integer)); integer res; dgetrf_ (&m,&n,rp,&m,auxipiv,&res); if(res>0) { res = 0; // FIXME } CHECK(res,res); int k; for (k=0; k=1 && n >=1 && ipivn == mn, BAD_SIZE); DEBUGMSG("lu_l_C"); integer* auxipiv = (integer*)malloc(mn*sizeof(integer)); integer res; zgetrf_ (&m,&n,rp,&m,auxipiv,&res); if(res>0) { res = 0; // FIXME } CHECK(res,res); int k; for (k=0; k=0 && x=0 && y