openjij
Framework for the Ising model and QUBO.
Loading...
Searching...
No Matches
dense.hpp
Go to the documentation of this file.
1// Copyright 2023 Jij Inc.
2
3// Licensed under the Apache License, Version 2.0 (the "License");
4// you may not use this file except in compliance with the License.
5// You may obtain a copy of the License at
6
7// http://www.apache.org/licenses/LICENSE-2.0
8
9// Unless required by applicable law or agreed to in writing, software
10// distributed under the License is distributed on an "AS IS" BASIS,
11// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
12// See the License for the specific language governing permissions and
13// limitations under the License.
14
15#pragma once
16
17#include <algorithm>
18#include <cassert>
19#include <cstddef>
20#include <exception>
21#include <map>
22#include <type_traits>
23#include <vector>
24
25#include <Eigen/Dense>
26
28
31
32namespace openjij {
33namespace graph {
34
44template <typename FloatType> class Dense : public Graph {
45 static_assert(std::is_floating_point<FloatType>::value,
46 "FloatType must be floating-point type.");
47
48public:
64 Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
65
70
71private:
76
77public:
83 explicit Dense(std::size_t num_spins)
84 : Graph(num_spins), _J(Interactions::Zero(num_spins + 1, num_spins + 1)) {
85 _J(num_spins, num_spins) = 1;
86 }
87
93 Dense(const json &j) : Dense(static_cast<size_t>(j["num_variables"])) {
94 // define bqm with ising variables
96 // interactions
97 // for(auto&& elem : bqm.get_quadratic()){
98 // const auto& key = elem.first;
99 // const auto& val = elem.second;
100 // J(key.first, key.second) += val;
101 //}
102 // local field
103 // for(auto&& elem : bqm.get_linear()){
104 // const auto& key = elem.first;
105 // const auto& val = elem.second;
106 // h(key) += val;
107 //}
108 // the above insertion is simplified as
109 this->_J = bqm.interaction_matrix();
110 }
111
117 void set_interaction_matrix(const Interactions &interaction) {
118 if (interaction.rows() != interaction.cols()) {
119 std::runtime_error("interaction.rows() != interaction.cols()");
120 }
121
122 if ((size_t)interaction.rows() != get_num_spins() + 1) {
123 throw std::runtime_error("invalid matrix size.");
124 }
125
126 // check if diagonal elements are zero
127 for (size_t i = 0; i < (size_t)(interaction.rows() - 1); i++) {
128 if (interaction(i, i) != 0) {
129 throw std::runtime_error(
130 "The diagonal elements of interaction matrix must be zero.");
131 }
132 }
133
134 if (interaction(interaction.rows() - 1, interaction.rows() - 1) != 1) {
135 throw std::runtime_error(
136 "The right bottom element of interaction matrix must be unity.");
137 }
138
139 _J = interaction.template selfadjointView<Eigen::Upper>();
140 }
141
145 Dense(const Dense<FloatType> &) = default;
146
150 Dense(Dense<FloatType> &&) = default;
151
161 return this->energy(spins);
162 }
163
164 FloatType calc_energy(const Eigen::Matrix<FloatType, Eigen::Dynamic, 1,
165 Eigen::ColMajor> &spins) const {
166 return this->energy(spins);
167 }
168
176 FloatType energy(const Spins &spins) const {
177 if (spins.size() != this->get_num_spins()) {
178 throw std::out_of_range("Out of range in energy in Dense graph.");
179 }
180
181 using Vec = Eigen::Matrix<FloatType, Eigen::Dynamic, 1, Eigen::ColMajor>;
182 Vec s(get_num_spins() + 1);
183 for (size_t i = 0; i < spins.size(); i++) {
184 s(i) = spins[i];
185 }
186 s(get_num_spins()) = 1;
187
188 // the energy must be consistent with BinaryQuadraticModel.
189 return (s.transpose() *
190 (_J.template triangularView<Eigen::Upper>() * s))(0, 0) -
191 1;
192 }
193
194 FloatType energy(const Eigen::Matrix<FloatType, Eigen::Dynamic, 1,
195 Eigen::ColMajor> &spins) const {
197 for (size_t i = 0; i < temp_spins.size(); i++) {
198 temp_spins[i] = spins(i);
199 }
200 return energy(temp_spins);
201 }
202
212 assert(i < get_num_spins());
214
215 if (i != j)
216 return _J(std::min(i, j), std::max(i, j));
217 else
218 return _J(std::min(i, j), get_num_spins());
219 }
220
229 const FloatType &J(Index i, Index j) const {
230 assert(i < get_num_spins());
232
233 if (i != j)
234 return _J(std::min(i, j), std::max(i, j));
235 else
236 return _J(std::min(i, j), get_num_spins());
237 }
238
247 assert(i < get_num_spins());
248 return J(i, i);
249 }
250
258 const FloatType &h(Index i) const {
259 assert(i < get_num_spins());
260 return J(i, i);
261 }
262
281 return this->_J.template selfadjointView<Eigen::Upper>();
282 }
283};
284} // namespace graph
285} // namespace openjij
two-body all-to-all interactions The Hamiltonian is like
Definition dense.hpp:44
Dense(const json &j)
Dense constructor (from nlohmann::json)
Definition dense.hpp:93
FloatType value_type
float type
Definition dense.hpp:69
FloatType & J(Index i, Index j)
access J_{ij}
Definition dense.hpp:211
FloatType energy(const Spins &spins) const
calculate total energy
Definition dense.hpp:176
Interactions _J
interactions
Definition dense.hpp:75
FloatType energy(const Eigen::Matrix< FloatType, Eigen::Dynamic, 1, Eigen::ColMajor > &spins) const
Definition dense.hpp:194
void set_interaction_matrix(const Interactions &interaction)
set interaction matrix from Eigen Matrix.
Definition dense.hpp:117
FloatType calc_energy(const Spins &spins) const
calculate total energy
Definition dense.hpp:160
Eigen::Matrix< FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > Interactions
interaction type (Eigen) The stored matrix has the following triangular form:
Definition dense.hpp:64
const FloatType & h(Index i) const
access h_{i} (local field)
Definition dense.hpp:258
const Interactions get_interactions() const
get interactions (Eigen Matrix)
Definition dense.hpp:280
Dense(std::size_t num_spins)
Dense constructor.
Definition dense.hpp:83
Dense(Dense< FloatType > &&)=default
Dense move constructor.
Dense(const Dense< FloatType > &)=default
Dense copy constructor.
FloatType calc_energy(const Eigen::Matrix< FloatType, Eigen::Dynamic, 1, Eigen::ColMajor > &spins) const
Definition dense.hpp:164
const FloatType & J(Index i, Index j) const
access J_{ij}
Definition dense.hpp:229
FloatType & h(Index i)
access h_{i} (local field)
Definition dense.hpp:246
Abstract graph class.
Definition graph.hpp:37
std::size_t get_num_spins() const noexcept
get number of spins
Definition graph.hpp:89
auto json_parse(const json &obj, bool relabel=true)
parse json object from bqm.to_serializable
Definition parse.hpp:50
std::vector< Spin > Spins
Definition graph.hpp:27
nlohmann::json json
Definition parse.hpp:33
std::size_t Index
Definition graph.hpp:30
Definition algorithm.hpp:24
double FloatType
Note:
Definition compile_config.hpp:37