openjij
Framework for the Ising model and QUBO.
Loading...
Searching...
No Matches
swendsen_wang.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 <cmath>
19#include <random>
20#include <unordered_map>
21
26
27namespace openjij {
28namespace updater {
29
35template <typename System> struct SwendsenWang;
36
42template <typename FloatType>
43struct SwendsenWang<system::ClassicalIsing<graph::Sparse<FloatType>>> {
44
46
47 template <typename RandomNumberEngine>
48 inline static void
49 update(ClIsing &system, RandomNumberEngine &random_number_engine,
50 const utility::ClassicalUpdaterParameter &parameter) {
51 auto urd = std::uniform_real_distribution<>(0, 1.0);
52
53 // num_spin = system size + additional spin
54 const size_t num_spin = system.spin.size();
55
56 // 1. update bonds
57 auto union_find_tree = utility::UnionFind(num_spin);
58 for (std::size_t node = 0; node < num_spin; ++node) {
59 for (typename ClIsing::SparseMatrixXx::InnerIterator it(
60 system.interaction, node);
61 it; ++it) {
62 // fetch adjacent node
63 std::size_t adj_node = it.index();
64 // fetch system.interaction(node, adj_node)
65 const FloatType &J = it.value();
66 if (node >= adj_node)
67 continue;
68 // check if bond can be connected
69 if (J * system.spin(node) * system.spin(adj_node) > 0)
70 continue;
71 const auto unite_rate =
72 std::max(static_cast<FloatType>(0.0),
73 static_cast<FloatType>(
74 1.0 - std::exp(-2.0 * parameter.beta * std::abs(J))));
75 if (urd(random_number_engine) < unite_rate)
76 union_find_tree.unite_sets(node, adj_node);
77 }
78 }
79
80 // 2. make clusters
81 const auto cluster_map = [num_spin, &union_find_tree]() {
82 auto cluster_map = std::unordered_multimap<utility::UnionFind::Node,
84 for (std::size_t node = 0; node < num_spin; ++node) {
85 cluster_map.insert({union_find_tree.find_set(node), node});
86 }
87 return cluster_map;
88 }();
89
90 // 3. update spin states in each cluster
91 for (auto &&c : union_find_tree.get_roots()) {
92 const auto range = cluster_map.equal_range(c);
93
94 // 3.1. decide spin state (flip with the probability 1/2)
95 const FloatType probability = 1.0 / 2.0;
96 if (urd(random_number_engine) < probability) {
97 // 3.2. update spin states
98 for (auto itr = range.first, last = range.second; itr != last; ++itr) {
99 const auto idx = itr->second;
100 system.spin(idx) *= -1;
101 }
102 }
103 }
104
105 return;
106 }
107};
108
114template <typename FloatType>
115struct SwendsenWang<system::ClassicalIsing<graph::CSRSparse<FloatType>>> {
116 //TODO: duplicate code of `SwendsenWang<system::ClassicalIsing<graph::Sparse<FloatType>>>`, need fix.
117
119
120 template <typename RandomNumberEngine>
121 inline static void
122 update(ClIsing &system, RandomNumberEngine &random_number_engine,
123 const utility::ClassicalUpdaterParameter &parameter) {
124 auto urd = std::uniform_real_distribution<>(0, 1.0);
125
126 // num_spin = system size + additional spin
127 const size_t num_spin = system.spin.size();
128
129 // 1. update bonds
130 auto union_find_tree = utility::UnionFind(num_spin);
131 for (std::size_t node = 0; node < num_spin; ++node) {
132 for (typename ClIsing::SparseMatrixXx::InnerIterator it(
133 system.interaction, node);
134 it; ++it) {
135 // fetch adjacent node
136 std::size_t adj_node = it.index();
137 // fetch system.interaction(node, adj_node)
138 const FloatType &J = it.value();
139 if (node >= adj_node)
140 continue;
141 // check if bond can be connected
142 if (J * system.spin(node) * system.spin(adj_node) > 0)
143 continue;
144 const auto unite_rate =
145 std::max(static_cast<FloatType>(0.0),
146 static_cast<FloatType>(
147 1.0 - std::exp(-2.0 * parameter.beta * std::abs(J))));
148 if (urd(random_number_engine) < unite_rate)
149 union_find_tree.unite_sets(node, adj_node);
150 }
151 }
152
153 // 2. make clusters
154 const auto cluster_map = [num_spin, &union_find_tree]() {
155 auto cluster_map = std::unordered_multimap<utility::UnionFind::Node,
157 for (std::size_t node = 0; node < num_spin; ++node) {
158 cluster_map.insert({union_find_tree.find_set(node), node});
159 }
160 return cluster_map;
161 }();
162
163 // 3. update spin states in each cluster
164 for (auto &&c : union_find_tree.get_roots()) {
165 const auto range = cluster_map.equal_range(c);
166
167 // 3.1. decide spin state (flip with the probability 1/2)
168 const FloatType probability = 1.0 / 2.0;
169 if (urd(random_number_engine) < probability) {
170 // 3.2. update spin states
171 for (auto itr = range.first, last = range.second; itr != last; ++itr) {
172 const auto idx = itr->second;
173 system.spin(idx) *= -1;
174 }
175 }
176 }
177
178 return;
179 }
180};
181} // namespace updater
182} // namespace openjij
Definition algorithm.hpp:24
double FloatType
Note:
Definition compile_config.hpp:37
ClassicalIsing structure for CSR Sparse graph (Eigen-based)
Definition classical_ising.hpp:190
VectorXx spin
spins (Eigen Vector)
Definition classical_ising.hpp:238
const SparseMatrixXx interaction
interaction (Eigen SparseMatrix)
Definition classical_ising.hpp:243
ClassicalIsing structure for Sparse graph (Eigen-based)
Definition classical_ising.hpp:118
VectorXx spin
spins (Eigen Vector)
Definition classical_ising.hpp:167
const SparseMatrixXx interaction
interaction (Eigen SparseMatrix)
Definition classical_ising.hpp:172
ClassicalIsing structure (system for classical Ising model)
Definition classical_ising.hpp:38
static void update(ClIsing &system, RandomNumberEngine &random_number_engine, const utility::ClassicalUpdaterParameter &parameter)
Definition swendsen_wang.hpp:122
static void update(ClIsing &system, RandomNumberEngine &random_number_engine, const utility::ClassicalUpdaterParameter &parameter)
Definition swendsen_wang.hpp:49
swendsen wang updater
Definition swendsen_wang.hpp:35
Definition union_find.hpp:24
std::size_t Node
Definition union_find.hpp:25
updater parameter for classical ising system
Definition schedule_list.hpp:37
double beta
inverse temperature
Definition schedule_list.hpp:46