32 static_assert(std::is_same_v<T, double> ||
33 std::is_same_v<T, std::complex<double>>,
34 "solve_Axeqb only available for Matrix<double> or "
35 "Matrix<complex<double>>");
42 auto x_gsl = x.as_gsl_view();
45 gsl_permutation *Am_perm =
46 gsl_permutation_alloc(std::max(Am.
rows(), Am.
cols()));
47 if constexpr (std::is_same_v<T, double>) {
48 gsl_linalg_LU_decomp(&Am_gsl.matrix, Am_perm, &sLU);
49 gsl_linalg_LU_solve(&Am_gsl.matrix, Am_perm, &b_gsl.vector, &x_gsl.vector);
50 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
51 gsl_linalg_complex_LU_decomp(&Am_gsl.matrix, Am_perm, &sLU);
52 gsl_linalg_complex_LU_solve(&Am_gsl.matrix, Am_perm, &b_gsl.vector,
55 gsl_permutation_free(Am_perm);
66 assert(A.rows() == A.cols());
69 auto eigen_vv = std::make_pair(
Vector<double>(A.rows()), std::move(A));
70 auto &[e_values, e_vectors] = eigen_vv;
72 int dim = (int)e_vectors.rows();
77 std::vector<T> work(3 * e_vectors.rows());
78 int workspace_size = (int)work.size();
80 if constexpr (std::is_same_v<T, double>) {
83 dsyev_(&jobz, &uplo, &dim, e_vectors.data(), &dim, e_values.data(),
84 work.data(), &workspace_size, &info);
86 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
91 std::vector<double> Rwork(3 * e_vectors.rows() - 2ul);
92 complex_double *evec_ptr =
93 reinterpret_cast<complex_double *
>(e_vectors.data());
94 complex_double *work_ptr =
reinterpret_cast<complex_double *
>(work.data());
95 zheev_(&jobz, &uplo, &dim, evec_ptr, &dim, e_values.data(), work_ptr,
96 &workspace_size, Rwork.data(), &info);
100 std::cerr <<
"\nError 91: symmhEigensystem " << info <<
" " << dim <<
" "
103 std::cerr <<
"The " << -info <<
"-th argument had an illegal value\n";
105 std::cerr <<
"The algorithm failed to converge; " << info
106 <<
" off-diagonal elements of an intermediate "
107 "tridiagonal form did not converge to zero.\n";
109 std::cerr << info <<
" " << A.size() <<
"\n";
118 assert(A.rows() == A.cols());
119 assert(number < (
int)A.rows());
124 auto &[e_values, e_vectors] = eigen_vv;
126 int dim = (int)A.rows();
134 double abstol{1.0e-12};
138 std::vector<T> work(8 * A.rows());
139 std::vector<int> iwork(5 * A.rows());
140 std::vector<int> ifail(A.rows());
141 int workspace_size = (int)work.size();
144 dsyevx_(&jobz, &range, &uplo, &dim, A.data(), &dim, &vl, &vu, &il, &iu,
145 &abstol, &m, e_values.data(), e_vectors.data(), &dim, work.data(),
146 &workspace_size, iwork.data(), ifail.data(), &info);
149 std::cerr <<
"\nError 135: symmhEigensystem " << info <<
" " << dim <<
" "
163 assert(A.rows() == A.cols());
164 assert(B.rows() == B.cols());
165 assert(A.rows() == B.rows());
168 auto eigen_vv = std::make_pair(
Vector<double>(A.rows()), std::move(A));
169 auto &[e_values, e_vectors] = eigen_vv;
172 int dim = (int)e_vectors.rows();
177 std::vector<T> work(3 * e_vectors.rows());
178 int workspace_size = (int)work.size();
180 if constexpr (std::is_same_v<T, double>) {
183 dsygv_(&itype, &jobz, &uplo, &dim, e_vectors.data(), &dim, B.data(), &dim,
184 e_values.data(), work.data(), &workspace_size, &info);
186 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
191 std::vector<double> Rwork(3 * e_vectors.rows() - 2ul);
192 complex_double *evec_ptr =
193 reinterpret_cast<complex_double *
>(e_vectors.data());
194 complex_double *Bmat_ptr =
reinterpret_cast<complex_double *
>(B.data());
195 complex_double *work_ptr =
reinterpret_cast<complex_double *
>(work.data());
196 zhegv_(&itype, &jobz, &uplo, &dim, evec_ptr, &dim, Bmat_ptr, &dim,
197 e_values.data(), work_ptr, &workspace_size, Rwork.data(), &info);
198 e_vectors.conj_in_place();
202 std::cerr <<
"\nError 198: symmhEigensystem " << info <<
" " << dim <<
" "
203 << B.size() << std::endl;
205 std::cerr <<
"The " << -info <<
"-th argument had an illegal value\n";
206 }
else if (info <= dim) {
207 std::cerr <<
"DSYEV failed to converge; " << info
208 <<
" off-diagonal elements of an intermediate tridiagonal form "
209 "did not converge to zero;";
211 std::cerr <<
"The leading minor of order " << info
212 <<
" of B is not positive definite. "
213 "The factorization of B could not be completed and no "
214 "eigenvalues or eigenvectors were computed.";
229 assert(A.rows() == A.cols());
230 static_assert(std::is_same_v<T, double>,
231 "genEigensystem only for Real matrix");
234 const auto dim = A.rows();
235 auto eigen_vv = std::make_pair(
Vector<std::complex<double>>(dim),
236 Matrix<std::complex<double>>(dim));
237 auto &[e_values, e_vectors] = eigen_vv;
241 auto A_gsl = A.as_gsl_view();
242 auto val_gsl = e_values.as_gsl_view();
243 auto vec_gsl = e_vectors.as_gsl_view();
244 gsl_eigen_nonsymmv_workspace *work = gsl_eigen_nonsymmv_alloc(dim);
245 gsl_eigen_nonsymmv(&A_gsl.matrix, &val_gsl.vector, &vec_gsl.matrix, work);
246 gsl_eigen_nonsymmv_free(work);
248 gsl_eigen_nonsymmv_sort(&val_gsl.vector, &vec_gsl.matrix,
249 GSL_EIGEN_SORT_ABS_ASC);
253 auto tmp = e_vectors.transpose();
254 e_vectors = std::move(tmp);
auto as_gsl_view()
Returns gsl_matrix_view (or _float_view, _complex_view, _complex_float_view). Call ....
Definition Matrix.ipp:310
std::pair< Vector< std::complex< double > >, Matrix< std::complex< double > > > genEigensystem(Matrix< T > A, bool sort)
Solves Av = ev for eigenvalues e and eigenvectors v of non-symmetric real matrix A....
Definition Solvers.ipp:228
Vector< T > solve_Axeqb(Matrix< T > Am, const Vector< T > &b)
Solves matrix equation Ax=b for x, for known square matrix A and vector b.
Definition Solvers.ipp:31
std::pair< Vector< double >, Matrix< T > > symmhEigensystem(Matrix< T > A)
Solves Av = ev for eigenvalues e and eigenvectors v of symmetric/Hermetian matrix A....
Definition Solvers.ipp:65