ampsci
c++ program for high-precision atomic structure calculations of single-valence systems
Loading...
Searching...
No Matches
Solvers.ipp
1//==============================================================================
2// Implementations:
3
4extern "C" {
5
6using complex_double = long double;
7// https://gcc.gnu.org/onlinedocs/gcc-4.5.4/gfortran/ISO_005fC_005fBINDING.html
8// on apple M1, long double is just double
9// nb: it actually doesn't matter, since all operations happen on fortran end
10
11void dsyev_(char *, char *, int *, double *, int *, double *, double *, int *,
12 int *);
13void zheev_(char *, char *, int *, complex_double *, int *, double *,
14 complex_double *, int *, double *, int *);
15void dsygv_(int *, char *, char *, int *, double *, int *, double *, int *,
16 double *, double *, int *, int *);
17void zhegv_(int *, char *, char *, int *, complex_double *, int *,
18 complex_double *, int *, double *, complex_double *, int *,
19 double *, int *);
20
21void dsyevx_(char *, char *, char *, int *, double *, int *, double *, double *,
22 int *, int *, double *, int *, double *, double *, int *, double *,
23 int *, int *, int *, int *);
24}
25
26namespace LinAlg {
27
28//==============================================================================
29// Solves matrix equation Ax=b for x, for known square matrix A and vector b.
30template <typename T>
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>>");
36 assert(Am.rows() == b.size());
37
38 Vector<T> x(Am.cols());
39
40 auto Am_gsl = Am.as_gsl_view();
41 const auto b_gsl = b.as_gsl_view();
42 auto x_gsl = x.as_gsl_view();
43
44 int sLU = 0;
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,
53 &x_gsl.vector);
54 }
55 gsl_permutation_free(Am_perm);
56
57 return x;
58}
59
60//============================================================================*
61// Solves Av = ev for eigenvalues e and eigenvectors v of symmetric/Hermetian
62// matrix A. Returns [e,v], where v(i,j) is the jth element of the ith
63// eigenvector corresponding to ith eigenvalue, e(i). e is always real.
64template <typename T>
65std::pair<Vector<double>, Matrix<T>> symmhEigensystem(Matrix<T> A) {
66 assert(A.rows() == A.cols());
67
68 // E. values of Hermetian complex matrix are _real_
69 auto eigen_vv = std::make_pair(Vector<double>(A.rows()), std::move(A));
70 auto &[e_values, e_vectors] = eigen_vv;
71
72 int dim = (int)e_vectors.rows();
73 char jobz{'V'};
74 char uplo{'U'};
75 int info = 0;
76
77 std::vector<T> work(3 * e_vectors.rows());
78 int workspace_size = (int)work.size();
79
80 if constexpr (std::is_same_v<T, double>) {
81
82 // For double (symmetric real) matrix
83 dsyev_(&jobz, &uplo, &dim, e_vectors.data(), &dim, e_values.data(),
84 work.data(), &workspace_size, &info);
85
86 } else if constexpr (std::is_same_v<T, std::complex<double>>) {
87 // static_assert(sizeof(complex_double) == 16); // LAPACK assumption
88 // nb: this isn't actually required
89
90 // For complex double (Hermetian) matrix
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);
97 }
98
99 if (info != 0) {
100 std::cerr << "\nError 91: symmhEigensystem " << info << " " << dim << " "
101 << std::endl;
102 if (info < 0) {
103 std::cerr << "The " << -info << "-th argument had an illegal value\n";
104 } else {
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";
108 }
109 std::cerr << info << " " << A.size() << "\n";
110 }
111
112 return eigen_vv;
113}
114
115//============================================================================*
116template <typename T>
117std::pair<Vector<double>, Matrix<T>> symmhEigensystem(Matrix<T> A, int number) {
118 assert(A.rows() == A.cols());
119 assert(number < (int)A.rows());
120
121 // E. values of Hermetian complex matrix are _real_
122 auto eigen_vv = std::make_pair(Vector<double>(A.rows()),
123 Matrix<T>(std::size_t(A.rows())));
124 auto &[e_values, e_vectors] = eigen_vv;
125
126 int dim = (int)A.rows();
127 char jobz{'V'};
128 char uplo{'U'};
129 char range{'I'};
130 double vl{0.0};
131 double vu{0.0};
132 int il{1};
133 int iu{number};
134 double abstol{1.0e-12}; //?
135 int m{0};
136 int info{0};
137
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();
142
143 // For double (symmetric real) matrix
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);
147
148 if (info != 0) {
149 std::cerr << "\nError 135: symmhEigensystem " << info << " " << dim << " "
150 << std::endl;
151 }
152
153 return eigen_vv;
154}
155
156//==============================================================================
157// Solves Av = eBv for eigenvalues e and eigenvectors v of symmetric/Hermetian
158// matrix pair A,B. Returns [e,v], where v(i,j) is the jth element of the ith
159// eigenvector corresponding to ith eigenvalue, e(i). e is always real.
160template <typename T>
161std::pair<Vector<double>, Matrix<T>> symmhEigensystem(Matrix<T> A,
162 Matrix<T> B) {
163 assert(A.rows() == A.cols());
164 assert(B.rows() == B.cols());
165 assert(A.rows() == B.rows());
166
167 // E. values of Hermetian complex matrix are _real_
168 auto eigen_vv = std::make_pair(Vector<double>(A.rows()), std::move(A));
169 auto &[e_values, e_vectors] = eigen_vv;
170
171 int itype = 1;
172 int dim = (int)e_vectors.rows();
173 char jobz{'V'};
174 char uplo{'U'};
175 int info = 0;
176
177 std::vector<T> work(3 * e_vectors.rows());
178 int workspace_size = (int)work.size();
179
180 if constexpr (std::is_same_v<T, double>) {
181
182 // For double (symmetric real) matrix
183 dsygv_(&itype, &jobz, &uplo, &dim, e_vectors.data(), &dim, B.data(), &dim,
184 e_values.data(), work.data(), &workspace_size, &info);
185
186 } else if constexpr (std::is_same_v<T, std::complex<double>>) {
187 // static_assert(sizeof(complex_double) == 16); // LAPACK assumption
188 // nb: this isn't actually required
189
190 // For complex double (Hermetian) matrix
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(); // ?? why?
199 }
200
201 if (info != 0) {
202 std::cerr << "\nError 198: symmhEigensystem " << info << " " << dim << " "
203 << B.size() << std::endl;
204 if (info < 0) {
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;";
210 } else {
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.";
215 }
216 }
217
218 return eigen_vv;
219}
220
221//==============================================================================
222// Solves Av = ev for eigenvalues e and eigenvectors v of non-symmetric real
223// matrix A. Returns [e,v], where v(i,j) is the jth element of the ith
224// eigenvector corresponding to ith eigenvalue, e(i). A must be real, while e
225// and v are complex.
226template <typename T>
227std::pair<Vector<std::complex<double>>, Matrix<std::complex<double>>>
229 assert(A.rows() == A.cols());
230 static_assert(std::is_same_v<T, double>,
231 "genEigensystem only for Real matrix");
232
233 // Solves Av = ev for eigenvalues e and eigenvectors v
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;
238
239 // GSL: This function computes eigenvalues and right eigenvectors of the
240 // n-by-n real nonsymmetric matrix A.
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);
247 if (sort)
248 gsl_eigen_nonsymmv_sort(&val_gsl.vector, &vec_gsl.matrix,
249 GSL_EIGEN_SORT_ABS_ASC);
250
251 // Eigensystems using GSL. NOTE: e-vectors are stored in COLUMNS (not rows) of
252 // matrix! Therefore, we transpose the matrix (duuumb)
253 auto tmp = e_vectors.transpose();
254 e_vectors = std::move(tmp);
255
256 return eigen_vv;
257}
258
259} // namespace LinAlg
Matrix class; row-major.
Definition Matrix.hpp:35
std::size_t size() const
Return rows*columns [total array size].
Definition Matrix.hpp:93
std::size_t rows() const
Return rows [major index size].
Definition Matrix.hpp:89
auto as_gsl_view()
Returns gsl_matrix_view (or _float_view, _complex_view, _complex_float_view). Call ....
Definition Matrix.ipp:310
std::size_t cols() const
Return columns [minor index size].
Definition Matrix.hpp:91
Vector class (inherits from Matrix)
Definition Vector.hpp:9
auto as_gsl_view()
Returns gsl_vector_view (or _float_view, _complex_view, _complex_float_view). Call ....
Definition Vector.ipp:79
Defines Matrix, Vector classes, and linear some algebra functions.
Definition Matrix.hpp:26
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