pairinteraction
A Rydberg Interaction Calculator
GreenTensor.hpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2025 Pairinteraction Developers
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
4#pragma once
5
9
10#include <Eigen/Dense>
11#include <complex>
12#include <map>
13#include <unsupported/Eigen/Splines>
14#include <variant>
15#include <vector>
16
17namespace pairinteraction {
18template <typename Scalar>
20public:
22
24 using complex_t = std::complex<real_t>;
25
27 friend class GreenTensor;
28
29 public:
30 Scalar val() const;
31 int row() const noexcept;
32 int col() const noexcept;
33
34 private:
35 ConstantEntry(int row, int col, Scalar val);
36 int row_;
37 int col_;
38 Scalar val_;
39 };
40
42 friend class GreenTensor;
43
44 public:
45 Scalar val(double omega) const;
46 int row() const noexcept;
47 int col() const noexcept;
48
49 private:
50 OmegaDependentEntry(int row, int col, Eigen::Spline<real_t, 1> real_spline,
51 Eigen::Spline<real_t, 1> imag_spline);
52 int row_;
53 int col_;
54 Eigen::Spline<real_t, 1> real_spline;
55 Eigen::Spline<real_t, 1> imag_spline;
56 };
57
58 using Entry = std::variant<ConstantEntry, OmegaDependentEntry>;
59
60 void
61 create_entries_from_cartesian(int kappa1, int kappa2,
62 const Eigen::MatrixX<Scalar> &tensor_in_cartesian_coordinates);
64 int kappa1, int kappa2,
65 const std::vector<Eigen::MatrixX<Scalar>> &tensors_in_cartesian_coordinates,
66 const std::vector<double> &omegas);
67 const std::vector<Entry> &get_spherical_entries(int kappa1, int kappa2) const;
68
69private:
70 std::map<std::pair<int, int>, std::vector<Entry>> entries_map;
71};
72
73extern template class GreenTensor<double>;
74extern template class GreenTensor<std::complex<double>>;
75} // namespace pairinteraction
typename traits::NumTraits< Scalar >::real_t real_t
Definition: GreenTensor.hpp:23
std::variant< ConstantEntry, OmegaDependentEntry > Entry
Definition: GreenTensor.hpp:58
void create_entries_from_cartesian(int kappa1, int kappa2, const Eigen::MatrixX< Scalar > &tensor_in_cartesian_coordinates)
Definition: GreenTensor.cpp:65
const std::vector< Entry > & get_spherical_entries(int kappa1, int kappa2) const
std::complex< real_t > complex_t
Definition: GreenTensor.hpp:24
Matrix< Type, Dynamic, Dynamic > MatrixX
Helper struct to extract types from a numerical type.
Definition: traits.hpp:35