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
25void dsyevr_(char *jobz, char *range, char *uplo, int *n, double *a, int *lda,
26 double *vl, double *vu, int *il, int *iu, double *abstol, int *m,
27 double *w, double *z, int *ldz, int *isuppz, double *work,
28 int *lwork, int *iwork, int *liwork, int *info);
29}
30
31namespace LinAlg {
32
33//==============================================================================
34// Solves matrix equation Ax=b for x, for known square matrix A and vector b.
35template <typename T>
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>>");
41 assert(Am.rows() == b.size());
42
43 Vector<T> x(Am.cols());
44
45 auto Am_gsl = Am.as_gsl_view();
46 const auto b_gsl = b.as_gsl_view();
47 auto x_gsl = x.as_gsl_view();
48
49 int sLU = 0;
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,
58 &x_gsl.vector);
59 }
60 gsl_permutation_free(Am_perm);
61
62 return x;
63}
64
65//============================================================================*
66// Solves Av = ev for eigenvalues e and eigenvectors v of symmetric/Hermetian
67// matrix A. Returns [e,v], where v(i,j) is the jth element of the ith
68// eigenvector corresponding to ith eigenvalue, e(i). e is always real.
69template <typename T>
70std::pair<Vector<double>, Matrix<T>> symmhEigensystem(Matrix<T> A) {
71 assert(A.rows() == A.cols());
72
73 // E. values of Hermetian complex matrix are _real_
74 auto eigen_vv = std::make_pair(Vector<double>(A.rows()), std::move(A));
75 auto &[e_values, e_vectors] = eigen_vv;
76
77 int dim = (int)e_vectors.rows();
78 char jobz{'V'};
79 char uplo{'U'};
80 int info = 0;
81
82 std::vector<T> work(3 * e_vectors.rows());
83 int workspace_size = (int)work.size();
84
85 if constexpr (std::is_same_v<T, double>) {
86
87 // For double (symmetric real) matrix
88 dsyev_(&jobz, &uplo, &dim, e_vectors.data(), &dim, e_values.data(),
89 work.data(), &workspace_size, &info);
90
91 } else if constexpr (std::is_same_v<T, std::complex<double>>) {
92 // static_assert(sizeof(complex_double) == 16); // LAPACK assumption
93 // nb: this isn't actually required
94
95 // For complex double (Hermetian) matrix
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);
102 }
103
104 if (info != 0) {
105 std::cerr << "\nError 91: symmhEigensystem " << info << " " << dim << " "
106 << std::endl;
107 if (info < 0) {
108 std::cerr << "The " << -info << "-th argument had an illegal value\n";
109 } else {
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";
113 }
114 std::cerr << info << " " << A.size() << "\n";
115 }
116
117 return eigen_vv;
118}
119
120//============================================================================*
121template <typename T>
122std::pair<Vector<double>, Matrix<T>> symmhEigensystem(Matrix<T> A, int number) {
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?");
128
129 // E. values of Hermetian complex matrix are _real_
130 auto eigen_vv = std::make_pair(Vector<double>(A.rows()), Matrix<T>(A.rows()));
131 auto &[e_values, e_vectors] = eigen_vv;
132
133 int dim = (int)A.rows();
134 char jobz{'V'};
135 char uplo{'U'};
136 char range{'I'};
137 double vl{0.0};
138 double vu{0.0};
139 int il{1};
140 int iu{number};
141 double abstol{1.0e-12}; //?
142 int m{0};
143 int info{0};
144
145 // std::vector<T> work(8 * A.rows());
146 // std::vector<int> iwork(5 * A.rows());
147 // std::vector<int> ifail(A.rows());
148
149 // ---- Workspace query (DSYEVX)
150 int lwork = -1;
151 double work_query = 0.0;
152 std::vector<int> iwork(5 * A.rows()); // DSYEVX requires at least 5*N
153 std::vector<int> ifail(A.rows()); // size N
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);
157
158 if (info != 0) {
159 std::cerr << "\nError: dsyevx workspace query info=" << info << "\n";
160 }
161
162 // ---- Allocate optimal WORK
163
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();
167
168 // For double (symmetric real) matrix
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);
172
173 if (info != 0) {
174 std::cerr << "\nError 135: symmhEigensystem " << info << " " << dim << " "
175 << std::endl;
176 }
177
178 return eigen_vv;
179}
180
181//============================================================================*
182template <typename T>
183std::tuple<int, Vector<double>, Matrix<T>> symmhEigensystem(Matrix<T> A,
184 double all_below) {
185 static_assert(std::is_same_v<T, double>,
186 "This overload uses dsyevr (real symmetric). "
187 "For complex Hermitian use zheevr/cheevr.");
188
189 assert(A.rows() == A.cols());
190 int n = static_cast<int>(A.rows());
191
192 // Return containers (same idea as your code: values size n, vectors n x n)
193 auto eigen_vv =
194 std::make_tuple(int{0}, Vector<double>(A.rows()), Matrix<T>(A.rows()));
195 auto &[num_solutions, e_values, e_vectors] = eigen_vv;
196
197 char jobz{'V'}; // eigenvectors
198 char range{'V'}; // value range
199 char uplo{'U'};
200
201 // dsyevr uses (vl,vu]; pick vl very negative so we get everything <= all_below
202 // Avoid +/-inf: many Fortran ABIs dislike it.
203 double vl = -std::numeric_limits<double>::max();
204 double vu = all_below;
205
206 int il{0}, iu{0}; // unused for RANGE='V'
207 double abstol{1.0e-12}; // your choice; could also use DLAMCH("S")
208 num_solutions = 0;
209 int info{0};
210
211 // Required by dsyevr
212 std::vector<int> isuppz(std::size_t(2 * std::max(1, n)));
213
214 // Workspace query
215 int lwork{-1}, liwork{-1};
216 double work_query{0.0};
217 int iwork_query{0};
218
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);
222
223 if (info != 0) {
224 std::cerr << "\nError: symmhEigensystem dsyevr(work query) info=" << info
225 << " n=" << n << "\n";
226 return eigen_vv;
227 }
228
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)));
233
234 // Actual call
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);
238
239 if (info != 0) {
240 std::cerr << "\nError: symmhEigensystem dsyevr info=" << info << " n=" << n
241 << " m=" << num_solutions << "\n";
242 }
243
244 return eigen_vv;
245}
246
247//==============================================================================
248// Solves Av = eBv for eigenvalues e and eigenvectors v of symmetric/Hermetian
249// matrix pair A,B. Returns [e,v], where v(i,j) is the jth element of the ith
250// eigenvector corresponding to ith eigenvalue, e(i). e is always real.
251template <typename T>
252std::pair<Vector<double>, Matrix<T>> symmhEigensystem(Matrix<T> A,
253 Matrix<T> B) {
254 assert(A.rows() == A.cols());
255 assert(B.rows() == B.cols());
256 assert(A.rows() == B.rows());
257
258 // E. values of Hermetian complex matrix are _real_
259 auto eigen_vv = std::make_pair(Vector<double>(A.rows()), std::move(A));
260 auto &[e_values, e_vectors] = eigen_vv;
261
262 int itype = 1;
263 int dim = (int)e_vectors.rows();
264 char jobz{'V'};
265 char uplo{'U'};
266 int info = 0;
267
268 std::vector<T> work(3 * e_vectors.rows());
269 int workspace_size = (int)work.size();
270
271 if constexpr (std::is_same_v<T, double>) {
272
273 // For double (symmetric real) matrix
274 dsygv_(&itype, &jobz, &uplo, &dim, e_vectors.data(), &dim, B.data(), &dim,
275 e_values.data(), work.data(), &workspace_size, &info);
276
277 } else if constexpr (std::is_same_v<T, std::complex<double>>) {
278 // static_assert(sizeof(complex_double) == 16); // LAPACK assumption
279 // nb: this isn't actually required
280
281 // For complex double (Hermetian) matrix
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(); // ?? why?
290 }
291
292 if (info != 0) {
293 std::cerr << "\nError 198: symmhEigensystem " << info << " " << dim << " "
294 << B.size() << std::endl;
295 if (info < 0) {
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;";
301 } else {
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.";
306 }
307 }
308
309 return eigen_vv;
310}
311
312//==============================================================================
313// Solves Av = ev for eigenvalues e and eigenvectors v of non-symmetric real
314// matrix A. Returns [e,v], where v(i,j) is the jth element of the ith
315// eigenvector corresponding to ith eigenvalue, e(i). A must be real, while e
316// and v are complex.
317template <typename T>
318std::pair<Vector<std::complex<double>>, Matrix<std::complex<double>>>
320 assert(A.rows() == A.cols());
321 static_assert(std::is_same_v<T, double>,
322 "genEigensystem only for Real matrix");
323
324 // Solves Av = ev for eigenvalues e and eigenvectors v
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;
329
330 // GSL: This function computes eigenvalues and right eigenvectors of the
331 // n-by-n real nonsymmetric matrix A.
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);
338 if (sort)
339 gsl_eigen_nonsymmv_sort(&val_gsl.vector, &vec_gsl.matrix,
340 GSL_EIGEN_SORT_ABS_ASC);
341
342 // Eigensystems using GSL. NOTE: e-vectors are stored in COLUMNS (not rows) of
343 // matrix! Therefore, we transpose the matrix (duuumb)
344 auto tmp = e_vectors.transpose();
345 e_vectors = std::move(tmp);
346
347 return eigen_vv;
348}
349
350} // namespace LinAlg
Matrix class; row-major.
Definition Matrix.hpp:35
std::size_t size() const
Return rows*columns [total array size].
Definition Matrix.hpp:112
std::size_t rows() const
Return rows [major index size].
Definition Matrix.hpp:108
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:110
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: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