diff --git a/test/triqs/arrays/eigenelements.cpp b/test/triqs/arrays/eigenelements.cpp index f20ab538..14cdc03d 100644 --- a/test/triqs/arrays/eigenelements.cpp +++ b/test/triqs/arrays/eigenelements.cpp @@ -19,57 +19,109 @@ * ******************************************************************************/ #include "./common.hpp" -#include -#include -#include #include +#include #include using namespace triqs::arrays; using namespace triqs::arrays::linalg; +using dcomplex = std::complex; + +template void check_eig(matrix M, matrix vectors, array values) { + auto _ = range(); + for (auto i : range(0,first_dim(M))) { + std::cerr << "check "<< i << std::endl; + std::cerr << (M -values(i))* vectors(i, _)< void test(M A) { + auto w = eigenelements(make_clone(A)); + std::cerr << "A = " << A << std::endl; + std::cerr << " values = " < -namespace triqs { namespace arrays { namespace linalg { +namespace triqs { +namespace arrays { + namespace linalg { - /** - * A worker to call lapack routine with the matrices. - * - * Handles both real and complex case. - * - * Usage : - * - construct from a VIEW of a matrix, that MUST be contiguous. - * - call invoke() - * - read the eigenvalues/vectors in values and vectors resp. - * NB : the content of the matrix is destroyed by the computation (it is .vectors() in fact, if Compute_Eigenvectors is true). - * For a one shot usage, prefer eigenelements, eigenvalues functions. - */ - template struct eigenelements_worker; - template struct eigenelements_worker_base; + extern "C" { + void TRIQS_FORTRAN_MANGLING(dsyev)(char *, char *, // JOBZ and UPLO + int &, // Matrix Size + double[], // matrix + int &, // LDA of the matrix + double[], // Eigenvalues array + double[], // WORK + int &, // LWORK + int & // INFO + ); - template - struct eigenelements_worker_base { - private: - void operator = ( eigenelements_worker_base const & x); - protected: - matrix_view mat; - triqs::arrays::vector ev; - triqs::arrays::vector work; - int dim,info,lwork; - char uplo,compz; - bool has_run; - eigenelements_worker_base ( matrix_view the_matrix) : mat(the_matrix) { - if (mat.is_empty()) TRIQS_RUNTIME_ERROR<<"eigenelements_worker : the matrix is empty : matrix = "< [], // matrix + int &, // LDA of the matrix + double[], // Eigenvalues array + std::complex [], // WORK + int &, // LWORK + double[], // WORK2 + int & // INFO + ); + } - public : - array values() const { - if (!has_run) TRIQS_RUNTIME_ERROR<<"eigenelements_worker has not been invoked !"; - return ev; - } - }; - - //-------------------------------- - - template - struct eigenelements_worker_base : eigenelements_worker_base { - protected: - eigenelements_worker_base ( matrix_view the_matrix) : eigenelements_worker_base (the_matrix) {this->compz='V'; } + /** + * A worker to call lapack routine with the matrices. Handles both real and complex case. + */ + template class eigenelements_worker { public: - typename matrix_view::regular_type vectors() const { - if (!this->has_run) TRIQS_RUNTIME_ERROR<<"eigenelements_worker has not been invoked !"; - return this->mat; - } + eigenelements_worker() = default; + + /// The eigenvalues + template array eigenvalues(M &mat) const { + _prepare(mat); + _invoke(is_complex(), 'N', mat); + return ev; + } + + /// The eigensystems + template std::pair, matrix> eigenelements(M &mat) const { + _prepare(mat); + _invoke(is_complex(), 'V', mat); + return {ev, _conj(mat, is_complex())}; + } + + private: + mutable array ev, work2; // work2 only used for T complex + mutable array work; + mutable int dim, lwork, info; + + // dispatch the implementation of invoke for T = double or complex + void _invoke(std::false_type, char compz, matrix_view mat) const { // the case double (is_complex = false) + char uplo = 'U'; + TRIQS_FORTRAN_MANGLING(dsyev)(&compz, &uplo, dim, mat.data_start(), dim, ev.data_start(), work.data_start(), lwork, info); + if (info) TRIQS_RUNTIME_ERROR << "eigenelements_worker :error code dsyev : " << info << " for matrix " << mat; + } + + void _invoke(std::true_type, char compz, matrix_view> mat) const { // the case complex (is_complex = true) + char uplo = 'U'; + TRIQS_FORTRAN_MANGLING(zheev)(&compz, &uplo, dim, mat.data_start(), dim, ev.data_start(), work.data_start(), lwork, + work2.data_start(), info); + if (info) TRIQS_RUNTIME_ERROR << "eigenelements_worker :error code zheev : " << info << " for matrix " << mat; + } + + template void _prepare(M const &mat) const { + if (mat.is_empty()) TRIQS_RUNTIME_ERROR << "eigenelements_worker : the matrix is empty : matrix = " << mat << " "; + if (!mat.is_square()) TRIQS_RUNTIME_ERROR << "eigenelements_worker : the matrix " << mat << " is not square "; + if (!mat.indexmap().is_contiguous()) + TRIQS_RUNTIME_ERROR << "eigenelements_worker : the matrix " << mat << " is not contiguous in memory"; + dim = first_dim(mat); + ev.resize(dim); + lwork = 64 * dim; + work.resize(lwork); + if (is_complex::value) work2.resize(lwork); + } + + template matrix _conj(M const &m, std::false_type) const { return m; } + + // impl : since we call fortran lapack, if the order is C (!), the matrix is transposed, or conjugated, so we obtain + // the conjugate of the eigenvectors... Fix #119. + // Do nothing if the order is fortran already... + template matrix> _conj(M const &m, std::true_type) const { + if (m.memory_layout_is_c()) + return conj(m); + else + return m.transpose(); // the matrix mat is understood as a fortran matrix. After the lapack, in memory, it contains the + // correct answer. + // but since it is a fortran matrix, the C will see its transpose. We need to compensate this transpose (!). + } }; - //-------------------------------- - extern "C" { - void TRIQS_FORTRAN_MANGLING(dsyev)(char*,char*, // JOBZ and UPLO - int &, // Matrix Size - double[], // matrix - int&, // LDA of the matrix - double[], // Eigenvalues array - double[], // WORK - int&, // LWORK - int & // INFO - ); + //-------------------------------- + /** + * Simple diagonalization call, return all eigenelements. + * Handles both real and complex case. + * @param M : the matrix or view. MUST be contiguous. It is modified by the call. + * If you wish not to modify it, call eigenelements(make_clone(A)) + */ + template + std::pair, matrix::value_type>> eigenelements(M &&m) { + return eigenelements_worker::value_type>().eigenelements(m); + } - void TRIQS_FORTRAN_MANGLING(zheev)(char*,char*, // JOBZ and UPLO - int &, // Matrix Size - std::complex [], // matrix - int&, // LDA of the matrix - double[], // Eigenvalues array - std::complex[], // WORK - int &, // LWORK - double[], // WORK2 - int & // INFO - ); + //-------------------------------- + + /** + * Simple diagonalization call, returning only the eigenvalues. + * Handles both real and complex case. + * @param M : the matrix VIEW : it MUST be contiguous + * @param take_copy : makes a copy of the matrix before calling lapack, so that the original is preserved. + * if false : no copy is made and the content of the matrix M is destroyed. + * if true : a copy is made, M is preserved, but of course it is slower... + */ + template array eigenvalues(M &&m) { + return eigenelements_worker::value_type>().eigenvalues(m); + } } - - //-------------------------------- - - template - struct eigenelements_worker< matrix_view ,Compute_Eigenvectors > :eigenelements_worker_base { - eigenelements_worker ( matrix_view the_matrix) : eigenelements_worker_base (the_matrix) {} - void invoke() { - int info; - //fortran_int_t info; - TRIQS_FORTRAN_MANGLING(dsyev) (&this->compz,&this->uplo,this->dim,this->mat.data_start(),this->dim,this->ev.data_start(),this->work.data_start(),this->lwork,info); - if (info) TRIQS_RUNTIME_ERROR<<"eigenelements_worker :error code dsyev : "<mat; - this->has_run = true; - } - }; - - //-------------------------------- - - template - struct eigenelements_worker< matrix_view, Opt>,Compute_Eigenvectors > :eigenelements_worker_base,Opt,Compute_Eigenvectors> { - triqs::arrays::vector work2; - public : - eigenelements_worker ( matrix_view ,Opt> the_matrix) : eigenelements_worker_base,Opt,Compute_Eigenvectors> (the_matrix) { work2.resize(this->lwork);} - void invoke() { - int info; - TRIQS_FORTRAN_MANGLING(zheev) (&this->compz,&this->uplo,this->dim,this->mat.data_start(), - this->dim,this->ev.data_start(),this->work.data_start(),this->lwork,this->work2.data_start(),info); - if (info) TRIQS_RUNTIME_ERROR<<"eigenelements_worker :error code zheev : "<mat; - this->has_run = true; - } - }; - - //-------------------------------- - - /** - * Simple diagonalization call, return all eigenelements. - * Handles both real and complex case. - * @param M : the matrix VIEW : it MUST be contiguous - * @param take_copy : makes a copy of the matrix before calling lapack, so that the original is preserved. - * if false : no copy is made and the content of the matrix M is destroyed (it is equal to vectors()). - * if true : a copy is made, M is preserved, but of course it is slower... - */ - template - std::pair, typename MatrixViewType::regular_type> eigenelements( MatrixViewType const & M, bool take_copy =false) { - eigenelements_worker W(take_copy ? MatrixViewType(make_clone(M)) : M); - W.invoke(); - return std::make_pair(W.values(),W.vectors()); - } - - //-------------------------------- - - /** - * Simple diagonalization call, returning only the eigenvalues. - * Handles both real and complex case. - * @param M : the matrix VIEW : it MUST be contiguous - * @param take_copy : makes a copy of the matrix before calling lapack, so that the original is preserved. - * if false : no copy is made and the content of the matrix M is destroyed. - * if true : a copy is made, M is preserved, but of course it is slower... - */ - template - triqs::arrays::vector_view eigenvalues( MatrixViewType const & M, bool take_copy = false) { - eigenelements_worker W(take_copy ? MatrixViewType(make_clone(M)) : M); W.invoke(); return W.values(); - } - -}}} // namespace triqs::arrays::linalg -#endif +} +} // namespace triqs::arrays::linalg diff --git a/triqs/lattice/tight_binding.cpp b/triqs/lattice/tight_binding.cpp index 5187529d..9f38042e 100644 --- a/triqs/lattice/tight_binding.cpp +++ b/triqs/lattice/tight_binding.cpp @@ -59,7 +59,7 @@ namespace lattice { array eval(norb, n_pts); k_t dk = (K2 - K1) / double(n_pts), k = K1; for (int i = 0; i < n_pts; ++i, k += dk) { - eval(range(), i) = linalg::eigenvalues(TK(k(range(0, ndim)))(), false); + eval(range(), i) = linalg::eigenvalues(TK(k(range(0, ndim)))()); } return eval; } @@ -86,7 +86,7 @@ namespace lattice { grid_generator grid(ndim, n_pts); array eval(norb, grid.size()); for (; grid; ++grid) { - eval(range(), grid.index()) = linalg::eigenvalues(TK((*grid)(range(0, ndim)))(), false); + eval(range(), grid.index()) = linalg::eigenvalues(TK((*grid)(range(0, ndim)))()); } return eval; }