pairinteraction
A Rydberg Interaction Calculator
SystemAtom.test.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2024 Pairinteraction Developers
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
5
17
18#include <Eigen/Eigenvalues>
19#include <doctest/doctest.h>
20#include <fmt/ranges.h>
21
22namespace pairinteraction {
23
24constexpr double VOLT_PER_CM_IN_ATOMIC_UNITS = 1 / 5.14220675112e9;
25constexpr double UM_IN_ATOMIC_UNITS = 1 / 5.29177210544e-5;
26constexpr double HARTREE_IN_GHZ = 6579683.920501762;
27
28DOCTEST_TEST_CASE("construct and diagonalize a small Hamiltonian") {
29 auto &database = Database::get_global_instance();
30 auto diagonalizer = DiagonalizerEigen<double>();
31
32 auto ket1 = KetAtomCreator()
33 .set_species("Rb")
38 .create(database);
39 auto ket2 = KetAtomCreator()
40 .set_species("Rb")
45 .create(database);
46 auto basis = BasisAtomCreator<double>().append_ket(ket1).append_ket(ket2).create(database);
47
48 auto system = SystemAtom<double>(basis);
49 system.set_electric_field({0, 0, 0.0001});
50
51 Eigen::MatrixXd tmp = Eigen::MatrixXd(1e5 * system.get_matrix()).array().round() / 1e5;
52 std::vector<double> matrix_vector(tmp.data(), tmp.data() + tmp.size());
53 DOCTEST_MESSAGE(fmt::format("Constructed: {}", fmt::join(matrix_vector, ", ")));
54
55 system.diagonalize(diagonalizer);
56 tmp = Eigen::MatrixXd(1e5 * system.get_matrix()).array().round() / 1e5;
57 matrix_vector = std::vector<double>(tmp.data(), tmp.data() + tmp.size());
58 DOCTEST_MESSAGE(fmt::format("Diagonalized: {}", fmt::join(matrix_vector, ", ")));
59
60 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eigensolver;
61 eigensolver.compute(system.get_matrix());
62 auto eigenenergies_eigen = eigensolver.eigenvalues();
63 auto eigenenergies_pairinteraction = system.get_eigenenergies();
64 for (int i = 0; i < eigenenergies_eigen.size(); ++i) {
65 DOCTEST_CHECK(std::abs(eigenenergies_eigen(i) - eigenenergies_pairinteraction(i)) < 1e-10);
66 }
67}
68
69DOCTEST_TEST_CASE("construct and diagonalize two Hamiltonians in parallel") {
70 auto &database = Database::get_global_instance();
71 auto diagonalizer = DiagonalizerEigen<std::complex<double>>();
72
74 .set_species("Rb")
76 .restrict_quantum_number_l(0, 1)
77 .create(database);
78
79 auto system1 = SystemAtom<std::complex<double>>(basis);
80 system1.set_electric_field({0, 0, 0.0001});
81
82 auto system2 = SystemAtom<std::complex<double>>(basis);
83 system2.set_electric_field({0, 0, 0.0002});
84
85 diagonalize<SystemAtom<std::complex<double>>>({system1, system2}, diagonalizer);
86
87 auto matrix1 = system1.get_matrix();
88 auto matrix2 = system2.get_matrix();
89 for (int i = 0; i < matrix1.rows(); ++i) {
90 for (int j = 0; j < matrix1.cols(); ++j) {
91 if (i != j) {
92 DOCTEST_CHECK(std::abs(matrix1.coeff(i, j)) < 1e-10);
93 DOCTEST_CHECK(std::abs(matrix2.coeff(i, j)) < 1e-10);
94 }
95 }
96 }
97}
98
99DOCTEST_TEST_CASE("construct and diagonalize multiple Hamiltonians in parallel" *
100 doctest::skip(true)) {
101 // TODO For a slow database, the fast parallelized construction of the tiny Hamiltonians seems
102 // to lead to the wrong Hamiltonians (visible in the warning "The floating point error (5e-324
103 // Hartree) is similar or larger than error estimated from the specified tolerance (0
104 // Hartree)."). This could be caused by an issue with thread safety or memory access.
105 int n = 10;
106
107 auto &database = Database::get_global_instance();
108 auto diagonalizer = DiagonalizerEigen<std::complex<double>>();
109
111 .set_species("Sr87_mqdt")
113 .restrict_quantum_number_l(0, 1)
114 .create(database);
115
116 std::vector<SystemAtom<std::complex<double>>> systems;
117 systems.reserve(n);
118 for (int i = 0; i < n; ++i) {
119 auto system = SystemAtom<std::complex<double>>(basis);
120 system.set_electric_field({0, 0, 0.0001 * i});
121 systems.push_back(std::move(system));
122 }
123
124 DOCTEST_MESSAGE("Basis size: ", basis->get_number_of_states());
125
126 diagonalize<SystemAtom<std::complex<double>>>(systems, diagonalizer);
127}
128
129DOCTEST_TEST_CASE("construct and diagonalize a Hamiltonian using different methods") {
130 auto &database = Database::get_global_instance();
131
133 .set_species("Rb")
135 .restrict_quantum_number_l(0, 1)
136 .create(database);
137
138 // Diagonalize using the Eigen library
139 auto system = SystemAtom<std::complex<double>>(basis);
142
143 Eigen::MatrixXcd matrix = system.get_matrix();
144 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXcd> eigensolver;
145 eigensolver.compute(matrix);
146 auto eigenenergies_eigen = eigensolver.eigenvalues();
147
148 // Create diagonalizers
149 std::vector<std::unique_ptr<DiagonalizerInterface<std::complex<double>>>> diagonalizers;
150 std::vector<double> rtols;
151 double eps{};
152 DOCTEST_SUBCASE("Double precision") {
153 diagonalizers.push_back(std::make_unique<DiagonalizerEigen<std::complex<double>>>());
154#ifdef WITH_LAPACKE
155 diagonalizers.push_back(std::make_unique<DiagonalizerLapackeEvd<std::complex<double>>>());
156 diagonalizers.push_back(std::make_unique<DiagonalizerLapackeEvr<std::complex<double>>>());
157#endif
158#ifdef WITH_MKL
159 diagonalizers.push_back(std::make_unique<DiagonalizerFeast<std::complex<double>>>(300));
160#endif
161 rtols = {1e-1, 1e-6, 1e-14};
162 eps = std::numeric_limits<double>::epsilon();
163 }
164
165 DOCTEST_SUBCASE("Single precision") {
166 diagonalizers.push_back(
167 std::make_unique<DiagonalizerEigen<std::complex<double>>>(FloatType::FLOAT32));
168#ifdef WITH_LAPACKE
169 diagonalizers.push_back(
170 std::make_unique<DiagonalizerLapackeEvd<std::complex<double>>>(FloatType::FLOAT32));
171 diagonalizers.push_back(
172 std::make_unique<DiagonalizerLapackeEvr<std::complex<double>>>(FloatType::FLOAT32));
173#endif
174#ifdef WITH_MKL
175 diagonalizers.push_back(
176 std::make_unique<DiagonalizerFeast<std::complex<double>>>(300, FloatType::FLOAT32));
177#endif
178 rtols = {1e-1, 1e-6};
179 eps = std::numeric_limits<float>::epsilon();
180 }
181
182 // Diagonalize using pairinteraction
183 for (double rtol_eigenenergies : rtols) {
184 double atol_eigenvectors =
185 std::max(0.5 * rtol_eigenenergies / std::sqrt(basis->get_number_of_states()), 4 * eps);
186 DOCTEST_MESSAGE("Precision: " << rtol_eigenenergies << " (rtol eigenenergies), "
187 << atol_eigenvectors << " (atol eigenvectors)");
188
189 for (const auto &diagonalizer : diagonalizers) {
190 auto system = SystemAtom<std::complex<double>>(basis);
194
195 // We specify a search interval because this is required if the FEAST routine is
196 // used. To avoid overflows, the interval ranges from half the smallest possible
197 // value to half the largest possible value.
198 system.diagonalize(*diagonalizer, std::numeric_limits<float>::lowest() / 2,
199 std::numeric_limits<float>::max() / 2, rtol_eigenenergies);
200 auto eigenenergies_pairinteraction = system.get_eigenenergies();
201 auto eigenvectors_pairinteraction = system.get_eigenbasis()->get_coefficients();
202
203 DOCTEST_CHECK(
204 (eigenenergies_eigen - eigenenergies_pairinteraction).array().abs().maxCoeff() <
205 rtol_eigenenergies * matrix.norm());
206
207 for (int i = 0; i < eigenvectors_pairinteraction.cols(); ++i) {
208 DOCTEST_CHECK(abs(1 - eigenvectors_pairinteraction.col(i).norm()) <
209 atol_eigenvectors * eigenvectors_pairinteraction.rows());
210 }
211 }
212 }
213}
214
215DOCTEST_TEST_CASE("construct and diagonalize a Hamiltonian with energy restrictions") {
216 double min_energy = 0.153355;
217 double max_energy = 0.153360;
218
219 auto &database = Database::get_global_instance();
220
221 auto basis = BasisAtomCreator<double>()
222 .set_species("Rb")
223 .restrict_quantum_number_n(58, 62)
224 .restrict_quantum_number_l(0, 1)
225 .create(database);
226
227 // Diagonalize using the Eigen library
228 auto system = SystemAtom<double>(basis);
229 system.set_electric_field(
231
232 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> eigensolver;
233 eigensolver.compute(system.get_matrix());
234 auto eigenenergies_all = eigensolver.eigenvalues();
235 std::vector<double> eigenenergies_eigen;
236 for (int i = 0; i < eigenenergies_all.size(); ++i) {
237 if (eigenenergies_all[i] > min_energy && eigenenergies_all[i] < max_energy) {
238 eigenenergies_eigen.push_back(eigenenergies_all[i]);
239 }
240 }
241
242 // Create diagonalizer
243 std::vector<std::unique_ptr<DiagonalizerInterface<double>>> diagonalizers;
244 diagonalizers.push_back(std::make_unique<DiagonalizerEigen<double>>(FloatType::FLOAT64));
245#ifdef WITH_LAPACKE
246 diagonalizers.push_back(std::make_unique<DiagonalizerLapackeEvd<double>>(FloatType::FLOAT64));
247 diagonalizers.push_back(std::make_unique<DiagonalizerLapackeEvr<double>>(FloatType::FLOAT64));
248#endif
249#ifdef WITH_MKL
250 diagonalizers.push_back(std::make_unique<DiagonalizerFeast<double>>(10, FloatType::FLOAT64));
251#endif
252
253 // Diagonalize using pairinteraction
254 for (const auto &diagonalizer : diagonalizers) {
255 auto system = SystemAtom<double>(basis);
256 system.set_electric_field(
258
259 system.diagonalize(*diagonalizer, min_energy, max_energy, 1e-6);
260 auto eigenenergies_pairinteraction = system.get_eigenenergies();
261
262 Eigen::MatrixXd tmp = (1e5 * eigenenergies_pairinteraction).array().round() / 1e5;
263 std::vector<double> eigenenergies_vector(tmp.data(), tmp.data() + tmp.size());
264 DOCTEST_MESSAGE(fmt::format("Eigenenergies: {}", fmt::join(eigenenergies_vector, ", ")));
265
266 DOCTEST_CHECK(eigenenergies_eigen.size() == 8);
267 DOCTEST_CHECK(eigenenergies_pairinteraction.size() == 8);
268 for (size_t i = 0; i < eigenenergies_eigen.size(); ++i) {
269 DOCTEST_CHECK(std::abs(eigenenergies_eigen[i] - eigenenergies_pairinteraction[i]) <
270 1e-10);
271 }
272 }
273}
274
275#ifdef WITH_MKL
276#include <Eigen/Dense>
277#include <mkl.h>
278DOCTEST_TEST_CASE("diagonalization with mkl") {
279 // We loop several times to check for memory errors
280 for (size_t i = 0; i < 10; ++i) {
281 // Create a symmetric matrix
282 int n = 100;
283 Eigen::MatrixXd matrix = Eigen::MatrixXd::Random(n, n);
284 matrix = (matrix + matrix.transpose()).eval();
285 Eigen::VectorXd eigenenergies(n);
286
287 // Diagonalize the matrix
288 int info =
289 LAPACKE_dsyev(LAPACK_COL_MAJOR, 'V', 'U', n, matrix.data(), n, eigenenergies.data());
290 DOCTEST_CHECK(info == 0);
291 }
292}
293#endif
294
295DOCTEST_TEST_CASE("handle it gracefully if no eigenenergies are within energy restrictions") {
296 double min_energy = -1;
297 double max_energy = -1;
298
299 auto &database = Database::get_global_instance();
300
301 auto basis = BasisAtomCreator<double>()
302 .set_species("Rb")
303 .restrict_quantum_number_n(58, 62)
304 .restrict_quantum_number_l(0, 1)
305 .create(database);
306
307 std::vector<std::unique_ptr<DiagonalizerInterface<double>>> diagonalizers;
308 diagonalizers.push_back(std::make_unique<DiagonalizerEigen<double>>());
309#ifdef WITH_LAPACKE
310 diagonalizers.push_back(std::make_unique<DiagonalizerLapackeEvd<double>>());
311#endif
312
313 for (const auto &diagonalizer : diagonalizers) {
314 auto system = SystemAtom<double>(basis);
315 system.set_electric_field(
317
318 system.diagonalize(*diagonalizer, min_energy, max_energy, 1e-6);
319 auto eigenenergies_pairinteraction = system.get_eigenenergies();
320
321 DOCTEST_CHECK(eigenenergies_pairinteraction.size() == 0);
322 }
323}
324
325DOCTEST_TEST_CASE("atom ion interaction") {
326 auto &database = Database::get_global_instance();
327 DiagonalizerEigen<double> diagonalizer;
328
329 auto ket = KetAtomCreator()
330 .set_species("Rb")
335 .create(database);
336 double energy = ket->get_energy();
337 double min_energy = energy - 50 / HARTREE_IN_GHZ;
338 double max_energy = energy + 50 / HARTREE_IN_GHZ;
339
340 auto basis = BasisAtomCreator<double>()
341 .set_species("Rb")
342 .restrict_quantum_number_n(58, 62)
343 .restrict_quantum_number_l(0, 3)
344 .restrict_quantum_number_m(0.5, 0.5)
345 .create(database);
346
347 auto system3 = SystemAtom<double>(basis);
348 system3.set_ion_interaction_order(3);
349 system3.set_ion_distance_vector({0, 0, 10 * UM_IN_ATOMIC_UNITS});
350 system3.diagonalize(diagonalizer, min_energy, max_energy, 1e-6);
351 auto energies3 = system3.get_eigenenergies();
352
353 auto system2 = SystemAtom<double>(basis);
354 system2.set_ion_interaction_order(2);
355 system2.set_ion_distance_vector({0, 0, 10 * UM_IN_ATOMIC_UNITS});
356 system2.diagonalize(diagonalizer, min_energy, max_energy, 1e-6);
357 auto energies2 = system2.get_eigenenergies();
358
359 // Ensure that the quadrupole order has a significant effect
360 size_t num_energies = std::min(energies2.size(), energies3.size());
361 for (size_t i = 0; i < num_energies; ++i) {
362 DOCTEST_CHECK(std::abs(energies3[i] - energies2[i]) * HARTREE_IN_GHZ > 1e-6);
363 }
364}
365
366} // namespace pairinteraction
Builder class for creating BasisAtom objects.
BasisAtomCreator< Scalar > & append_ket(const std::shared_ptr< const ket_t > &ket)
BasisAtomCreator< Scalar > & set_species(const std::string &value)
BasisAtomCreator< Scalar > & restrict_quantum_number_n(int min, int max)
BasisAtomCreator< Scalar > & restrict_quantum_number_nu(real_t min, real_t max)
static Database & get_global_instance()
Definition: Database.cpp:1068
Builder class for creating KetAtom objects.
KetAtomCreator & set_quantum_number_n(int value)
std::shared_ptr< const KetAtom > create(Database &database) const
KetAtomCreator & set_quantum_number_j(double value)
KetAtomCreator & set_quantum_number_l(double value)
KetAtomCreator & set_species(const std::string &value)
KetAtomCreator & set_quantum_number_m(double value)
Type & set_electric_field(const std::array< real_t, 3 > &field)
Definition: SystemAtom.cpp:26
const Eigen::SparseMatrix< scalar_t, Eigen::RowMajor > & get_matrix() const
Definition: System.cpp:110
System< Derived > & diagonalize(const DiagonalizerInterface< scalar_t > &diagonalizer, std::optional< real_t > min_eigenenergy={}, std::optional< real_t > max_eigenenergy={}, double rtol=1e-6)
Definition: System.cpp:184
Eigen::VectorX< real_t > get_eigenenergies() const
Definition: System.cpp:97
constexpr double UM_IN_ATOMIC_UNITS
constexpr double HARTREE_IN_GHZ
DOCTEST_TEST_CASE("create a basis for strontium 88")
constexpr double VOLT_PER_CM_IN_ATOMIC_UNITS