37 static_assert(std::is_same_v<T, double> ||
38 std::is_same_v<T, std::complex<double>>,
39 "solve_Axeqb only available for Matrix<double> or "
40 "Matrix<complex<double>>");
47 auto x_gsl = x.as_gsl_view();
50 gsl_permutation *Am_perm =
51 gsl_permutation_alloc(std::max(Am.
rows(), Am.
cols()));
52 if constexpr (std::is_same_v<T, double>) {
53 gsl_linalg_LU_decomp(&Am_gsl.matrix, Am_perm, &sLU);
54 gsl_linalg_LU_solve(&Am_gsl.matrix, Am_perm, &b_gsl.vector, &x_gsl.vector);
55 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
56 gsl_linalg_complex_LU_decomp(&Am_gsl.matrix, Am_perm, &sLU);
57 gsl_linalg_complex_LU_solve(&Am_gsl.matrix, Am_perm, &b_gsl.vector,
60 gsl_permutation_free(Am_perm);
71 assert(A.rows() == A.cols());
74 auto eigen_vv = std::make_pair(
Vector<double>(A.rows()), std::move(A));
75 auto &[e_values, e_vectors] = eigen_vv;
77 int dim = (int)e_vectors.rows();
82 std::vector<T> work(3 * e_vectors.rows());
83 int workspace_size = (int)work.size();
85 if constexpr (std::is_same_v<T, double>) {
88 dsyev_(&jobz, &uplo, &dim, e_vectors.data(), &dim, e_values.data(),
89 work.data(), &workspace_size, &info);
91 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
96 std::vector<double> Rwork(3 * e_vectors.rows() - 2ul);
97 complex_double *evec_ptr =
98 reinterpret_cast<complex_double *
>(e_vectors.data());
99 complex_double *work_ptr =
reinterpret_cast<complex_double *
>(work.data());
100 zheev_(&jobz, &uplo, &dim, evec_ptr, &dim, e_values.data(), work_ptr,
101 &workspace_size, Rwork.data(), &info);
105 std::cerr <<
"\nError 91: symmhEigensystem " << info <<
" " << dim <<
" "
108 std::cerr <<
"The " << -info <<
"-th argument had an illegal value\n";
110 std::cerr <<
"The algorithm failed to converge; " << info
111 <<
" off-diagonal elements of an intermediate "
112 "tridiagonal form did not converge to zero.\n";
114 std::cerr << info <<
" " << A.size() <<
"\n";
123 assert(A.rows() == A.cols());
124 assert(number < (
int)A.rows());
125 static_assert(std::is_same_v<T, double>,
126 "This overload uses dsyevx_ (real symmetric). "
127 "For complex Hermitian use zheevx/cheevx?");
131 auto &[e_values, e_vectors] = eigen_vv;
133 int dim = (int)A.rows();
141 double abstol{1.0e-12};
151 double work_query = 0.0;
152 std::vector<int> iwork(5 * A.rows());
153 std::vector<int> ifail(A.rows());
154 dsyevx_(&jobz, &range, &uplo, &dim, A.data(), &dim, &vl, &vu, &il, &iu,
155 &abstol, &m, e_values.data(), e_vectors.data(), &dim, &work_query,
156 &lwork, iwork.data(), ifail.data(), &info);
159 std::cerr <<
"\nError: dsyevx workspace query info=" << info <<
"\n";
164 lwork =
static_cast<int>(work_query);
165 std::vector<double> work(std::size_t(std::max(1, lwork)));
166 int workspace_size = (int)work.size();
169 dsyevx_(&jobz, &range, &uplo, &dim, A.data(), &dim, &vl, &vu, &il, &iu,
170 &abstol, &m, e_values.data(), e_vectors.data(), &dim, work.data(),
171 &workspace_size, iwork.data(), ifail.data(), &info);
174 std::cerr <<
"\nError 135: symmhEigensystem " << info <<
" " << dim <<
" "
185 static_assert(std::is_same_v<T, double>,
186 "This overload uses dsyevr (real symmetric). "
187 "For complex Hermitian use zheevr/cheevr.");
189 assert(A.rows() == A.cols());
190 int n =
static_cast<int>(A.rows());
195 auto &[num_solutions, e_values, e_vectors] = eigen_vv;
203 double vl = -std::numeric_limits<double>::max();
204 double vu = all_below;
207 double abstol{1.0e-12};
212 std::vector<int> isuppz(std::size_t(2 * std::max(1, n)));
215 int lwork{-1}, liwork{-1};
216 double work_query{0.0};
219 dsyevr_(&jobz, &range, &uplo, &n, A.data(), &n, &vl, &vu, &il, &iu, &abstol,
220 &num_solutions, e_values.data(), e_vectors.data(), &n, isuppz.data(),
221 &work_query, &lwork, &iwork_query, &liwork, &info);
224 std::cerr <<
"\nError: symmhEigensystem dsyevr(work query) info=" << info
225 <<
" n=" << n <<
"\n";
229 lwork =
static_cast<int>(work_query);
230 liwork = iwork_query;
231 std::vector<double> work(std::size_t(std::max(1, lwork)));
232 std::vector<int> iwork(std::size_t(std::max(1, liwork)));
235 dsyevr_(&jobz, &range, &uplo, &n, A.data(), &n, &vl, &vu, &il, &iu, &abstol,
236 &num_solutions, e_values.data(), e_vectors.data(), &n, isuppz.data(),
237 work.data(), &lwork, iwork.data(), &liwork, &info);
240 std::cerr <<
"\nError: symmhEigensystem dsyevr info=" << info <<
" n=" << n
241 <<
" m=" << num_solutions <<
"\n";
254 assert(A.rows() == A.cols());
255 assert(B.rows() == B.cols());
256 assert(A.rows() == B.rows());
259 auto eigen_vv = std::make_pair(
Vector<double>(A.rows()), std::move(A));
260 auto &[e_values, e_vectors] = eigen_vv;
263 int dim = (int)e_vectors.rows();
268 std::vector<T> work(3 * e_vectors.rows());
269 int workspace_size = (int)work.size();
271 if constexpr (std::is_same_v<T, double>) {
274 dsygv_(&itype, &jobz, &uplo, &dim, e_vectors.data(), &dim, B.data(), &dim,
275 e_values.data(), work.data(), &workspace_size, &info);
277 }
else if constexpr (std::is_same_v<T, std::complex<double>>) {
282 std::vector<double> Rwork(3 * e_vectors.rows() - 2ul);
283 complex_double *evec_ptr =
284 reinterpret_cast<complex_double *
>(e_vectors.data());
285 complex_double *Bmat_ptr =
reinterpret_cast<complex_double *
>(B.data());
286 complex_double *work_ptr =
reinterpret_cast<complex_double *
>(work.data());
287 zhegv_(&itype, &jobz, &uplo, &dim, evec_ptr, &dim, Bmat_ptr, &dim,
288 e_values.data(), work_ptr, &workspace_size, Rwork.data(), &info);
289 e_vectors.conj_in_place();
293 std::cerr <<
"\nError 198: symmhEigensystem " << info <<
" " << dim <<
" "
294 << B.size() << std::endl;
296 std::cerr <<
"The " << -info <<
"-th argument had an illegal value\n";
297 }
else if (info <= dim) {
298 std::cerr <<
"DSYEV failed to converge; " << info
299 <<
" off-diagonal elements of an intermediate tridiagonal form "
300 "did not converge to zero;";
302 std::cerr <<
"The leading minor of order " << info
303 <<
" of B is not positive definite. "
304 "The factorization of B could not be completed and no "
305 "eigenvalues or eigenvectors were computed.";
320 assert(A.rows() == A.cols());
321 static_assert(std::is_same_v<T, double>,
322 "genEigensystem only for Real matrix");
325 const auto dim = A.rows();
326 auto eigen_vv = std::make_pair(
Vector<std::complex<double>>(dim),
327 Matrix<std::complex<double>>(dim));
328 auto &[e_values, e_vectors] = eigen_vv;
332 auto A_gsl = A.as_gsl_view();
333 auto val_gsl = e_values.as_gsl_view();
334 auto vec_gsl = e_vectors.as_gsl_view();
335 gsl_eigen_nonsymmv_workspace *work = gsl_eigen_nonsymmv_alloc(dim);
336 gsl_eigen_nonsymmv(&A_gsl.matrix, &val_gsl.vector, &vec_gsl.matrix, work);
337 gsl_eigen_nonsymmv_free(work);
339 gsl_eigen_nonsymmv_sort(&val_gsl.vector, &vec_gsl.matrix,
340 GSL_EIGEN_SORT_ABS_ASC);
344 auto tmp = e_vectors.transpose();
345 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:319
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:36
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:70