pairinteraction
A Rydberg Interaction Calculator
System.py.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2024 Pairinteraction Developers
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
4#include "./System.py.hpp"
5
13
14#include <nanobind/eigen/dense.h>
15#include <nanobind/eigen/sparse.h>
16#include <nanobind/nanobind.h>
17#include <nanobind/stl/array.h>
18#include <nanobind/stl/complex.h>
19#include <nanobind/stl/optional.h>
20#include <nanobind/stl/shared_ptr.h>
21#include <nanobind/stl/variant.h>
22#include <nanobind/stl/vector.h>
23
24namespace nb = nanobind;
25using namespace nb::literals;
26using namespace pairinteraction;
27
28template <typename T>
29static void declare_system(nb::module_ &m, std::string const &type_name) {
30 using S = System<T>;
31 using scalar_t = typename System<T>::scalar_t;
32
33 std::string pyclass_name = "System" + type_name;
34
35 nb::class_<System<T>, TransformationBuilderInterface<scalar_t>> pyclass(m,
36 pyclass_name.c_str());
37 pyclass.def("get_basis", &S::get_basis)
38 .def("get_eigenbasis", &S::get_eigenbasis)
39 .def("get_eigenenergies", &S::get_eigenenergies)
40 .def("get_matrix", &S::get_matrix)
41 .def("get_transformation", &S::get_transformation)
42 .def("get_rotator", &S::get_rotator)
43 .def("get_sorter", &S::get_sorter)
44 .def("get_indices_of_blocks", &S::get_indices_of_blocks)
45 .def("transform", nb::overload_cast<const Transformation<scalar_t> &>(&S::transform))
46 .def("transform", nb::overload_cast<const Sorting &>(&S::transform))
47 .def("diagonalize", &S::diagonalize, "diagonalizer"_a, "min_eigenenergy"_a = nb::none(),
48 "max_eigenenergy"_a = nb::none(), "rtol"_a = 1e-6)
49 .def("is_diagonal", &S::is_diagonal);
50}
51
52template <typename T>
53static void declare_system_atom(nb::module_ &m, std::string const &type_name) {
54 using S = SystemAtom<T>;
55 using basis_t = typename SystemAtom<T>::basis_t;
56
57 std::string pyclass_name = "SystemAtom" + type_name;
58
59 nb::class_<S, System<S>> pyclass(m, pyclass_name.c_str());
60 pyclass.def(nb::init<std::shared_ptr<const basis_t>>())
61 .def("set_electric_field", &S::set_electric_field)
62 .def("set_magnetic_field", &S::set_magnetic_field)
63 .def("set_diamagnetism_enabled", &S::set_diamagnetism_enabled)
64 .def("set_ion_distance_vector", &S::set_ion_distance_vector)
65 .def("set_ion_charge", &S::set_ion_charge)
66 .def("set_ion_interaction_order", &S::set_ion_interaction_order);
67}
68
69template <typename T>
70static void declare_system_pair(nb::module_ &m, std::string const &type_name) {
71 using S = SystemPair<T>;
72 using basis_t = typename SystemPair<T>::basis_t;
73
74 std::string pyclass_name = "SystemPair" + type_name;
75
76 nb::class_<S, System<S>> pyclass(m, pyclass_name.c_str());
77 pyclass.def(nb::init<std::shared_ptr<const basis_t>>())
78 .def("set_interaction_order", &S::set_interaction_order)
79 .def("set_distance_vector", &S::set_distance_vector)
80 .def("set_green_tensor", &S::set_green_tensor);
81}
82
83template <typename T>
84static void declare_green_tensor(nb::module_ &m, std::string const &type_name) {
85 using CE = typename GreenTensor<T>::ConstantEntry;
86 using OE = typename GreenTensor<T>::OmegaDependentEntry;
87 using GT = GreenTensor<T>;
88
89 std::string ce_name = "ConstantEntry" + type_name;
90 std::string oe_name = "OmegaDependentEntry" + type_name;
91 std::string gt_name = "GreenTensor" + type_name;
92
93 nb::class_<CE>(m, ce_name.c_str())
94 .def("row", &CE::row)
95 .def("col", &CE::col)
96 .def("val", &CE::val);
97
98 nb::class_<OE>(m, oe_name.c_str())
99 .def("row", &OE::row)
100 .def("col", &OE::col)
101 .def("val", &OE::val);
102
103 nb::class_<GT>(m, gt_name.c_str())
104 .def(nb::init<>())
105 .def("create_entries_from_cartesian",
106 nb::overload_cast<int, int, const Eigen::MatrixX<T> &>(
107 &GT::create_entries_from_cartesian))
108 .def("create_entries_from_cartesian",
109 nb::overload_cast<int, int, const std::vector<Eigen::MatrixX<T>> &,
110 const std::vector<double> &>(&GT::create_entries_from_cartesian))
111 .def("get_spherical_entries",
112 nb::overload_cast<int, int>(&GT::get_spherical_entries, nb::const_));
113}
114
115void bind_system(nb::module_ &m) {
116 declare_system<SystemAtom<double>>(m, "SystemAtomReal");
117 declare_system<SystemAtom<std::complex<double>>>(m, "SystemAtomComplex");
118 declare_system_atom<double>(m, "Real");
119 declare_system_atom<std::complex<double>>(m, "Complex");
120
121 declare_system<SystemPair<double>>(m, "SystemPairReal");
122 declare_system<SystemPair<std::complex<double>>>(m, "SystemPairComplex");
123 declare_system_pair<double>(m, "Real");
124 declare_system_pair<std::complex<double>>(m, "Complex");
125
126 declare_green_tensor<double>(m, "Real");
127 declare_green_tensor<std::complex<double>>(m, "Complex");
128}
void bind_system(nb::module_ &m)
Definition: System.py.cpp:115
typename traits::CrtpTraits< Type >::basis_t basis_t
Definition: SystemAtom.hpp:43
typename traits::CrtpTraits< Type >::basis_t basis_t
Definition: SystemPair.hpp:53
typename traits::CrtpTraits< Derived >::scalar_t scalar_t
Definition: System.hpp:27
Matrix< Type, Dynamic, Dynamic > MatrixX
void diagonalize(std::initializer_list< std::reference_wrapper< Derived > > systems, const DiagonalizerInterface< typename Derived::scalar_t > &diagonalizer, std::optional< typename Derived::real_t > min_eigenenergy={}, std::optional< typename Derived::real_t > max_eigenenergy={}, double rtol=1e-6)
Definition: diagonalize.cpp:17