pairinteraction
A Rydberg Interaction Calculator
Basis.hpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2024 Pairinteraction Developers
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
4#pragma once
5
10
11#include <Eigen/Dense>
12#include <Eigen/SparseCore>
13#include <memory>
14#include <set>
15#include <unordered_map>
16#include <vector>
17
18namespace pairinteraction {
19enum class Parity : int;
20enum class TransformationType : unsigned char;
21enum class OperatorType;
22
23/**
24 * @class Basis
25 *
26 * @brief Base class for a basis
27 *
28 * This base class represents a basis. It comprises a list of ket states and a matrix of
29 * coefficients. The rows of the coefficient matrix correspond to indices of ket states and
30 * the columns to indices of basis vectors.
31 * Using CRPT, it is a base class for specific basis implementations. Its
32 * constructor is protected to indicate that derived classes should not allow direct instantiation.
33 * Instead, a factory class should be provided that is a friend of the derived class and can create
34 * instances of it.
35 *
36 * @tparam Derived Derived class.
37 */
38
39template <typename Derived>
40class Basis
41 : public TransformationBuilderInterface<typename traits::CrtpTraits<Derived>::scalar_t> {
42public:
47
48 Basis() = delete;
49 virtual ~Basis() = default;
50
51 bool has_quantum_number_f() const;
52 bool has_quantum_number_m() const;
53 bool has_parity() const;
54
55 const ketvec_t &get_kets() const;
56 size_t get_number_of_states() const;
57 size_t get_number_of_kets() const;
58 real_t get_quantum_number_f(size_t state_index) const;
59 real_t get_quantum_number_m(size_t state_index) const;
60 Parity get_parity(size_t state_index) const;
61 std::shared_ptr<const Derived> get_state(size_t state_index) const;
62 std::shared_ptr<const ket_t> get_ket(size_t ket_index) const;
63 std::shared_ptr<const ket_t> get_corresponding_ket(size_t state_index) const;
64 std::shared_ptr<const ket_t> get_corresponding_ket(std::shared_ptr<const Derived> state) const;
65 size_t get_corresponding_ket_index(size_t state_index) const;
66 size_t get_corresponding_ket_index(std::shared_ptr<const Derived> state) const;
67 std::shared_ptr<const Derived> get_corresponding_state(size_t ket_index) const;
68 std::shared_ptr<const Derived> get_corresponding_state(std::shared_ptr<const ket_t> ket) const;
69 size_t get_corresponding_state_index(size_t ket_index) const;
70 size_t get_corresponding_state_index(std::shared_ptr<const ket_t> ket) const;
71 std::shared_ptr<const Derived> get_canonical_state_from_ket(size_t ket_index) const;
72 std::shared_ptr<const Derived>
73 get_canonical_state_from_ket(std::shared_ptr<const ket_t> ket) const;
74 const Eigen::SparseMatrix<scalar_t, Eigen::RowMajor> &get_coefficients() const;
75 Eigen::SparseMatrix<scalar_t, Eigen::RowMajor> &get_coefficients();
76 Eigen::VectorX<scalar_t> get_amplitudes(std::shared_ptr<const ket_t> ket) const;
77 Eigen::SparseMatrix<scalar_t, Eigen::RowMajor>
78 get_amplitudes(std::shared_ptr<const Derived> other) const;
79 Eigen::VectorX<real_t> get_overlaps(std::shared_ptr<const ket_t> ket) const;
80 Eigen::SparseMatrix<real_t, Eigen::RowMajor>
81 get_overlaps(std::shared_ptr<const Derived> other) const;
82 virtual Eigen::VectorX<scalar_t> get_matrix_elements(std::shared_ptr<const ket_t> ket,
83 OperatorType type, int q = 0) const = 0;
84 Eigen::SparseMatrix<scalar_t, Eigen::RowMajor> virtual get_matrix_elements(
85 std::shared_ptr<const Derived> other, OperatorType type, int q = 0) const = 0;
86
87 class Iterator {
88 public:
89 Iterator(typename ketvec_t::const_iterator it);
90 bool operator!=(const Iterator &other) const;
91 std::shared_ptr<const ket_t> operator*() const;
93
94 private:
95 typename ketvec_t::const_iterator it;
96 };
97
98 Iterator begin() const;
99 Iterator end() const;
100
101 const Transformation<scalar_t> &get_transformation() const override;
102 Transformation<scalar_t> get_rotator(real_t alpha, real_t beta, real_t gamma) const override;
103 Sorting get_sorter(const std::vector<TransformationType> &labels) const override;
104 std::vector<IndicesOfBlock>
105 get_indices_of_blocks(const std::vector<TransformationType> &labels) const override;
106
107 void perform_sorter_checks(const std::vector<TransformationType> &labels) const;
108 void perform_blocks_checks(const std::set<TransformationType> &unique_labels) const;
109 void get_sorter_without_checks(const std::vector<TransformationType> &labels,
110 Sorting &transformation) const;
111 void get_indices_of_blocks_without_checks(const std::set<TransformationType> &unique_labels,
112 IndicesOfBlocksCreator &blocks) const;
113
114 std::shared_ptr<const Derived>
115 transformed(const Transformation<scalar_t> &transformation) const;
116 std::shared_ptr<const Derived> transformed(const Sorting &transformation) const;
117
118protected:
120 int get_ket_index_from_ket(std::shared_ptr<const ket_t> ket) const;
122
123private:
124 const Derived &derived() const;
125
126 struct hash {
127 std::size_t operator()(const std::shared_ptr<const ket_t> &k) const;
128 };
129
130 struct equal_to {
131 bool operator()(const std::shared_ptr<const ket_t> &lhs,
132 const std::shared_ptr<const ket_t> &rhs) const;
133 };
134
135 Transformation<scalar_t> coefficients;
136
137 std::unordered_map<std::shared_ptr<const ket_t>, size_t, hash, equal_to> ket_to_ket_index;
138 std::vector<size_t> ket_index_to_state_index;
139
140 std::vector<real_t> state_index_to_quantum_number_f;
141 std::vector<real_t> state_index_to_quantum_number_m;
142 std::vector<Parity> state_index_to_parity;
143 std::vector<size_t> state_index_to_ket_index;
144
145 bool _has_quantum_number_f{true};
146 bool _has_quantum_number_m{true};
147 bool _has_parity{true};
148};
149} // namespace pairinteraction
Iterator(typename ketvec_t::const_iterator it)
Definition: Basis.cpp:345
std::shared_ptr< const ket_t > operator*() const
Definition: Basis.cpp:353
bool operator!=(const Iterator &other) const
Definition: Basis.cpp:348
Base class for a basis.
Definition: Basis.hpp:41
Eigen::VectorX< real_t > get_overlaps(std::shared_ptr< const ket_t > ket) const
Definition: Basis.cpp:158
const Transformation< scalar_t > & get_transformation() const override
Definition: Basis.cpp:375
virtual Eigen::SparseMatrix< scalar_t, Eigen::RowMajor > get_matrix_elements(std::shared_ptr< const Derived > other, OperatorType type, int q=0) const =0
Iterator begin() const
Definition: Basis.cpp:335
void get_indices_of_blocks_without_checks(const std::set< TransformationType > &unique_labels, IndicesOfBlocksCreator &blocks) const
Definition: Basis.cpp:535
const ketvec_t & get_kets() const
Definition: Basis.cpp:114
bool has_parity() const
Definition: Basis.cpp:104
void get_sorter_without_checks(const std::vector< TransformationType > &labels, Sorting &transformation) const
Definition: Basis.cpp:454
int get_ket_index_from_ket(std::shared_ptr< const ket_t > ket) const
Definition: Basis.cpp:131
std::shared_ptr< const Derived > get_canonical_state_from_ket(size_t ket_index) const
Definition: Basis.cpp:296
typename traits::CrtpTraits< Derived >::ketvec_t ketvec_t
Definition: Basis.hpp:46
size_t get_corresponding_ket_index(size_t state_index) const
Definition: Basis.cpp:281
Iterator end() const
Definition: Basis.cpp:340
size_t get_number_of_states() const
Definition: Basis.cpp:364
virtual Eigen::VectorX< scalar_t > get_matrix_elements(std::shared_ptr< const ket_t > ket, OperatorType type, int q=0) const =0
Eigen::VectorX< scalar_t > get_amplitudes(std::shared_ptr< const ket_t > ket) const
Definition: Basis.cpp:140
typename traits::CrtpTraits< Derived >::scalar_t scalar_t
Definition: Basis.hpp:43
std::shared_ptr< const Derived > get_corresponding_state(size_t ket_index) const
Definition: Basis.cpp:244
typename traits::CrtpTraits< Derived >::ket_t ket_t
Definition: Basis.hpp:45
Transformation< scalar_t > get_rotator(real_t alpha, real_t beta, real_t gamma) const override
Definition: Basis.cpp:381
std::shared_ptr< const Derived > get_state(size_t state_index) const
Definition: Basis.cpp:212
void perform_sorter_checks(const std::vector< TransformationType > &labels) const
Definition: Basis.cpp:26
real_t get_quantum_number_m(size_t state_index) const
Definition: Basis.cpp:178
void perform_blocks_checks(const std::set< TransformationType > &unique_labels) const
Definition: Basis.cpp:36
std::shared_ptr< const ket_t > get_ket(size_t ket_index) const
Definition: Basis.cpp:239
Parity get_parity(size_t state_index) const
Definition: Basis.cpp:187
std::vector< IndicesOfBlock > get_indices_of_blocks(const std::vector< TransformationType > &labels) const override
Definition: Basis.cpp:440
real_t get_quantum_number_f(size_t state_index) const
Definition: Basis.cpp:169
bool has_quantum_number_m() const
Definition: Basis.cpp:99
std::shared_ptr< const ket_t > get_corresponding_ket(size_t state_index) const
Definition: Basis.cpp:197
typename traits::CrtpTraits< Derived >::real_t real_t
Definition: Basis.hpp:44
bool has_quantum_number_f() const
Definition: Basis.cpp:94
std::shared_ptr< const Derived > transformed(const Transformation< scalar_t > &transformation) const
Definition: Basis.cpp:614
virtual ~Basis()=default
Sorting get_sorter(const std::vector< TransformationType > &labels) const override
Definition: Basis.cpp:412
size_t get_corresponding_state_index(size_t ket_index) const
Definition: Basis.cpp:263
const Eigen::SparseMatrix< scalar_t, Eigen::RowMajor > & get_coefficients() const
Definition: Basis.cpp:120
size_t get_number_of_kets() const
Definition: Basis.cpp:369
Matrix< Type, Dynamic, 1 > VectorX
Helper struct to extract types from a derived basis type. Must be specialized for each derived basis ...
Definition: traits.hpp:24