/* LA: linear algebra C++ interface library Copyright (C) 2008 Jiri Pittner or This program is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program. If not, see . */ //this can be safely included since it contains ifdefs NONCBLAS and NONCLAPACK anyway #include "noncblas.h" #include "vec.h" #include "smat.h" #include "mat.h" #include "nonclass.h" #include "qsort.h" namespace LA { #ifdef FORTRAN_ #define FORNAME(x) x##_ #else #define FORNAME(x) x #endif #define INSTANTIZE(T) \ template void lawritemat(FILE *file,const T *a,int r,int c,const char *form0, \ int nodim,int modulo, int issym); INSTANTIZE(double) INSTANTIZE(complex) INSTANTIZE(int) INSTANTIZE(short) INSTANTIZE(char) INSTANTIZE(long) INSTANTIZE(long long) INSTANTIZE(unsigned char) INSTANTIZE(unsigned short) INSTANTIZE(unsigned int) INSTANTIZE(unsigned long) INSTANTIZE(unsigned long long) #define EPSDET 1e-300 template void lawritemat(FILE *file,const T *a,int r,int c,const char *form0, int nodim,int modulo, int issym) { int i,j; const char *f; /*print out title before %*/ f=form0; skiptext: while (*f && *f !='%' ) {fputc(*f++,file);} if (*f=='%' && f[1]=='%') { fputc(*f,file); f+=2; goto skiptext; } /* this has to be avoided when const arguments should be allowed *f=0; */ /*use the rest as a format for numbers*/ if (modulo) nodim=0; if (nodim==2) fprintf(file,"%d %d\n",r,c); if (nodim==1) fprintf(file,"%d\n",c); if (modulo) { int n1, n2, l, m; char ff[32]; /* prepare integer format for column numbering */ if (sscanf(f+1,"%d",&l) != 1) l=128/modulo; l -= 2; m = l/2; l = l-m; sprintf(ff,"%%%ds%%3d%%%ds", l, m); n1 = 1; while(n1 <= c) { n2=n1+modulo-1; if (n2 > c) n2 = c; /*write block between columns n1 and n2 */ fprintf(file,"\n "); for (i=n1; i<=n2; i++) fprintf(file,ff," ",i," "); fprintf(file,"\n\n"); for (i=1; i<=r; i++) { fprintf(file, "%3d ", i); for (j=n1; j<=n2; j++) { if(issym) { int ii,jj; if (i >= j) { ii=i; jj=j; } else { ii=j; jj=i; } fprintf(file, f, ((complex)a[ii*(ii+1)/2+jj]).real(), ((complex)a[ii*(ii+1)/2+jj]).imag()); } else fprintf(file, f, ((complex)a[(i-1)*c+j-1]).real(), ((complex)a[(i-1)*c+j-1]).imag()); if (j < n2) fputc(' ',file); } fprintf(file, "\n"); } n1 = n2+1; } } else { for (i=1; i<=r; i++) { for (j=1; j<=c; j++) { if (issym) { int ii,jj; if (i >= j) { ii=i; jj=j; } else { ii=j; jj=i; } fprintf(file, f, ((complex)a[ii*(ii+1)/2+jj]).real(), ((complex)a[ii*(ii+1)/2+jj]).imag()); } else fprintf(file,f,((complex)a[(i-1)*c+j-1]).real(), ((complex)a[(i-1)*c+j-1]).imag()); putc(j &A, double *B, const int nrhs, const int ldb, double *det, int n) { int r, *ipiv; if (n==A.nrows() && A.nrows() != A.ncols()) laerror("linear_solve() call for non-square matrix"); A.copyonwrite(); ipiv = new int[A.nrows()]; r = clapack_dgesv(CblasRowMajor, n, nrhs, A[0], A.ncols(), ipiv, B , ldb); if (r < 0) { delete[] ipiv; laerror("illegal argument in lapack_gesv"); } if (det && r==0) { *det = 1.; //take into account some numerical instabilities in dgesv for singular matrices for (int i=0; i0) *det = 0; /* std::cout <<"ipiv = "; for (int i=0; i0 && B) laerror("singular matrix in lapack_gesv"); } void linear_solve(NRMat &A, NRMat *B, double *det, int n) { if(n<=0) n=A.nrows(); //default - whole matrix if (n==A.nrows() && B && A.nrows() != B->ncols() || B && n>B->ncols() ||n>A.nrows()) laerror("incompatible matrices in linear_solve()"); if(B) B->copyonwrite(); linear_solve_do(A,B?(*B)[0]:NULL,B?B->nrows() : 0, B?B->ncols():A.nrows(), det,n); } void linear_solve(NRMat &A, NRVec &B, double *det, int n) { if(n<=0) n=A.nrows(); //default - whole matrix if(n==A.nrows() && A.nrows() != B.size() || n>B.size()||n>A.nrows() ) laerror("incompatible matrices in linear_solve()"); B.copyonwrite(); linear_solve_do(A,&B[0],1,A.nrows(),det,n); } // Next routines are not available in clapack, fotran ones will be used with an // additional swap/transpose of outputs when needed extern "C" void FORNAME(dspsv)(const char *UPLO, const int *N, const int *NRHS, double *AP, int *IPIV, double *B, const int *LDB, int *INFO); static void linear_solve_do(NRSMat &a, double *b, const int nrhs, const int ldb, double *det, int n) { int r, *ipiv; a.copyonwrite(); ipiv = new int[n]; char U = 'U'; FORNAME(dspsv)(&U, &n, &nrhs, a, ipiv, b, &ldb,&r); if (r < 0) { delete[] ipiv; laerror("illegal argument in spsv() call of linear_solve()"); } if (det && r == 0) { *det = 1.; for (int i=1; i0) *det = 0; delete[] ipiv; if (r > 0 && b) laerror("singular matrix in linear_solve(SMat&, Mat*, double*"); } void linear_solve(NRSMat &a, NRMat *B, double *det, int n) { if(n<=0) n=a.nrows(); if (n==a.nrows() && B && a.nrows() != B->ncols() || B && n>B->ncols() || n>a.nrows()) laerror("incompatible matrices in symmetric linear_solve()"); if (B) B->copyonwrite(); linear_solve_do(a,B?(*B)[0]:NULL,B?B->nrows() : 0, B?B->ncols():a.nrows(),det,n); } void linear_solve(NRSMat &a, NRVec &B, double *det, int n) { if(n<=0) n=a.nrows(); if (n==a.nrows() && a.nrows()!= B.size() || n>B.size() || n>a.nrows()) laerror("incompatible matrices in symmetric linear_solve()"); B.copyonwrite(); linear_solve_do(a,&B[0],1,a.nrows(),det,n); } //other version of linear solver based on gesvx extern "C" void FORNAME(zgesvx)(const char *fact, const char *trans, const int *n, const int *nrhs, complex *A, const int *lda, complex *AF, const int *ldaf, const int *ipiv, char *equed, double *R,double *C, complex *B, const int *ldb, complex *X, const int *ldx, double *rcond, double *ferr, double *berr, complex *work, double *rwork, int *info); extern "C" void FORNAME(dgesvx)(const char *fact, const char *trans, const int *n, const int *nrhs, double *A, const int *lda, double *AF, const int *ldaf, const int *ipiv, char *equed, double *R,double *C, double *B, const int *ldb, double *X, const int *ldx, double *rcond, double *ferr, double *berr, double *work, double *rwork, int *info); int linear_solve_x_(NRMat > &A, complex *B, const bool eq, const int nrhs, const int ldb, const char trans) { const int n_rows = A.nrows(); const int n_cols = A.ncols(); if(n_rows != n_cols)laerror("non-squre matrix in linear_solve_x"); const int n = n_rows; const char fact = eq?'E':'N'; char equed = 'B';//fact = 'N' => equed is an output argument int info, lwork; double rcond, ferr[nrhs], berr[nrhs], rwork[2*n]; double R[n], C[n]; complex *AF = new complex[n*n]; complex *work = new complex[2*n]; NRMat > X(n, nrhs); int ipiv[n]; A.copyonwrite(); FORNAME(zgesvx)(&fact, &trans, &n_rows, &nrhs, \ A[0], &n_rows, &AF[0], &n_rows, &ipiv[0], &equed, &R[0], &C[0], \ &B[0], &ldb, X[0], &n_rows, &rcond, &ferr[0], &berr[0], &work[0], &rwork[0], &info); delete[] work; delete[] AF; memcpy(B, X[0], sizeof(complex)*n*nrhs); return info; } int linear_solve_x_(NRMat &A, double *B, const bool eq, const int nrhs, const int ldb, const char trans) { const int n_rows = A.nrows(); const int n_cols = A.ncols(); if(n_rows != n_cols)laerror("non-squre matrix in linear_solve_x"); const int n = n_rows; const char fact = eq?'E':'N'; char equed = 'B';//fact = 'N' => equed is an output argument int info, lwork; double rcond, ferr[nrhs], berr[nrhs], rwork[2*n]; double R[n], C[n]; double *AF = new double[n*n]; double *work = new double[2*n]; NRMat X(n, nrhs); int ipiv[n]; A.copyonwrite(); FORNAME(dgesvx)(&fact, &trans, &n_rows, &nrhs, \ A[0], &n_rows, &AF[0], &n_rows, &ipiv[0], &equed, &R[0], &C[0], \ &B[0], &ldb, X[0], &n_rows, &rcond, &ferr[0], &berr[0], &work[0], &rwork[0], &info); delete[] work; delete[] AF; memcpy(B, X[0], sizeof(double)*n*nrhs); return info; } extern "C" void FORNAME(dsyev)(const char *JOBZ, const char *UPLO, const int *N, double *A, const int *LDA, double *W, double *WORK, const int *LWORK, int *INFO); extern "C" void FORNAME(dsygv)(const int *ITYPE, const char *JOBZ, const char *UPLO, const int *N, double *A, const int *LDA, double *B, const int *LDB, double *W, double *WORK, const int *LWORK, int *INFO); // a will contain eigenvectors (columns if corder==1), w eigenvalues void diagonalize(NRMat &a, NRVec &w, const bool eivec, const bool corder, int n, NRMat *b, const int itype) { int m = a.nrows(); if (m != a.ncols()) laerror("diagonalize() call with non-square matrix"); if (a.nrows() != w.size()) laerror("inconsistent dimension of eigenvalue vector in diagonalize()"); if(n==0) n=m; if(n<0||n>m) laerror("actual dimension out of range in diagonalize"); if(b) if(n>b->nrows() || n> b->ncols()) laerror("wrong B matrix dimension in diagonalize"); a.copyonwrite(); w.copyonwrite(); if(b) b->copyonwrite(); int r = 0; char U ='U'; char vectors = 'V'; if (!eivec) vectors = 'N'; int LWORK = -1; double WORKX; int ldb=0; if(b) ldb=b->ncols(); // First call is to determine size of workspace if(b) FORNAME(dsygv)(&itype,&vectors, &U, &n, a, &m, *b, &ldb, w, &WORKX, &LWORK, &r ); else FORNAME(dsyev)(&vectors, &U, &n, a, &m, w, &WORKX, &LWORK, &r ); LWORK = (int)WORKX; double *WORK = new double[LWORK]; if(b) FORNAME(dsygv)(&itype,&vectors, &U, &n, a, &m, *b,&ldb, w, WORK, &LWORK, &r ); else FORNAME(dsyev)(&vectors, &U, &n, a, &m, w, WORK, &LWORK, &r ); delete[] WORK; if (vectors == 'V' && corder) a.transposeme(n); if (r < 0) laerror("illegal argument in sygv/syev in diagonalize()"); if (r > 0) laerror("convergence problem in sygv/syev in diagonalize()"); } extern "C" void FORNAME(zheev)(const char *JOBZ, const char *UPLO, const int *N, complex *A, const int *LDA, double *W, complex *WORK, const int *LWORK, double *RWORK, int *INFO); extern "C" void FORNAME(zhegv)(const int *ITYPE, const char *JOBZ, const char *UPLO, const int *N, complex *A, const int *LDA, complex *B, const int *LDB, double *W, complex *WORK, const int *LWORK, double *RWORK, int *INFO); // a will contain eigenvectors (columns if corder==1), w eigenvalues void diagonalize(NRMat > &a, NRVec &w, const bool eivec, const bool corder, int n, NRMat > *b, const int itype) { int m = a.nrows(); if (m != a.ncols()) laerror("diagonalize() call with non-square matrix"); if (a.nrows() != w.size()) laerror("inconsistent dimension of eigenvalue vector in diagonalize()"); if(n==0) n=m; if(n<0||n>m) laerror("actual dimension out of range in diagonalize"); if(b) if(n>b->nrows() || n> b->ncols()) laerror("wrong B matrix dimension in diagonalize"); a.copyonwrite(); w.copyonwrite(); if(b) b->copyonwrite(); int r = 0; char U ='U'; char vectors = 'V'; if (!eivec) vectors = 'N'; int LWORK = -1; complex WORKX; int ldb=0; if(b) ldb=b->ncols(); // First call is to determine size of workspace double *RWORK = new double[3*n+2]; if(b) FORNAME(zhegv)(&itype,&vectors, &U, &n, a, &m, *b, &ldb, w, &WORKX, &LWORK, RWORK, &r ); else FORNAME(zheev)(&vectors, &U, &n, a, &m, w, &WORKX, &LWORK, RWORK, &r ); LWORK = (int)WORKX.real(); complex *WORK = new complex[LWORK]; if(b) FORNAME(zhegv)(&itype,&vectors, &U, &n, a, &m, *b,&ldb, w, WORK, &LWORK, RWORK, &r ); else FORNAME(zheev)(&vectors, &U, &n, a, &m, w, WORK, &LWORK, RWORK, &r ); delete[] WORK; delete[] RWORK; if (vectors == 'V' && corder) a.transposeme(n); if (r < 0) laerror("illegal argument in hegv/heev in diagonalize()"); if (r > 0) laerror("convergence problem in hegv/heev in diagonalize()"); } extern "C" void FORNAME(dspev)(const char *JOBZ, const char *UPLO, const int *N, double *AP, double *W, double *Z, const int *LDZ, double *WORK, int *INFO); extern "C" void FORNAME(dspgv)(const int *ITYPE, const char *JOBZ, const char *UPLO, const int *N, double *AP, double *BP, double *W, double *Z, const int *LDZ, double *WORK, int *INFO); // v will contain eigenvectors, w eigenvalues void diagonalize(NRSMat &a, NRVec &w, NRMat *v, const bool corder, int n, NRSMat *b, const int itype) { if(n<=0) n = a.nrows(); if (v) if (v->nrows() != v ->ncols() || n > v->nrows() || n > a.nrows()) laerror("diagonalize() call with inconsistent dimensions"); if (n==a.nrows() && n != w.size() || n>w.size()) laerror("inconsistent dimension of eigenvalue vector"); if(b) if(n>b->nrows() || n> b->ncols()) laerror("wrong B matrix dimension in diagonalize"); a.copyonwrite(); w.copyonwrite(); if(v) v->copyonwrite(); if(b) b->copyonwrite(); int r = 0; char U = 'U'; char job = v ? 'v' : 'n'; double *WORK = new double[3*n]; int ldv=v?v->ncols():n; if(b) FORNAME(dspgv)(&itype,&job, &U, &n, a, *b, w, v?(*v)[0]:(double *)0, &ldv, WORK, &r ); else FORNAME(dspev)(&job, &U, &n, a, w, v?(*v)[0]:(double *)0, &ldv, WORK, &r ); delete[] WORK; if (v && corder) v->transposeme(n); if (r < 0) laerror("illegal argument in spgv/spev in diagonalize()"); if (r > 0) laerror("convergence problem in spgv/spev in diagonalize()"); } extern "C" void FORNAME(zhpev)(const char *JOBZ, const char *UPLO, const int *N, complex *AP, double *W, complex *Z, const int *LDZ, complex *WORK, double *RWORK, int *INFO); extern "C" void FORNAME(zhpgv)(const int *ITYPE, const char *JOBZ, const char *UPLO, const int *N, complex *AP, complex *BP, double *W, complex *Z, const int *LDZ, complex *WORK, double *RWORK, int *INFO); // v will contain eigenvectors, w eigenvalues void diagonalize(NRSMat > &a, NRVec &w, NRMat > *v, const bool corder, int n, NRSMat > *b, const int itype) { if(n<=0) n = a.nrows(); if (v) if (v->nrows() != v ->ncols() || n > v->nrows() || n > a.nrows()) laerror("diagonalize() call with inconsistent dimensions"); if (n==a.nrows() && n != w.size() || n>w.size()) laerror("inconsistent dimension of eigenvalue vector"); if(b) if(n>b->nrows() || n> b->ncols()) laerror("wrong B matrix dimension in diagonalize"); a.copyonwrite(); w.copyonwrite(); if(v) v->copyonwrite(); if(b) b->copyonwrite(); int r = 0; char U = 'U'; char job = v ? 'v' : 'n'; complex *WORK = new complex[2*n]; double *RWORK = new double[3*n]; int ldv=v?v->ncols():n; if(b) FORNAME(zhpgv)(&itype,&job, &U, &n, a, *b, w, v?(*v)[0]:(complex *)0, &ldv, WORK, RWORK, &r ); else FORNAME(zhpev)(&job, &U, &n, a, w, v?(*v)[0]:(complex *)0, &ldv, WORK, RWORK, &r ); delete[] WORK; delete[] RWORK; if (v && corder) v->transposeme(n); if (r < 0) laerror("illegal argument in hpgv/hpev in diagonalize()"); if (r > 0) laerror("convergence problem in hpgv/hpev in diagonalize()"); } extern "C" void FORNAME(dgesvd)(const char *JOBU, const char *JOBVT, const int *M, const int *N, double *A, const int *LDA, double *S, double *U, const int *LDU, double *VT, const int *LDVT, double *WORK, const int *LWORK, int *INFO ); void singular_decomposition(NRMat &a, NRMat *u, NRVec &s, NRMat *v, const bool corder, int m, int n) { int m0 = a.nrows(); int n0 = a.ncols(); if(m<=0) m=m0; if(n<=0) n=n0; if(n>n0 || m>m0) laerror("bad dimension in singular_decomposition"); if (u) if (m > u->nrows() || m> u->ncols()) laerror("inconsistent dimension of U Mat in singular_decomposition()"); if (s.size() < m && s.size() < n) laerror("inconsistent dimension of S Vec in singular_decomposition()"); if (v) if (n > v->nrows() || n > v->ncols()) laerror("inconsistent dimension of V Mat in singular_decomposition()"); a.copyonwrite(); s.copyonwrite(); if (u) u->copyonwrite(); if (v) v->copyonwrite(); // C-order (transposed) input and swap u,v matrices, // v should be transposed at the end char jobu = u ? 'A' : 'N'; char jobv = v ? 'A' : 'N'; double work0; int lwork = -1; int r; FORNAME(dgesvd)(&jobv, &jobu, &n, &m, a, &n0, s, v?(*v)[0]:0, &n0, u?(*u)[0]:0, &m0, &work0, &lwork, &r); lwork = (int) work0; double *work = new double[lwork]; FORNAME(dgesvd)(&jobv, &jobu, &n, &m, a, &n0, s, v?(*v)[0]:0, &n0, u?(*u)[0]:0, &m0, work, &lwork, &r); delete[] work; if (v && corder) v->transposeme(n); if (r < 0) laerror("illegal argument in gesvd() of singular_decomposition()"); if (r > 0) laerror("convergence problem in gesvd() of ingular_decomposition()"); } //QR decomposition //extern "C" void FORNAME(dgeqrf)(const int *M, const int *N, double *A, const int *LDA, double *TAU, double *WORK, int *LWORK, int *INFO); extern "C" void FORNAME(dgeev)(const char *JOBVL, const char *JOBVR, const int *N, double *A, const int *LDA, double *WR, double *WI, double *VL, const int *LDVL, double *VR, const int *LDVR, double *WORK, const int *LWORK, int *INFO ); extern "C" void FORNAME(dggev)(const char *JOBVL, const char *JOBVR, const int *N, double *A, const int *LDA, double *B, const int *LDB, double *WR, double *WI, double *WBETA, double *VL, const int *LDVL, double *VR, const int *LDVR, double *WORK, const int *LWORK, int *INFO ); //statics for sorting static int *gdperm; static double *gdwr, *gdwi, *gdbeta; //compare methods static double realonly(const int i, const int j) { if(gdbeta) { if(gdbeta[i]==0. && gdbeta[j]!=0) return 1.; if(gdbeta[j]==0. && gdbeta[i]!=0) return -1.; if(gdbeta[i]==0. && gdbeta[j]==0) return 0.; double tmp = gdwr[i]/gdbeta[i]-gdwr[j]/gdbeta[j]; if(tmp) return tmp; return gdwi[j]/gdbeta[j]-gdwi[i]/gdbeta[i]; } //else double tmp = gdwr[i]-gdwr[j]; if(tmp) return tmp; return gdwi[j]-gdwi[i]; } static double realfirst(const int i, const int j) { if(gdwi[i] && ! gdwi[j]) return 1.; if(!gdwi[i] && gdwi[j]) return -1.; return realonly(i,j); } static double (* gdcompar[2])(const int, const int) = {&realonly, &realfirst}; //swap method static void gdswap(const int i, const int j) { double tmp; int itmp; itmp=gdperm[i]; gdperm[i]=gdperm[j]; gdperm[j]=itmp; tmp=gdwr[i]; gdwr[i]=gdwr[j]; gdwr[j]=tmp; tmp=gdwi[i]; gdwi[i]=gdwi[j]; gdwi[j]=tmp; if(gdbeta) {tmp=gdbeta[i]; gdbeta[i]=gdbeta[j]; gdbeta[j]=tmp;} } void gdiagonalize(NRMat &a, NRVec &wr, NRVec &wi, NRMat *vl, NRMat *vr, const bool corder, int n, const int sorttype, const int biorthonormalize, NRMat *b, NRVec *beta) { if(n<=0) n = a.nrows(); if (n > a.ncols() || n>a.nrows() ) laerror("gdiagonalize() call for a non-square matrix"); if (n > wr.size()) laerror("inconsistent dimension of eigen vector in gdiagonalize()"); if (vl) if (n > vl->nrows() || n > vl->ncols()) laerror("inconsistent dimension of vl in gdiagonalize()"); if (vr) if (n > vr->nrows() || n > vr->ncols()) laerror("inconsistent dimension of vr in gdiagonalize()"); if (beta) if(n > beta ->size()) laerror("inconsistent dimension of beta in gdiagonalize()"); if(b) if(n > b->nrows() || n > b->ncols()) laerror("inconsistent dimension of b in gdiagonalize()"); if(b && !beta || beta && !b) laerror("missing array for generalized diagonalization"); a.copyonwrite(); wr.copyonwrite(); wi.copyonwrite(); if (vl) vl->copyonwrite(); if (vr) vr->copyonwrite(); if (beta) beta->copyonwrite(); if (b) b->copyonwrite(); char jobvl = vl ? 'V' : 'N'; char jobvr = vr ? 'V' : 'N'; double work0; int lwork = -1; int r; int lda=a.ncols(); int ldb=0; if(b) ldb=b->ncols(); int ldvl= vl?vl->ncols():lda; int ldvr= vr?vr->ncols():lda; if(b) FORNAME(dggev)(&jobvr, &jobvl, &n, a, &lda, *b, &ldb, wr, wi, *beta, vr?vr[0]:(double *)0, &ldvr, vl?vl[0]:(double *)0, &ldvl, &work0, &lwork, &r); else FORNAME(dgeev)(&jobvr, &jobvl, &n, a, &lda, wr, wi, vr?vr[0]:(double *)0, &ldvr, vl?vl[0]:(double *)0, &ldvl, &work0, &lwork, &r); lwork = (int) work0; double *work = new double[lwork]; if(b) FORNAME(dggev)(&jobvr, &jobvl, &n, a, &lda, *b, &ldb, wr, wi, *beta, vr?vr[0]:(double *)0, &ldvr, vl?vl[0]:(double *)0, &ldvl, work, &lwork, &r); else FORNAME(dgeev)(&jobvr, &jobvl, &n, a, &lda, wr, wi, vr?vr[0]:(double *)0, &ldvr, vl?vl[0]:(double *)0, &ldvl, work, &lwork, &r); delete[] work; //std::cout <<"TEST dgeev\n"< 0) laerror("convergence problem in ggev/geev in gdiagonalize()"); if(biorthonormalize && vl && vr) { if(b || beta) laerror("@@@ biorthonormalize not implemented yet for generalized non-symmetric eigenproblem");//metric b would be needed int i=0; while(i0) { NRVec perm(n); for(int i=0; itransposeme(n); if (vr) vr->transposeme(n); } } void gdiagonalize(NRMat &a, NRVec< complex > &w, NRMat< complex >*vl, NRMat< complex > *vr, const bool corder, int n, const int sorttype, const int biorthonormalize, NRMat *b, NRVec *beta) { if(!corder) laerror("gdiagonalize() corder 0 not implemented"); if(n<=0) n = a.nrows(); if(n> a.nrows() || n == a.nrows() && n != a.ncols()) laerror("gdiagonalize() call for a non-square matrix"); NRVec wr(n), wi(n); NRMat *rvl = 0; NRMat *rvr = 0; if (vl) rvl = new NRMat(n, n); if (vr) rvr = new NRMat(n, n); gdiagonalize(a, wr, wi, rvl, rvr, 0, n, sorttype, biorthonormalize, b, beta); //process the results into complex matrices int i; for (i=0; i(wr[i], wi[i]); if (rvl || rvr) { i = 0; while (i < n) { if (wi[i] == 0) { if (vl) for (int j=0; j((*rvl)[i][j], (*rvl)[i+1][j]); (*vl)[i+1][j] = complex((*rvl)[i][j], -(*rvl)[i+1][j]); } if (vr) for (int j=0; j((*rvr)[i][j], (*rvr)[i+1][j]); (*vr)[i+1][j] = complex((*rvr)[i][j], -(*rvr)[i+1][j]); } i += 2; } } } if (rvl) delete rvl; if (rvr) delete rvr; } const NRMat realpart(const NRMat< complex > &a) { NRMat result(a.nrows(), a.ncols()); cblas_dcopy(a.nrows()*a.ncols(), (const double *)a[0], 2, result, 1); return result; } const NRMat imagpart(const NRMat< complex > &a) { NRMat result(a.nrows(), a.ncols()); cblas_dcopy(a.nrows()*a.ncols(), (const double *)a[0]+1, 2, result, 1); return result; } const NRMat< complex > realmatrix (const NRMat &a) { NRMat > result(a.nrows(), a.ncols()); cblas_dcopy(a.nrows()*a.ncols(), a, 1, (double *)result[0], 2); return result; } const NRMat< complex > imagmatrix (const NRMat &a) { NRMat< complex > result(a.nrows(), a.ncols()); cblas_dcopy(a.nrows()*a.ncols(), a, 1, (double *)result[0]+1, 2); return result; } const NRMat< complex > complexmatrix (const NRMat &re, const NRMat &im) { if(re.nrows()!=im.nrows() || re.ncols() != im.ncols()) laerror("incompatible sizes of real and imaginary parts"); NRMat< complex > result(re.nrows(), re.ncols()); cblas_dcopy(re.nrows()*re.ncols(), re, 1, (double *)result[0], 2); cblas_dcopy(re.nrows()*re.ncols(), im, 1, (double *)result[0]+1, 2); return result; } NRMat matrixfunction(NRMat a, complex (*f)(const complex &), const bool adjust) { int n = a.nrows(); NRMat< complex > u(n, n), v(n, n); NRVec< complex > w(n); /* NRMat > a0=complexify(a); */ gdiagonalize(a, w, &u, &v);//a gets destroyed, eigenvectors are rows NRVec< complex > z = diagofproduct(u, v, 1, 1); /* std::cout <<"TEST matrixfunction\n"< > wz(n); for (int i=0; i > r(n, n); r.gemm(0.0, v, 'c', u, 'n', 1.0); double inorm = cblas_dnrm2(n*n, (double *)r[0]+1, 2); if (inorm > 1e-10) { std::cout << "norm = " << inorm << std::endl; laerror("nonzero norm of imaginary part of real matrixfunction"); } return realpart(r); } NRMat matrixfunction(NRSMat a, double (*f) (double)) { int n = a.nrows(); NRVec w(n); NRMat v(n, n); diagonalize(a, w, &v, 0); for (int i=0; i u = v; v.diagmultl(w); NRMat r(n, n); r.gemm(0.0, u, 't', v, 'n', 1.0); return r; } NRMat realmatrixfunction(NRMat a, double (*f) (const double)) { int n = a.nrows(); NRVec w(n); diagonalize(a, w, true, false); for (int i=0; i u = a; a.diagmultl(w); NRMat r(n, n); r.gemm(0.0, u, 't', a, 'n', 1.0); return r; } NRMat > complexmatrixfunction(NRMat a, double (*fre) (const double), double (*fim) (const double)) { int n = a.nrows(); NRVec wre(n),wim(n); diagonalize(a, wre, true, false); for (int i=0; i u = a; NRMat b = a; a.diagmultl(wre); b.diagmultl(wim); NRMat t(n,n),tt(n,n); t.gemm(0.0, u, 't', a, 'n', 1.0); tt.gemm(0.0, u, 't', b, 'n', 1.0); NRMat > r(n, n); for (int i=0; i(t(i,j),tt(i,j)); return r; } // instantize template to an addresable function complex myclog (const complex &x) { return log(x); } complex mycexp (const complex &x) { return std::exp(x); } complex sqrtinv (const complex &x) { return 1./std::sqrt(x); } double sqrtinv (const double x) { return 1./std::sqrt(x); } NRMat log(const NRMat &a) { return matrixfunction(a, &myclog, 1); } NRMat exp0(const NRMat &a) { return matrixfunction(a, &mycexp, 1); } const NRVec diagofproduct(const NRMat &a, const NRMat &b, bool trb, bool conjb) { if (trb && (a.nrows() != b.nrows() || a.ncols() != b.ncols()) || !trb && (a.nrows() != b.ncols() || a.ncols() != b.nrows())) laerror("incompatible Mats in diagofproduct()"); NRVec result(a.nrows()); if (trb) for(int i=0; i > diagofproduct(const NRMat< complex > &a, const NRMat< complex > &b, bool trb, bool conjb) { if (trb && (a.nrows() != b.nrows() || a.ncols() != b.ncols()) || !trb && (a.nrows() != b.ncols() || a.ncols() != b.nrows())) laerror("incompatible Mats in diagofproduct()"); NRVec< complex > result(a.nrows()); if (trb) { if (conjb) { for(int i=0; i &a, const NRMat &b, bool trb) { if (trb && (a.nrows() != b.nrows() || a.ncols() != b.ncols()) || !trb && (a.nrows() != b.ncols() || a.ncols() != b.nrows())) laerror("incompatible Mats in trace2()"); if (trb) return cblas_ddot(a.nrows()*a.ncols(), a, 1, b, 1); double sum = 0.0; for (int i=0; i &a, const NRSMat &b, const bool diagscaled) { if (a.nrows() != b.nrows()) laerror("incompatible SMats in trace2()"); //double r = 0; for (int i=0; i &a, const NRMat &b, const bool diagscaled) { if (a.nrows() != b.nrows()||b.nrows()!=b.ncols()) laerror("incompatible SMats in trace2()"); double r=0; int i, j, k=0; for (i=0; i *A, const int *LDA, int *INFO); void cholesky(NRMat &a, bool upper) { if(a.nrows()!=a.ncols()) laerror("matrix must be square in Cholesky"); int lda=a.ncols(); int n=a.nrows(); char uplo=upper?'U':'L'; int info; a.copyonwrite(); FORNAME(dpotrf)(&uplo, &n, a, &lda, &info); if(info) {std::cerr << "Lapack error "< > &a, bool upper) { if(a.nrows()!=a.ncols()) laerror("matrix must be square in Cholesky"); int lda=a.ncols(); int n=a.nrows(); char uplo=upper?'U':'L'; int info; a.copyonwrite(); a.transposeme();//switch to Fortran order FORNAME(zpotrf)(&uplo, &n, a, &lda, &info); if(info) {std::cerr << "Lapack error "< *A, const int *LDA, double *WORK); extern "C" double FORNAME(dlange)( const char *NORM, const int *M, const int *N, double *A, const int *LDA, double *WORK); double MatrixNorm(NRMat > &A, const char norm) { const char TypNorm = (tolower(norm) == 'o')?'I':'O'; //switch c-order/fortran-order const int M = A.nrows(); const int N = A.ncols(); double work[M]; const double ret = FORNAME(zlange)(&TypNorm, &M, &N, A[0], &M, &work[0]); return ret; } double MatrixNorm(NRMat &A, const char norm) { const char TypNorm = (tolower(norm) == 'o')?'I':'O'; //switch c-order/fortran-order const int M = A.nrows(); const int N = A.ncols(); double work[M]; const double ret = FORNAME(dlange)(&TypNorm, &M, &N, A[0], &M, &work[0]); return ret; } //condition number extern "C" void FORNAME(zgecon)( const char *norm, const int *n, complex *A, const int *LDA, const double *anorm, double *rcond, complex *work, double *rwork, int *info); extern "C" void FORNAME(dgecon)( const char *norm, const int *n, double *A, const int *LDA, const double *anorm, double *rcond, double *work, double *rwork, int *info); double CondNumber(NRMat > &A, const char norm) { const char TypNorm = (tolower(norm) == 'o')?'I':'O'; //switch c-order/fortran-order const int N = A.nrows(); double Norma(0.0), ret(0.0); int info; complex *work; double *rwork; if(N != A.ncols()){ laerror("nonsquare matrix in zgecon"); return 0.0; } work = new complex[2*N]; rwork = new double[2*N]; Norma = MatrixNorm(A, norm); FORNAME(zgecon)(&TypNorm, &N, A[0], &N, &Norma, &ret, &work[0], &rwork[0], &info); delete[] work; delete[] rwork; return ret; } double CondNumber(NRMat &A, const char norm) { const char TypNorm = (tolower(norm) == 'o')?'I':'O'; //switch c-order/fortran-order const int N = A.nrows(); double Norma(0.0), ret(0.0); int info; double *work; double *rwork; if(N != A.ncols()){ laerror("nonsquare matrix in zgecon"); return 0.0; } work = new double[2*N]; rwork = new double[2*N]; Norma = MatrixNorm(A, norm); FORNAME(dgecon)(&TypNorm, &N, A[0], &N, &Norma, &ret, &work[0], &rwork[0], &info); delete[] work; delete[] rwork; return ret; } #ifdef obsolete void gendiagonalize(NRMat &a, NRVec &w, NRMat b, int n) { if(a.nrows()!=a.ncols() || a.nrows()!=w.size() || a.nrows()!=b.nrows() || b.nrows()!=b.ncols() ) laerror("incompatible Mats in gendiagonalize"); a.copyonwrite(); w.copyonwrite(); b.copyonwrite(); int m=w.size(); NRVec dl(m); int i,j; double x; if(n==0) n=m; if(n<0 || n>m) laerror("actual dimension in gendiagonalize out of range"); //transform the problem to usual diagonalization //cholesky decompose in b and dl for(i=0; i=0; --i)//component loop { if(i &v) { int n=v.nrows(); double det=determinant(v); int nchange=0; for(int i=0; i