18#include <Eigen/Eigenvalues>
19#include <doctest/doctest.h>
20#include <fmt/ranges.h>
49 system.set_electric_field({0, 0, 0.0001});
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,
", ")));
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,
", ")));
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);
76 .restrict_quantum_number_l(0, 1)
85 diagonalize<SystemAtom<std::complex<double>>>({system1, system2}, diagonalizer);
88 auto matrix2 = system2.get_matrix();
89 for (
int i = 0; i < matrix1.rows(); ++i) {
90 for (
int j = 0; j < matrix1.cols(); ++j) {
92 DOCTEST_CHECK(std::abs(matrix1.coeff(i, j)) < 1e-10);
93 DOCTEST_CHECK(std::abs(matrix2.coeff(i, j)) < 1e-10);
100 doctest::skip(
true)) {
111 .set_species(
"Sr87_mqdt")
113 .restrict_quantum_number_l(0, 1)
116 std::vector<SystemAtom<std::complex<double>>> systems;
118 for (
int i = 0; i < n; ++i) {
121 systems.push_back(std::move(system));
124 DOCTEST_MESSAGE(
"Basis size: ", basis->get_number_of_states());
126 diagonalize<SystemAtom<std::complex<double>>>(systems, diagonalizer);
135 .restrict_quantum_number_l(0, 1)
143 Eigen::MatrixXcd matrix = system.
get_matrix();
144 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXcd> eigensolver;
145 eigensolver.compute(matrix);
146 auto eigenenergies_eigen = eigensolver.eigenvalues();
149 std::vector<std::unique_ptr<DiagonalizerInterface<std::complex<double>>>> diagonalizers;
150 std::vector<double> rtols;
152 DOCTEST_SUBCASE(
"Double precision") {
153 diagonalizers.push_back(std::make_unique<
DiagonalizerEigen<std::complex<double>>>());
159 diagonalizers.push_back(std::make_unique<
DiagonalizerFeast<std::complex<double>>>(300));
161 rtols = {1e-1, 1e-6, 1e-14};
162 eps = std::numeric_limits<double>::epsilon();
165 DOCTEST_SUBCASE(
"Single precision") {
166 diagonalizers.push_back(
169 diagonalizers.push_back(
171 diagonalizers.push_back(
175 diagonalizers.push_back(
178 rtols = {1e-1, 1e-6};
179 eps = std::numeric_limits<float>::epsilon();
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)");
189 for (
const auto &diagonalizer : diagonalizers) {
198 system.
diagonalize(*diagonalizer, std::numeric_limits<float>::lowest() / 2,
199 std::numeric_limits<float>::max() / 2, rtol_eigenenergies);
201 auto eigenvectors_pairinteraction = system.get_eigenbasis()->get_coefficients();
204 (eigenenergies_eigen - eigenenergies_pairinteraction).array().abs().maxCoeff() <
205 rtol_eigenenergies * matrix.norm());
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());
216 double min_energy = 0.153355;
217 double max_energy = 0.153360;
223 .restrict_quantum_number_n(58, 62)
224 .restrict_quantum_number_l(0, 1)
229 system.set_electric_field(
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]);
243 std::vector<std::unique_ptr<DiagonalizerInterface<double>>> diagonalizers;
254 for (
const auto &diagonalizer : diagonalizers) {
256 system.set_electric_field(
259 system.diagonalize(*diagonalizer, min_energy, max_energy, 1e-6);
260 auto eigenenergies_pairinteraction = system.get_eigenenergies();
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,
", ")));
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]) <
276#include <Eigen/Dense>
280 for (
size_t i = 0; i < 10; ++i) {
283 Eigen::MatrixXd matrix = Eigen::MatrixXd::Random(n, n);
284 matrix = (matrix + matrix.transpose()).eval();
285 Eigen::VectorXd eigenenergies(n);
289 LAPACKE_dsyev(LAPACK_COL_MAJOR,
'V',
'U', n, matrix.data(), n, eigenenergies.data());
290 DOCTEST_CHECK(info == 0);
296 double min_energy = -1;
297 double max_energy = -1;
303 .restrict_quantum_number_n(58, 62)
304 .restrict_quantum_number_l(0, 1)
307 std::vector<std::unique_ptr<DiagonalizerInterface<double>>> diagonalizers;
313 for (
const auto &diagonalizer : diagonalizers) {
315 system.set_electric_field(
318 system.diagonalize(*diagonalizer, min_energy, max_energy, 1e-6);
319 auto eigenenergies_pairinteraction = system.get_eigenenergies();
321 DOCTEST_CHECK(eigenenergies_pairinteraction.size() == 0);
336 double energy = ket->get_energy();
342 .restrict_quantum_number_n(58, 62)
343 .restrict_quantum_number_l(0, 3)
344 .restrict_quantum_number_m(0.5, 0.5)
348 system3.set_ion_interaction_order(3);
350 system3.diagonalize(diagonalizer, min_energy, max_energy, 1e-6);
351 auto energies3 = system3.get_eigenenergies();
354 system2.set_ion_interaction_order(2);
356 system2.diagonalize(diagonalizer, min_energy, max_energy, 1e-6);
357 auto energies2 = system2.get_eigenenergies();
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);
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()
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)
const Eigen::SparseMatrix< scalar_t, Eigen::RowMajor > & get_matrix() const
System< Derived > & diagonalize(const DiagonalizerInterface< scalar_t > &diagonalizer, std::optional< real_t > min_eigenenergy={}, std::optional< real_t > max_eigenenergy={}, double rtol=1e-6)
Eigen::VectorX< real_t > get_eigenenergies() const
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