pairinteraction
A Rydberg Interaction Calculator
System.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/SparseCore>
12#include <memory>
13#include <optional>
14#include <set>
15#include <vector>
16
17namespace pairinteraction {
18enum class TransformationType : unsigned char;
19
20template <typename Scalar>
21class DiagonalizerInterface;
22
23template <typename Derived>
24class System
25 : public TransformationBuilderInterface<typename traits::CrtpTraits<Derived>::scalar_t> {
26public:
32
33 System(std::shared_ptr<const basis_t> basis);
34 System(const System &other);
35 System(System &&other) noexcept;
36 System<Derived> &operator=(const System &other);
37 System<Derived> &operator=(System &&other) noexcept;
38 virtual ~System();
39
40 std::shared_ptr<const basis_t> get_basis() const;
41 std::shared_ptr<const basis_t> get_eigenbasis() const;
43
44 const Eigen::SparseMatrix<scalar_t, Eigen::RowMajor> &get_matrix() const;
45
46 const Transformation<scalar_t> &get_transformation() const override;
47 Transformation<scalar_t> get_rotator(real_t alpha, real_t beta, real_t gamma) const override;
48 Sorting get_sorter(const std::vector<TransformationType> &labels) const override;
49 std::vector<IndicesOfBlock>
50 get_indices_of_blocks(const std::vector<TransformationType> &labels) const override;
51
52 void get_indices_of_blocks_without_checks(const std::set<TransformationType> &unique_labels,
53 IndicesOfBlocksCreator &blocks) const;
54
56 System<Derived> &transform(const Sorting &transformation);
57
59 std::optional<real_t> min_eigenenergy = {},
60 std::optional<real_t> max_eigenenergy = {}, double rtol = 1e-6);
61 bool is_diagonal() const;
62
63protected:
64 mutable std::unique_ptr<operator_t> hamiltonian;
66 mutable bool hamiltonian_is_diagonal{false};
67 mutable std::vector<TransformationType> blockdiagonalizing_labels;
68
69 virtual void construct_hamiltonian() const = 0;
70
71private:
72 const Derived &derived() const;
73};
74} // namespace pairinteraction
const Eigen::SparseMatrix< scalar_t, Eigen::RowMajor > & get_matrix() const
Definition: System.cpp:110
std::vector< IndicesOfBlock > get_indices_of_blocks(const std::vector< TransformationType > &labels) const override
Definition: System.cpp:149
std::vector< TransformationType > blockdiagonalizing_labels
Definition: System.hpp:67
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
Transformation< scalar_t > get_rotator(real_t alpha, real_t beta, real_t gamma) const override
Definition: System.cpp:130
const Transformation< scalar_t > & get_transformation() const override
Definition: System.cpp:120
typename traits::CrtpTraits< Derived >::operator_t operator_t
Definition: System.hpp:31
void get_indices_of_blocks_without_checks(const std::set< TransformationType > &unique_labels, IndicesOfBlocksCreator &blocks) const
bool is_diagonal() const
Definition: System.cpp:321
std::shared_ptr< const basis_t > get_basis() const
Definition: System.cpp:76
System< Derived > & operator=(const System &other)
Definition: System.cpp:46
bool hamiltonian_requires_construction
Definition: System.hpp:65
virtual void construct_hamiltonian() const =0
Sorting get_sorter(const std::vector< TransformationType > &labels) const override
Definition: System.cpp:139
typename traits::CrtpTraits< Derived >::ketvec_t ketvec_t
Definition: System.hpp:29
typename traits::CrtpTraits< Derived >::real_t real_t
Definition: System.hpp:28
System(std::shared_ptr< const basis_t > basis)
Definition: System.cpp:28
Eigen::VectorX< real_t > get_eigenenergies() const
Definition: System.cpp:97
System< Derived > & transform(const Transformation< scalar_t > &transformation)
Definition: System.cpp:158
typename traits::CrtpTraits< Derived >::basis_t basis_t
Definition: System.hpp:30
std::shared_ptr< const basis_t > get_eigenbasis() const
Definition: System.cpp:85
std::unique_ptr< operator_t> hamiltonian
Definition: System.hpp:64
typename traits::CrtpTraits< Derived >::scalar_t scalar_t
Definition: System.hpp:27
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