extern "C" { #include "atlas_enum.h" #include "clapack.h" } #include "la.h" #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(char) 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, NRMat *B, double *det) { int r, *ipiv; if (A.nrows() != A.ncols()) laerror("linear_solve() call for non-square matrix"); if (B && A.nrows() != B->ncols()) laerror("incompatible matrices in linear_solve()"); A.copyonwrite(); if (B) B->copyonwrite(); ipiv = new int[A.nrows()]; r = clapack_dgesv(CblasRowMajor, A.nrows(), B ? B->nrows() : 0, A[0], A.ncols(), ipiv, B ? B[0] : (double *)0, B ? B->ncols() : A.nrows()); if (r < 0) { delete[] ipiv; laerror("illegal argument in lapack_gesv"); } if (det && r>=0) { *det = A[0][0]; for (int i=1; i0 && B) laerror("singular matrix in lapack_gesv"); } // Next routines are not available in clapack, fotran ones will b 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); void linear_solve(NRSMat &a, NRMat *b, double *det) { int r, *ipiv; if (det) cerr << "@@@ sign of the determinant not implemented correctly yet\n"; if (b && a.nrows() != b->ncols()) laerror("incompatible matrices in symmetric linear_solve()"); a.copyonwrite(); if (b) b->copyonwrite(); ipiv = new int[a.nrows()]; char U = 'U'; int n = a.nrows(); int nrhs = 0; if (b) nrhs = b->nrows(); int ldb = b ? b->ncols() : a.nrows(); FORNAME(dspsv)(&U, &n, &nrhs, a, ipiv, b?(*b)[0]:0, &ldb,&r); if (r < 0) { delete[] ipiv; laerror("illegal argument in spsv() call of linear_solve()"); } if (det && r >= 0) { *det = a(0,0); for (int i=1; i 0 && b) laerror("singular matrix in linear_solve(SMat&, Mat*, double*"); } 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); // a will contain eigenvectors (columns if corder==1), w eigenvalues void diagonalize(NRMat &a, NRVec &w, const bool eivec, const bool corder) { int n = a.nrows(); if (n != a.ncols()) laerror("diagonalize() call with non-square matrix"); if (a.nrows() != w.size()) laerror("inconsistent dimension of eigenvalue vector in diagonalize()"); a.copyonwrite(); w.copyonwrite(); int r = 0; char U ='U'; char vectors = 'V'; if (!eivec) vectors = 'N'; int LWORK = -1; double WORKX; // First call is to determine size of workspace FORNAME(dsyev)(&vectors, &U, &n, a, &n, w, (double *)&WORKX, &LWORK, &r ); LWORK = (int)WORKX; double *WORK = new double[LWORK]; FORNAME(dsyev)(&vectors, &U, &n, a, &n, w, WORK, &LWORK, &r ); delete[] WORK; if (vectors == 'V' && corder) a.transposeme(); if (r < 0) laerror("illegal argument in syev() of diagonalize()"); if (r > 0) laerror("convergence problem in syev() of 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); // v will contain eigenvectors, w eigenvalues void diagonalize(NRSMat &a, NRVec &w, NRMat *v, const bool corder) { int n = a.nrows(); if (v) if (v->nrows() != v ->ncols() || n != v->nrows()) laerror("diagonalize() call with inconsistent dimensions"); if (n != w.size()) laerror("inconsistent dimension of eigenvalue vector"); a.copyonwrite(); w.copyonwrite(); int r = 0; char U = 'U'; char job = v ? 'v' : 'n'; double *WORK = new double[3*n]; FORNAME(dspev)(&job, &U, &n, a, w, v?(*v)[0]:(double *)0, &n, WORK, &r ); delete[] WORK; if (v && corder) v->transposeme(); if (r < 0) laerror("illegal argument in spev() of diagonalize()"); if (r > 0) laerror("convergence problem in spev() of 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 = a.nrows(); int n = a.ncols(); 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, &n, s, v?(*v)[0]:0, &n, u?(*u)[0]:0, &m, &work0, &lwork, &r); lwork = (int) work0; double *work = new double[lwork]; FORNAME(dgesvd)(&jobv, &jobu, &n, &m, a, &n, s, v?(*v)[0]:0, &n, u?(*u)[0]:0, &m, &work0, &lwork, &r); delete[] work; if (v && corder) v->transposeme(); if (r < 0) laerror("illegal argument in gesvd() of singular_decomposition()"); if (r > 0) laerror("convergence problem in gesvd() of ingular_decomposition()"); } 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 ); void gdiagonalize(NRMat &a, NRVec &wr, NRVec &wi, NRMat *vl, NRMat *vr, const bool corder) { int n = a.nrows(); if (n != a.ncols()) 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()"); a.copyonwrite(); wr.copyonwrite(); wi.copyonwrite(); if (vl) vl->copyonwrite(); if (vr) vr->copyonwrite(); char jobvl = vl ? 'V' : 'N'; char jobvr = vr ? 'V' : 'N'; double work0; int lwork = -1; int r; FORNAME(dgeev)(&jobvr, &jobvl, &n, a, &n, wr, wi, vr?vr[0]:(double *)0, &n, vl?vl[0]:(double *)0, &n, &work0, &lwork, &r); lwork = (int) work0; double *work = new double[lwork]; FORNAME(dgeev)(&jobvr, &jobvl, &n, a, &n, wr, wi, vr?vr[0]:(double *)0, &n, vl?vl[0]:(double *)0, &n, &work0, &lwork, &r); delete[] work; if (corder) { if (vl) vl->transposeme(); if (vr) vr->transposeme(); } if (r < 0) laerror("illegal argument in geev() of gdiagonalize()"); if (r > 0) laerror("convergence problem in geev() of gdiagonalize()"); } void gdiagonalize(NRMat &a, NRVec< complex > &w, NRMat< complex >*vl, NRMat< complex > *vr) { int n = a.nrows(); if(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); //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; } 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); gdiagonalize(a, w, &u, &v); NRVec< complex > z = diagofproduct(u, v, 1, 1); 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) { cout << "norm = " << inorm << 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; } // instantize template to an addresable function complex myclog (const complex &x) { return log(x); } NRMat log(const NRMat &a) { return matrixfunction(a, &myclog, 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 diagofproduct()"); 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 = 2.0*cblas_ddot(a.nrows()*(a.nrows()+1)/2, a, 1, b, 1); if (diagscaled) return r; for (int i=0; i &a, NRVec &w, NRMat b) { 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 n=w.size(); NRVec dl(n); int i,j; double x; //transform the problem to usual diagonalization /* c c this routine is a translation of the algol procedure reduc1, c num. math. 11, 99-110(1968) by martin and wilkinson. c handbook for auto. comp., vol.ii-linear algebra, 303-314(1971). c c this routine reduces the generalized symmetric eigenproblem c ax=(lambda)bx, where b is positive definite, to the standard c symmetric eigenproblem using the cholesky factorization of b. c c on input c c nm must be set to the row dimension of two-dimensional c array parameters as declared in the calling program c dimension statement. c c n is the order of the matrices a and b. if the cholesky c factor l of b is already available, n should be prefixed c with a minus sign. c c a and b contain the real symmetric input matrices. only the c full upper triangles of the matrices need be supplied. if c n is negative, the strict lower triangle of b contains, c instead, the strict lower triangle of its cholesky factor l. c c dl contains, if n is negative, the diagonal elements of l. c c on output c c a contains in its full lower triangle the full lower triangle c of the symmetric matrix derived from the reduction to the c standard form. the strict upper triangle of a is unaltered. c c b contains in its strict lower triangle the strict lower c triangle of its cholesky factor l. the full upper c triangle of b is unaltered. c c dl contains the diagonal elements of l. c c ierr is set to c zero for normal return, c 7*n+1 if b is not positive definite. c c questions and comments should be directed to burton s. garbow, c mathematics and computer science div, argonne national laboratory c c this version dated august 1983. c */ // .......... form l in the arrays b and dl .......... for(i=0; i=0; --i)//component loop { if(i