openjij
Framework for the Ising model and QUBO.
Loading...
Searching...
No Matches
get_solution.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 <limits>
20
21#include "openjij/graph/all.hpp"
23
24#ifdef USE_CUDA
26#endif
27
28namespace openjij {
29namespace result {
30
39template <typename GraphType>
40const graph::Spins
42 // convert from Eigen::Vector to std::vector
43 graph::Spins ret_spins(system.num_spins);
44 for (std::size_t i = 0; i < system.num_spins; i++) {
45 ret_spins[i] = static_cast<graph::Spin>(system.spin(i) *
46 system.spin(system.num_spins));
47 }
48 return ret_spins;
49}
50
59template <typename GraphType>
60const graph::Spins
62 std::size_t minimum_trotter = 0;
63 // aliases
64 auto &spins = system.trotter_spins;
65 double energy = 0.0;
66 double min_energy = std::numeric_limits<double>::max();
67 // get number of trotter slices
68 std::size_t num_trotter_slices = system.trotter_spins.cols();
69 for (std::size_t t = 0; t < num_trotter_slices; t++) {
70 // calculate classical energy in each classical spin
71 energy = spins.col(t).transpose() * system.interaction * spins.col(t);
72 if (energy < min_energy) {
73 minimum_trotter = t;
74 min_energy = energy;
75 }
76 }
77
78 // convert from Eigen::Vector to std::vector
79 graph::Spins ret_spins(system.num_classical_spins);
80 for (std::size_t i = 0; i < system.num_classical_spins; i++) {
81 ret_spins[i] = static_cast<graph::Spin>(spins(i, minimum_trotter));
82 }
83 return ret_spins;
84}
85
86template <typename GraphType>
87const graph::Spins
89 return system.variables;
90}
91
92template <typename GraphType>
93const graph::Spins
95 return system.binaries;
96}
97
106template <typename GraphType>
109 auto spins = system.get_slice_at(0.0);
110
111 if (system.get_auxiliary_spin(0.0) < 0) {
112 /*if auxiliary spin is negative, flip all spins*/
113 for (auto &spin : spins) {
114 spin *= -1;
115 }
116 }
117
118 return spins;
119}
120
121#ifdef USE_CUDA
122
134template <typename FloatType, std::size_t rows_per_block,
135 std::size_t cols_per_block, std::size_t trotters_per_block>
136const graph::Spins
137get_solution(const system::ChimeraTransverseGPU<FloatType, rows_per_block,
138 cols_per_block,
139 trotters_per_block> &system) {
140
141 std::size_t localsize =
142 system.info.rows * system.info.cols * system.info.chimera_unitsize;
143
144 graph::Spins ret_spins(localsize);
145
146 size_t select_t = system.info.trotters / 2;
147 HANDLE_ERROR_CUDA(
148 cudaMemcpy(ret_spins.data(), system.spin.get() + (localsize * select_t),
149 localsize * sizeof(int), cudaMemcpyDeviceToHost));
150
151 return ret_spins;
152}
153
164template <typename FloatType, std::size_t rows_per_block,
165 std::size_t cols_per_block>
166const graph::Spins
167get_solution(const system::ChimeraClassicalGPU<FloatType, rows_per_block,
168 cols_per_block> &system) {
169 using Base = typename system::ChimeraClassicalGPU<FloatType, rows_per_block,
170 cols_per_block>::Base;
171
172 return get_solution(static_cast<const Base &>(system));
173}
174#endif
175
176} // namespace result
177} // namespace openjij
ClassicalIsingPolynomial class, which is a system to solve higher order unconstrained binary optimiza...
Definition classical_ising_polynomial.hpp:36
KLocalPolynomial class, which is a system to solve higher order unconstrained binary optimization (HU...
Definition k_local_polynomial.hpp:32
std::vector< Spin > Spins
Definition graph.hpp:27
int Spin
Definition graph.hpp:26
const graph::Spins get_solution(const system::ClassicalIsing< GraphType > &system)
get solution of classical ising system
Definition get_solution.hpp:41
Definition algorithm.hpp:24
double FloatType
Note:
Definition compile_config.hpp:37
ClassicalIsing structure (system for classical Ising model)
Definition classical_ising.hpp:38
Continuous Time Quantum Ising system.
Definition continuous_time_ising.hpp:34
TransverseIsing structure with discrete-time trotter spins.
Definition transverse_ising.hpp:39