18#include <Eigen/Sparse>
31template <
typename GraphType,
int Options = Eigen::ColMajor>
40template <
typename FloatType,
int Options>
43 Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Options>;
52template <
typename FloatType,
int Options>
54 using type = Eigen::SparseMatrix<FloatType, Options>;
66template <
typename FloatType,
int Options = Eigen::ColMajor>
67inline static Eigen::Matrix<FloatType, Eigen::Dynamic, 1, Options>
69 Eigen::Matrix<FloatType, Eigen::Dynamic, 1, Options> ret_vec(
70 init_spin.size() + 1);
73 for (
size_t i = 0; i < init_spin.size(); i++) {
74 ret_vec(i) = init_spin[i];
78 ret_vec[init_spin.size()] = 1;
93template <
typename FloatType,
int Options = Eigen::ColMajor>
94inline static Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Options>
96 Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Options> ret_mat(
97 trotter_spins[0].size() + 1, trotter_spins.size());
100 for (
size_t j = 0; j < trotter_spins.size(); j++) {
101 for (
size_t i = 0; i < trotter_spins[j].size(); i++) {
102 ret_mat(i, j) = trotter_spins[j][i];
107 for (
size_t j = 0; j < trotter_spins.size(); j++) {
108 ret_mat(trotter_spins[0].size(), j) = 1;
124template <
int Options = Eigen::ColMajor,
typename FloatType>
125inline static Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Options>
128 Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Options> ret_mat(
135 ret_mat(i, j) = graph.
J(i, j);
136 ret_mat(j, i) = graph.
J(i, j);
162template <
int Options = Eigen::ColMajor,
typename FloatType>
163inline static Eigen::SparseMatrix<FloatType, Options>
166 Eigen::SparseMatrix<FloatType, Options> ret_mat(graph.
get_num_spins() + 1,
172 using T = std::vector<Eigen::Triplet<FloatType>>;
176 for (
size_t adj_ind : graph.
adj_nodes(ind)) {
177 if (ind != adj_ind) {
178 t_list.emplace_back(ind, adj_ind, graph.
J(ind, adj_ind));
188 ret_mat.setFromTriplets(t_list.begin(), t_list.end());
two-body all-to-all interactions The Hamiltonian is like
Definition dense.hpp:44
FloatType & J(Index i, Index j)
access J_{ij}
Definition dense.hpp:211
FloatType & h(Index i)
access h_{i} (local field)
Definition dense.hpp:246
std::size_t get_num_spins() const noexcept
get number of spins
Definition graph.hpp:89
Sparse graph: two-body intereactions with O(1) connectivity The Hamiltonian is like.
Definition sparse.hpp:40
FloatType & h(Index i)
access h_{i} (local field)
Definition sparse.hpp:300
FloatType & J(Index i, Index j)
access J_{ij}
Definition sparse.hpp:269
const Nodes & adj_nodes(Index ind) const
list of adjacent nodes
Definition sparse.hpp:201
std::vector< Spin > Spins
Definition graph.hpp:27
static Eigen::Matrix< FloatType, Eigen::Dynamic, Eigen::Dynamic, Options > gen_matrix_from_graph(const graph::Dense< FloatType > &graph)
generate Eigen Dense Matrix from Dense graph
Definition eigen.hpp:126
static Eigen::Matrix< FloatType, Eigen::Dynamic, 1, Options > gen_vector_from_std_vector(const graph::Spins &init_spin)
generate Eigen Vector from std::vector
Definition eigen.hpp:68
static Eigen::Matrix< FloatType, Eigen::Dynamic, Eigen::Dynamic, Options > gen_matrix_from_trotter_spins(const std::vector< graph::Spins > &trotter_spins)
generate Eigen Matrix from TrotterSpins
Definition eigen.hpp:95
Definition algorithm.hpp:24
Eigen::Matrix< FloatType, Eigen::Dynamic, Eigen::Dynamic, Options > type
Definition eigen.hpp:43
Eigen::SparseMatrix< FloatType, Options > type
Definition eigen.hpp:54
get Eigen Matrix type from Graph Type
Definition eigen.hpp:32