cimod
C++ library for a binary (and polynomial) quadratic model.
Loading...
Searching...
No Matches
binary_quadratic_model_dict.hpp
Go to the documentation of this file.
1// Copyright 2020-2025 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 <cstdint>
19#include <functional>
20#include <iostream>
21#include <set>
22#include <string>
23#include <tuple>
24#include <typeinfo>
25#include <unordered_map>
26#include <unordered_set>
27#include <utility>
28#include <vector>
29
30#include <Eigen/Dense>
31
34#include "cimod/hash.hpp"
35#include "cimod/json.hpp"
36#include "cimod/utilities.hpp"
37#include "cimod/vartypes.hpp"
38
39namespace cimod {
40
41 struct Dict { };
42
47 template<typename IndexType, typename FloatType>
49 using DataType = Dict;
50
51 public:
52 using DenseMatrix = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
53 using SparseMatrix = Eigen::SparseMatrix<FloatType, Eigen::RowMajor>;
54
58 using Matrix = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>;
59
60 protected:
66
72
78
84
85 public:
95 const Linear<IndexType, FloatType> &linear,
96 const Quadratic<IndexType, FloatType> &quadratic,
97 const FloatType &offset,
98 const Vartype vartype ) :
99 m_offset( offset ),
100 m_vartype( vartype ) {
101 add_variables_from( linear );
102 add_interactions_from( quadratic );
103 }
104
113 const Linear<IndexType, FloatType> &linear,
114 const Quadratic<IndexType, FloatType> &quadratic,
115 const Vartype vartype ) :
116 BinaryQuadraticModel( linear, quadratic, 0.0, vartype ) {
117 }
118
130 const Eigen::Ref<const DenseMatrix> &,
131 const std::vector<IndexType> &,
132 const FloatType &,
133 const Vartype,
134 bool ) {
135 throw std::runtime_error( "Initialization from matrix is not implemented on dict-type BQM" );
136 }
137
147 BinaryQuadraticModel( const Eigen::Ref<const DenseMatrix> &, const std::vector<IndexType> &, const Vartype, bool ) {
148 throw std::runtime_error( "Initialization from matrix is not implemented on dict-type BQM" );
149 }
150
161 BinaryQuadraticModel( const SparseMatrix &, const std::vector<IndexType> &, const FloatType &, const Vartype ) {
162 throw std::runtime_error( "Initialization from matrix is not implemented on dict-type BQM" );
163 }
164
174 BinaryQuadraticModel( const SparseMatrix &, const std::vector<IndexType> &, const Vartype ) {
175 throw std::runtime_error( "Initialization from matrix is not implemented on dict-type BQM" );
176 }
177
184
190 std::vector<IndexType> _generate_indices() const {
191 std::unordered_set<IndexType> index_set;
192 for ( auto elem : m_linear ) {
193 index_set.insert( elem.first );
194 }
195
196 for ( auto elem : m_quadratic ) {
197 index_set.insert( elem.first.first );
198 index_set.insert( elem.first.second );
199 }
200
201 auto ret = std::vector<IndexType>( index_set.begin(), index_set.end() );
202 std::sort( ret.begin(), ret.end() );
203 return ret;
204 }
205
212 size_t length() const {
213 return get_num_variables();
214 }
215
222 size_t get_num_variables() const {
223 return m_linear.size();
224 }
225
232 bool contains( const IndexType &v ) const {
233 if ( m_linear.count( v ) != 0 ) {
234 return true;
235 } else {
236 return false;
237 }
238 }
239
248 return this->m_linear.at( label_i );
249 }
250
257 return this->m_linear;
258 }
259
266 return this->m_quadratic.at( std::make_pair( std::min( label_i, label_j ), std::max( label_i, label_j ) ) );
267 }
268
275 return this->m_quadratic;
276 }
277
283 const FloatType &get_offset() const {
284 return this->m_offset;
285 }
286
292 const Vartype &get_vartype() const {
293 return this->m_vartype;
294 }
295
301 const std::vector<IndexType> get_variables() const {
302 std::vector<IndexType> variables;
303 variables.reserve( m_linear.size() );
304 for ( auto &&elem : m_linear ) {
305 variables.push_back( elem.first );
306 }
307
308 std::sort( variables.begin(), variables.end() );
309
310 return variables;
311 }
312
321
322 /* Update methods */
323
330 void add_variable( const IndexType &v, const FloatType &bias ) {
331 FloatType b = bias;
332
333 Vartype vartype = this->get_vartype();
334
335 // handle the case that a different vartype is provided
336 if ( ( vartype != Vartype::NONE ) && ( vartype != m_vartype ) ) {
337 if ( ( m_vartype == Vartype::SPIN ) && ( vartype == Vartype::BINARY ) ) {
338 b /= 2;
339 m_offset += b;
340 } else if ( ( m_vartype == Vartype::BINARY ) && ( vartype == Vartype::SPIN ) ) {
341 m_offset -= b;
342 b *= 2;
343 } else {
344 throw std::runtime_error( "Unknown vartype" );
345 }
346 }
347
348 // Insert or assign the bias
349 FloatType value = 0;
350 if ( m_linear.count( v ) != 0 ) {
351 value = m_linear[ v ];
352 }
353 insert_or_assign( m_linear, v, value + b );
354 }
355
362 for ( auto &it : linear ) {
363 add_variable( it.first, it.second );
364 }
365 }
366
374 void add_interaction( const IndexType &arg_u, const IndexType &arg_v, const FloatType &bias ) {
375 IndexType u = std::min( arg_u, arg_v );
376 IndexType v = std::max( arg_u, arg_v );
377
378 Vartype vartype = this->get_vartype();
379
380 if ( u == v ) {
381 throw std::runtime_error( "No self-loops allowed" );
382 } else {
383 // Check vartype when m_linear.empty() is true
384 if ( m_linear.empty() && m_vartype == Vartype::NONE ) {
385 if ( vartype != Vartype::NONE ) {
386 m_vartype = vartype;
387 } else {
388 throw std::runtime_error(
389 "Binary quadratic model is empty. Please set vartype to Vartype::SPIN or Vartype::BINARY" );
390 }
391 }
392
393 FloatType b = bias;
394 if ( ( vartype != Vartype::NONE ) && ( vartype != m_vartype ) ) {
395 if ( ( m_vartype == Vartype::SPIN ) && ( vartype == Vartype::BINARY ) ) {
396 // convert from binary to spin
397 b /= 4;
398 add_offset( b );
399 add_variable( u, b );
400 add_variable( v, b );
401 } else if ( ( m_vartype == Vartype::BINARY ) && ( vartype == Vartype::SPIN ) ) {
402 // convert from spin to binary
403 add_offset( b );
404 add_variable( u, -2 * b );
405 add_variable( v, -2 * b );
406 b *= 4;
407 } else {
408 throw std::runtime_error( "Unknown vartype" );
409 }
410 } else {
411 if ( m_linear.count( u ) == 0 ) {
412 add_variable( u, 0 );
413 }
414 if ( m_linear.count( v ) == 0 ) {
415 add_variable( v, 0 );
416 }
417 }
418
419 FloatType value = 0;
420 std::pair<IndexType, IndexType> p1 = std::make_pair( u, v );
421 if ( m_quadratic.count( p1 ) != 0 ) {
422 value = m_quadratic[ p1 ];
423 }
424 insert_or_assign( m_quadratic, p1, value + b );
425 }
426 }
427
434 for ( auto &it : quadratic ) {
435 add_interaction( it.first.first, it.first.second, it.second );
436 }
437 }
438
444 void remove_variable( const IndexType &v ) {
445 std::vector<std::pair<IndexType, IndexType>> interactions;
446 for ( auto &it : m_quadratic ) {
447 if ( it.first.first == v || it.first.second == v ) {
448 interactions.push_back( it.first );
449 }
450 }
452 m_linear.erase( v );
453 }
454
460 void remove_variables_from( const std::vector<IndexType> &variables ) {
461 for ( auto &it : variables ) {
463 }
464 }
465
473 IndexType u = std::min( arg_u, arg_v );
474 IndexType v = std::max( arg_u, arg_v );
475 auto p = std::make_pair( u, v );
476 if ( m_quadratic.count( p ) != 0 ) {
477 m_quadratic.erase( p );
478 }
479
480 bool u_flag = true;
481 for ( auto &it : m_quadratic ) {
482 if ( ( it.first.first == u ) || ( it.first.second == u ) ) {
483 u_flag = false;
484 break;
485 }
486 }
487
488 if ( u_flag && ( m_linear[ u ] == 0 ) ) {
490 }
491
492 bool v_flag = true;
493 for ( auto &it : m_quadratic ) {
494 if ( ( it.first.first == v ) || ( it.first.second == v ) ) {
495 v_flag = false;
496 break;
497 }
498 }
499
500 if ( v_flag && ( m_linear[ v ] == 0 ) ) {
502 }
503 }
504
510 void remove_interactions_from( const std::vector<std::pair<IndexType, IndexType>> &interactions ) {
511 for ( auto &it : interactions ) {
512 remove_interaction( it.first, it.second );
513 }
514 }
515
521 void add_offset( const FloatType &offset ) {
522 m_offset += offset;
523 }
524
530 }
531
540 void scale(
541 const FloatType &scalar,
542 const std::vector<IndexType> &ignored_variables = {},
543 const std::vector<std::pair<IndexType, IndexType>> &ignored_interactions = {},
544 const bool ignored_offset = false ) {
545 // scaling linear
546 for ( auto &it : m_linear ) {
547 if ( std::find( ignored_variables.begin(), ignored_variables.end(), it.first ) != ignored_variables.end()
548 || ignored_variables.empty() ) {
549 it.second *= scalar;
550 }
551 }
552
553 // scaling quadratic
554 for ( auto &it : m_quadratic ) {
555 if ( std::find( ignored_interactions.begin(), ignored_interactions.end(), it.first ) != ignored_interactions.end()
556 || ignored_variables.empty() ) {
557 it.second *= scalar;
558 }
559 }
560
561 // scaling offset
562 if ( !ignored_offset ) {
563 m_offset *= scalar;
564 }
565 }
566
580 const std::pair<FloatType, FloatType> &bias_range = { 1.0, 1.0 },
581 const bool use_quadratic_range = false,
582 const std::pair<FloatType, FloatType> &quadratic_range = { 1.0, 1.0 },
583 const std::vector<IndexType> &ignored_variables = {},
584 const std::vector<std::pair<IndexType, IndexType>> &ignored_interactions = {},
585 const bool ignored_offset = false ) {
586 if ( m_linear.empty() ) {
587 return;
588 }
589 // parse range
590 std::pair<FloatType, FloatType> l_range = bias_range;
591 std::pair<FloatType, FloatType> q_range;
592 if ( !use_quadratic_range ) {
594 } else {
596 }
597
598 // calculate scaling value
599 auto comp = []( const auto &a, const auto &b ) { return a.second < b.second; };
600 auto it_lin_min = std::min_element( m_linear.begin(), m_linear.end(), comp );
601 auto it_lin_max = std::max_element( m_linear.begin(), m_linear.end(), comp );
602 auto it_quad_min = std::min_element( m_quadratic.begin(), m_quadratic.end(), comp );
603 auto it_quad_max = std::max_element( m_quadratic.begin(), m_quadratic.end(), comp );
604
605 std::vector<FloatType> v_scale
606 = { it_lin_min->second / l_range.first,
607 it_lin_max->second / l_range.second,
608 it_quad_min->second / q_range.first,
609 it_quad_max->second / q_range.second };
610 FloatType inv_scale = *std::max_element( v_scale.begin(), v_scale.end() );
611
612 // scaling
613 if ( inv_scale != 0.0 ) {
615 }
616 }
617
624 void fix_variable( const IndexType &v, const int32_t &value ) {
625 std::vector<std::pair<IndexType, IndexType>> interactions;
626 for ( auto &it : m_quadratic ) {
627 if ( it.first.first == v ) {
628 add_variable( it.first.second, value * it.second );
629 interactions.push_back( it.first );
630 } else if ( it.first.second == v ) {
631 add_variable( it.first.first, value * it.second );
632 interactions.push_back( it.first );
633 }
634 }
636 add_offset( m_linear[ v ] * value );
638 }
639
645 void fix_variables( const std::vector<std::pair<IndexType, int32_t>> &fixed ) {
646 for ( auto &it : fixed ) {
647 fix_variable( it.first, it.second );
648 }
649 }
650
656 void flip_variable( const IndexType &v ) {
657 // check variable
658 if ( m_linear.count( v ) == 0 ) {
659 throw std::runtime_error( "not a variable in the binary quadratic model." );
660 }
661
662 if ( m_vartype == Vartype::SPIN ) {
663 m_linear[ v ] *= -1.0;
664 for ( auto &it : m_quadratic ) {
665 if ( it.first.first == v || it.first.second == v ) {
666 it.second *= -1.0;
667 }
668 }
669 } else if ( m_vartype == Vartype::BINARY ) {
670 add_offset( m_linear[ v ] );
671 m_linear[ v ] *= -1.0;
672
673 for ( auto &it : m_quadratic ) {
674 if ( it.first.first == v ) {
675 m_linear[ it.first.second ] += it.second;
676 it.second *= -1.0;
677 } else if ( it.first.second == v ) {
678 m_linear[ it.first.first ] += it.second;
679 it.second *= -1.0;
680 }
681 }
682 }
683 }
684
691 void contract_variables( const IndexType &u, const IndexType &v ) {
692 // check variable
693 if ( m_linear.count( v ) == 0 ) {
694 throw std::runtime_error( "not a variable in the binary quadratic model." );
695 }
696 if ( m_linear.count( u ) == 0 ) {
697 throw std::runtime_error( "not a variable in the binary quadratic model." );
698 }
699
700 auto p1 = std::make_pair( u, v );
701 auto p2 = std::make_pair( v, u );
702 if ( m_quadratic.count( p1 ) != 0 ) {
703 if ( m_vartype == Vartype::BINARY ) {
704 add_variable( u, m_quadratic[ p1 ] );
705 } else if ( m_vartype == Vartype::SPIN ) {
706 add_offset( m_quadratic[ p1 ] );
707 }
709 }
710 if ( m_quadratic.count( p2 ) != 0 ) {
711 if ( m_vartype == Vartype::BINARY ) {
712 add_variable( u, m_quadratic[ p2 ] );
713 } else if ( m_vartype == Vartype::SPIN ) {
714 add_offset( m_quadratic[ p2 ] );
715 }
717 }
718
719 std::vector<std::pair<IndexType, IndexType>> interactions;
720 for ( auto &it : m_quadratic ) {
721 if ( it.first.first == v ) {
722 add_interaction( u, it.first.second, it.second );
723 interactions.push_back( it.first );
724 } else if ( it.first.second == v ) {
725 add_interaction( it.first.first, u, it.second );
726 interactions.push_back( it.first );
727 }
728 }
730
731 add_variable( u, m_linear[ v ] );
733 }
734
735 /* Transformations */
736
743 void change_vartype( const Vartype &vartype ) {
746 FloatType offset = 0.0;
747
748 if ( m_vartype == Vartype::BINARY && vartype == Vartype::SPIN ) // binary -> spin
749 {
750 std::tie( linear, quadratic, offset ) = binary_to_spin( m_linear, m_quadratic, m_offset );
751 } else if ( m_vartype == Vartype::SPIN && vartype == Vartype::BINARY ) // spin -> binary
752 {
753 std::tie( linear, quadratic, offset ) = spin_to_binary( m_linear, m_quadratic, m_offset );
754 } else {
755 std::tie( linear, quadratic, offset ) = std::tie( m_linear, m_quadratic, m_offset );
756 }
757
758 // inplace
759 m_linear = linear;
760 m_quadratic = quadratic;
761 m_offset = offset;
762 m_vartype = vartype;
763 }
764
775 FloatType offset = 0.0;
776
777 if ( m_vartype == Vartype::BINARY && vartype == Vartype::SPIN ) // binary -> spin
778 {
779 std::tie( linear, quadratic, offset ) = binary_to_spin( m_linear, m_quadratic, m_offset );
780 } else if ( m_vartype == Vartype::SPIN && vartype == Vartype::BINARY ) // spin -> binary
781 {
782 std::tie( linear, quadratic, offset ) = spin_to_binary( m_linear, m_quadratic, m_offset );
783 } else {
784 std::tie( linear, quadratic, offset ) = std::tie( m_linear, m_quadratic, m_offset );
785 }
786
787 BinaryQuadraticModel bqm( linear, quadratic, offset, vartype );
788
789 if ( inplace == true ) {
790 // inplace
791 m_linear = bqm.get_linear();
792 m_quadratic = bqm.get_quadratic();
793 m_offset = bqm.get_offset();
794 m_vartype = bqm.get_vartype();
795 }
796
797 return bqm;
798 }
799
800 /* Static methods */
811 static std::tuple<Linear<IndexType, FloatType>, Quadratic<IndexType, FloatType>, FloatType> spin_to_binary(
812 const Linear<IndexType, FloatType> &linear,
813 const Quadratic<IndexType, FloatType> &quadratic,
814 const FloatType &offset ) {
817 FloatType new_offset = offset;
820
821 for ( auto &it : linear ) {
822 insert_or_assign( new_linear, it.first, static_cast<FloatType>( 2.0 * it.second ) );
823 linear_offset += it.second;
824 }
825
826 for ( auto &it : quadratic ) {
827 insert_or_assign( new_quadratic, it.first, static_cast<FloatType>( 4.0 * it.second ) );
828 new_linear[ it.first.first ] -= 2.0 * it.second;
829 new_linear[ it.first.second ] -= 2.0 * it.second;
830 quadratic_offset += it.second;
831 }
832
834
835 return std::make_tuple( new_linear, new_quadratic, new_offset );
836 }
837
848 static std::tuple<Linear<IndexType, FloatType>, Quadratic<IndexType, FloatType>, FloatType> binary_to_spin(
849 const Linear<IndexType, FloatType> &linear,
850 const Quadratic<IndexType, FloatType> &quadratic,
851 const FloatType &offset ) {
854 FloatType new_offset = offset;
857
858 for ( auto &it : linear ) {
859 insert_or_assign( h, it.first, static_cast<FloatType>( 0.5 * it.second ) );
860 linear_offset += it.second;
861 }
862
863 for ( auto &it : quadratic ) {
864 insert_or_assign( J, it.first, static_cast<FloatType>( 0.25 * it.second ) );
865 h[ it.first.first ] += 0.25 * it.second;
866 h[ it.first.second ] += 0.25 * it.second;
867 quadratic_offset += it.second;
868 }
869
870 new_offset += 0.5 * linear_offset + 0.25 * quadratic_offset;
871
872 return std::make_tuple( h, J, new_offset );
873 }
874
875 /* Methods */
876
885 for ( auto &&it : m_linear ) {
886 if ( check_vartype( sample.at( it.first ), m_vartype ) ) {
887 en += static_cast<FloatType>( sample.at( it.first ) ) * it.second;
888 }
889 }
890 for ( auto &it : m_quadratic ) {
891 if ( check_vartype( sample.at( it.first.first ), m_vartype )
892 && check_vartype( sample.at( it.first.second ), m_vartype ) ) {
893 en += static_cast<FloatType>( sample.at( it.first.first ) )
894 * static_cast<FloatType>( sample.at( it.first.second ) ) * it.second;
895 }
896 }
897 return en;
898 }
899
906 std::vector<FloatType> energies( const std::vector<Sample<IndexType>> &samples_like ) const {
907 std::vector<FloatType> en_vec;
908 for ( auto &it : samples_like ) {
909 en_vec.push_back( energy( it ) );
910 }
911 return en_vec;
912 }
913
914 /* Conversions */
920 std::tuple<Quadratic<IndexType, FloatType>, FloatType> to_qubo() {
921 // change vartype to binary
923
924 Linear<IndexType, FloatType> linear = bqm.get_linear();
925 Quadratic<IndexType, FloatType> Q = bqm.get_quadratic();
926 FloatType offset = bqm.get_offset();
927 for ( auto &it : linear ) {
928 insert_or_assign( Q, std::make_pair( it.first, it.first ), it.second );
929 }
930 return std::make_tuple( Q, offset );
931 }
932
944
945 for ( auto &&elem : Q ) {
946 const auto &key = elem.first;
947 const auto &value = elem.second;
948 if ( key.first == key.second ) {
949 linear[ key.first ] = value;
950 } else {
951 quadratic[ std::make_pair( key.first, key.second ) ] = value;
952 }
953 }
954
956 }
957
963 std::tuple<Linear<IndexType, FloatType>, Quadratic<IndexType, FloatType>, FloatType> to_ising() {
964 // change vartype to spin
966
967 Linear<IndexType, FloatType> linear = bqm.get_linear();
968 Quadratic<IndexType, FloatType> quadratic = bqm.get_quadratic();
969 FloatType offset = bqm.get_offset();
970 return std::make_tuple( linear, quadratic, offset );
971 }
972
983 const Linear<IndexType, FloatType> &linear,
984 const Quadratic<IndexType, FloatType> &quadratic,
985 FloatType offset = 0.0 ) {
986 return BinaryQuadraticModel<IndexType, FloatType, DataType>( linear, quadratic, offset, Vartype::SPIN );
987 }
988
1026 Matrix interaction_matrix( const std::vector<IndexType> &indices ) const {
1027 // generate matrix
1028 size_t system_size = indices.size();
1029 Matrix _interaction_matrix = Matrix::Zero( system_size, system_size );
1030 const Linear<IndexType, FloatType> &linear = m_linear;
1031 const Quadratic<IndexType, FloatType> &quadratic = m_quadratic;
1032
1033 for ( size_t i = 0; i < indices.size(); i++ ) {
1034 const IndexType &i_index = indices[ i ];
1035 _interaction_matrix( i, i ) = ( linear.find( i_index ) != linear.end() ) ? linear.at( i_index ) : 0;
1036 for ( size_t j = i + 1; j < indices.size(); j++ ) {
1037 const IndexType &j_index = indices[ j ];
1038 FloatType jval = 0.0;
1039
1040 if ( quadratic.find( std::make_pair( i_index, j_index ) ) != quadratic.end() ) {
1041 jval += quadratic.at( std::make_pair( i_index, j_index ) );
1042 }
1043 if ( quadratic.find( std::make_pair( j_index, i_index ) ) != quadratic.end() ) {
1044 jval += quadratic.at( std::make_pair( j_index, i_index ) );
1045 }
1046
1047 _interaction_matrix( i, j ) = jval;
1048 _interaction_matrix( j, i ) = jval;
1049 }
1050 }
1051
1052 return _interaction_matrix;
1053 }
1054
1070 // generate matrix
1071 auto indices = this->get_variables();
1072 size_t system_size = this->get_num_variables();
1073 Matrix _interaction_matrix = Matrix::Zero( system_size + 1, system_size + 1 );
1074 const Linear<IndexType, FloatType> &linear = m_linear;
1075 const Quadratic<IndexType, FloatType> &quadratic = m_quadratic;
1076
1077 for ( size_t i = 0; i < indices.size(); i++ ) {
1078 const IndexType &i_index = indices[ i ];
1079 _interaction_matrix( i, system_size ) = ( linear.find( i_index ) != linear.end() ) ? linear.at( i_index ) : 0;
1080 for ( size_t j = i + 1; j < indices.size(); j++ ) {
1081 const IndexType &j_index = indices[ j ];
1082 FloatType jval = 0.0;
1083
1084 if ( quadratic.find( std::make_pair( i_index, j_index ) ) != quadratic.end() ) {
1085 jval += quadratic.at( std::make_pair( i_index, j_index ) );
1086 }
1087 if ( quadratic.find( std::make_pair( j_index, i_index ) ) != quadratic.end() ) {
1088 jval += quadratic.at( std::make_pair( j_index, i_index ) );
1089 }
1090
1091 _interaction_matrix( i, j ) = jval;
1092 }
1093 }
1094
1095 return _interaction_matrix;
1096 }
1097
1098 using json = nlohmann::json;
1099
1107 std::string schema_version = "3.0.0";
1108 // all variables are contained in the keys of m_linear.
1109 // All we need is to traverse the keys of m_linear.
1110 /*
1111 * output sample
1112 * >>> bqm = dimod.BinaryQuadraticModel({'c': -1, 'd': 1}, {('a', 'd'): 2, ('b', 'e'): 5, ('a', 'c'): 3}, 0.0,
1113 * dimod.BINARY)
1114 *
1115 * >>> bqm.to_serializable()
1116 * {'type': 'BinaryQuadraticModel', 'version': {'bqm_schema': '3.0.0'}, 'use_bytes': False, 'index_type': 'uint16',
1117 * 'bias_type': 'float32', 'num_variables': 5, 'num_interactions': 3, 'variable_labels': ['a', 'b', 'c', 'd', 'e'],
1118 * 'variable_type': 'BINARY', 'offset': 0.0, 'info': {}, 'linear_biases': [0.0, 0.0, -1.0, 1.0, 0.0],
1119 * 'quadratic_biases': [3.0, 2.0, 5.0], 'quadratic_head': [0, 0, 1], 'quadratic_tail': [2, 3, 4]}
1120 */
1121
1122 // set variables (sorted)
1123 std::vector<IndexType> variables;
1124 variables.reserve( m_linear.size() );
1125 for ( auto &&elem : m_linear ) {
1126 variables.push_back( elem.first );
1127 }
1128
1129 std::sort( variables.begin(), variables.end() );
1130
1131 size_t num_variables = variables.size();
1132
1133 // set sorted linear biases
1134 std::vector<FloatType> l_bias;
1135 for ( auto &&key : variables ) {
1136 l_bias.push_back( m_linear.at( key ) );
1137 }
1138
1139 // set quadratic head, tail and biases
1140 std::vector<size_t> q_head, q_tail;
1141 std::vector<FloatType> q_bias;
1142 for ( auto &&elem : m_quadratic ) {
1143 auto it_head = std::find( variables.begin(), variables.end(), elem.first.first );
1144 auto it_tail = std::find( variables.begin(), variables.end(), elem.first.second );
1145 size_t idx_head = std::distance( variables.begin(), it_head );
1146 size_t idx_tail = std::distance( variables.begin(), it_tail );
1147 q_head.push_back( idx_head );
1148 q_tail.push_back( idx_tail );
1149 q_bias.push_back( elem.second );
1150 }
1151
1152 size_t num_interactions = m_quadratic.size();
1153
1154 // set index_dtype
1155 std::string index_dtype = num_variables <= 65536UL ? "uint16" : "uint32";
1156
1157 // set bias_type
1158 std::string bias_type;
1159 if ( typeid( m_offset ) == typeid( float ) ) {
1160 bias_type = "float32";
1161 } else if ( typeid( m_offset ) == typeid( double ) ) {
1162 bias_type = "float64";
1163 } else {
1164 throw std::runtime_error( "FloatType must be float or double." );
1165 }
1166
1167 // set variable type
1168 std::string variable_type;
1169 if ( m_vartype == Vartype::SPIN ) {
1170 variable_type = "SPIN";
1171 } else if ( m_vartype == Vartype::BINARY ) {
1172 variable_type = "BINARY";
1173 } else {
1174 throw std::runtime_error( "Variable type must be SPIN or BINARY." );
1175 }
1176
1177 json output;
1178 output[ "type" ] = "BinaryQuadraticModel";
1179 output[ "version" ] = { { "bqm_schema", "3.0.0" } };
1180 output[ "variable_labels" ] = variables;
1181 output[ "use_bytes" ] = false;
1182 output[ "index_type" ] = index_dtype;
1183 output[ "bias_type" ] = bias_type;
1184 output[ "num_variables" ] = num_variables;
1185 output[ "num_interactions" ] = num_interactions;
1186 output[ "variable_labels" ] = variables;
1187 output[ "variable_type" ] = variable_type;
1188 output[ "offset" ] = m_offset;
1189 output[ "info" ] = "";
1190 output[ "linear_biases" ] = l_bias;
1191 output[ "quadratic_biases" ] = q_bias;
1192 output[ "quadratic_head" ] = q_head;
1193 output[ "quadratic_tail" ] = q_tail;
1194
1195 return output;
1196 }
1197
1206 template<typename IndexType_serial = IndexType, typename FloatType_serial = FloatType>
1208 // extract type and version
1209 std::string type = input[ "type" ];
1210 if ( type != "BinaryQuadraticModel" ) {
1211 throw std::runtime_error( "Type must be \"BinaryQuadraticModel\".\n" );
1212 }
1213 std::string version = input[ "version" ][ "bqm_schema" ];
1214 if ( version != "3.0.0" ) {
1215 throw std::runtime_error( "bqm_schema must be 3.0.0.\n" );
1216 }
1217
1218 // extract variable_type
1219 Vartype vartype;
1220 std::string variable_type = input[ "variable_type" ];
1221 if ( variable_type == "SPIN" ) {
1222 vartype = Vartype::SPIN;
1223 } else if ( variable_type == "BINARY" ) {
1224 vartype = Vartype::BINARY;
1225 } else {
1226 throw std::runtime_error( "variable_type must be SPIN or BINARY." );
1227 }
1228
1229 // extract linear biases
1231 std::vector<IndexType_serial> variables = input[ "variable_labels" ];
1232 std::vector<FloatType_serial> l_bias = input[ "linear_biases" ];
1233 for ( size_t i = 0; i < variables.size(); ++i ) {
1234 insert_or_assign( linear, variables[ i ], l_bias[ i ] );
1235 }
1236
1237 // extract quadratic biases
1239 std::vector<size_t> q_head = input[ "quadratic_head" ];
1240 std::vector<size_t> q_tail = input[ "quadratic_tail" ];
1241 std::vector<FloatType_serial> q_bias = input[ "quadratic_biases" ];
1242 for ( size_t i = 0; i < q_head.size(); ++i ) {
1243 insert_or_assign( quadratic, std::make_pair( variables[ q_head[ i ] ], variables[ q_tail[ i ] ] ), q_bias[ i ] );
1244 }
1245
1246 // extract offset and info
1247 FloatType_serial offset = input[ "offset" ];
1248
1249 // construct a BinaryQuadraticModel instance
1250 BinaryQuadraticModel<IndexType_serial, FloatType_serial, DataType> bqm( linear, quadratic, offset, vartype );
1251 return bqm;
1252 }
1253 };
1254} // namespace cimod
Linear< IndexType, FloatType > m_linear
Linear biases as a dictionary.
Definition binary_quadratic_model_dict.hpp:65
void add_interaction(const IndexType &arg_u, const IndexType &arg_v, const FloatType &bias)
Add an interaction and/or quadratic bias to a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:374
void remove_variables_from(const std::vector< IndexType > &variables)
Remove specified variables and all of their interactions from a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:460
BinaryQuadraticModel(const Linear< IndexType, FloatType > &linear, const Quadratic< IndexType, FloatType > &quadratic, const FloatType &offset, const Vartype vartype)
BinaryQuadraticModel constructor.
Definition binary_quadratic_model_dict.hpp:94
std::tuple< Linear< IndexType, FloatType >, Quadratic< IndexType, FloatType >, FloatType > to_ising()
Convert a binary quadratic model to Ising format.
Definition binary_quadratic_model_dict.hpp:963
FloatType energy(const Sample< IndexType > &sample) const
Determine the energy of the specified sample of a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:883
Quadratic< IndexType, FloatType > m_quadratic
Quadratic biases as a dictionary.
Definition binary_quadratic_model_dict.hpp:71
const FloatType & get_offset() const
Get the offset.
Definition binary_quadratic_model_dict.hpp:283
FloatType m_offset
The energy offset associated with the model.
Definition binary_quadratic_model_dict.hpp:77
BinaryQuadraticModel(const Eigen::Ref< const DenseMatrix > &, const std::vector< IndexType > &, const Vartype, bool)
BinaryQuadraticModel constructor (with matrix); This constructor is not be implemented.
Definition binary_quadratic_model_dict.hpp:147
Matrix interaction_matrix(const std::vector< IndexType > &indices) const
generate interaction matrix with given list of indices The generated matrix will be the following sym...
Definition binary_quadratic_model_dict.hpp:1026
FloatType get_quadratic(IndexType label_i, IndexType label_j) const
Get the element of quadratic object.
Definition binary_quadratic_model_dict.hpp:265
void normalize(const std::pair< FloatType, FloatType > &bias_range={ 1.0, 1.0 }, const bool use_quadratic_range=false, const std::pair< FloatType, FloatType > &quadratic_range={ 1.0, 1.0 }, const std::vector< IndexType > &ignored_variables={}, const std::vector< std::pair< IndexType, IndexType > > &ignored_interactions={}, const bool ignored_offset=false)
Normalizes the biases of the binary quadratic model such that they fall in the provided range(s),...
Definition binary_quadratic_model_dict.hpp:579
void remove_interactions_from(const std::vector< std::pair< IndexType, IndexType > > &interactions)
Remove all specified interactions from the binary quadratic model.
Definition binary_quadratic_model_dict.hpp:510
void fix_variable(const IndexType &v, const int32_t &value)
Fix the value of a variable and remove it from a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:624
void remove_variable(const IndexType &v)
Remove variable v and all its interactions from a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:444
void flip_variable(const IndexType &v)
Flip variable v in a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:656
BinaryQuadraticModel(const BinaryQuadraticModel &)=default
Copy constructor of BinaryQuadraticModel.
const Linear< IndexType, FloatType > & get_linear() const
Get the linear object.
Definition binary_quadratic_model_dict.hpp:256
std::vector< FloatType > energies(const std::vector< Sample< IndexType > > &samples_like) const
Determine the energies of the given samples.
Definition binary_quadratic_model_dict.hpp:906
void scale(const FloatType &scalar, const std::vector< IndexType > &ignored_variables={}, const std::vector< std::pair< IndexType, IndexType > > &ignored_interactions={}, const bool ignored_offset=false)
Multiply by the specified scalar all the biases and offset of a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:540
size_t length() const
Return the number of variables.
Definition binary_quadratic_model_dict.hpp:212
std::vector< IndexType > _generate_indices() const
generate indices
Definition binary_quadratic_model_dict.hpp:190
void fix_variables(const std::vector< std::pair< IndexType, int32_t > > &fixed)
Fix the value of the variables and remove it from a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:645
static BinaryQuadraticModel from_ising(const Linear< IndexType, FloatType > &linear, const Quadratic< IndexType, FloatType > &quadratic, FloatType offset=0.0)
Create a binary quadratic model from an Ising problem.
Definition binary_quadratic_model_dict.hpp:982
void add_offset(const FloatType &offset)
Add specified value to the offset of a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:521
Matrix interaction_matrix() const
generate interaction matrix with given list of indices The generated matrix will be the following sym...
Definition binary_quadratic_model_dict.hpp:1069
BinaryQuadraticModel(const Eigen::Ref< const DenseMatrix > &, const std::vector< IndexType > &, const FloatType &, const Vartype, bool)
BinaryQuadraticModel constructor (with matrix); This constructor is not be implemented.
Definition binary_quadratic_model_dict.hpp:129
Eigen::Matrix< FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > Matrix
Eigen Matrix.
Definition binary_quadratic_model_dict.hpp:58
std::tuple< Quadratic< IndexType, FloatType >, FloatType > to_qubo()
Convert a binary quadratic model to QUBO format.
Definition binary_quadratic_model_dict.hpp:920
void change_vartype(const Vartype &vartype)
Create a binary quadratic model with the specified vartype.
Definition binary_quadratic_model_dict.hpp:743
nlohmann::json json
Definition binary_quadratic_model_dict.hpp:1098
BinaryQuadraticModel(const SparseMatrix &, const std::vector< IndexType > &, const Vartype)
BinaryQuadraticModel constructor (with sparse matrix); this constructor is for developers.
Definition binary_quadratic_model_dict.hpp:174
static BinaryQuadraticModel< IndexType_serial, FloatType_serial, DataType > from_serializable(const json &input)
Create a BinaryQuadraticModel instance from a serializable object.
Definition binary_quadratic_model_dict.hpp:1207
BinaryQuadraticModel(const Linear< IndexType, FloatType > &linear, const Quadratic< IndexType, FloatType > &quadratic, const Vartype vartype)
BinaryQuadraticModel constructor.
Definition binary_quadratic_model_dict.hpp:112
const std::vector< IndexType > get_variables() const
Get variables.
Definition binary_quadratic_model_dict.hpp:301
json to_serializable() const
Convert the binary quadratic model to a serializable object user_bytes is assume to be set to False.
Definition binary_quadratic_model_dict.hpp:1106
static std::tuple< Linear< IndexType, FloatType >, Quadratic< IndexType, FloatType >, FloatType > spin_to_binary(const Linear< IndexType, FloatType > &linear, const Quadratic< IndexType, FloatType > &quadratic, const FloatType &offset)
Convert linear, quadratic, and offset from spin to binary.
Definition binary_quadratic_model_dict.hpp:811
void add_variable(const IndexType &v, const FloatType &bias)
Add variable v and/or its bias to a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:330
void remove_interaction(const IndexType &arg_u, const IndexType &arg_v)
Remove interaction of variables u, v from a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:472
size_t get_num_variables() const
Return the number of variables.
Definition binary_quadratic_model_dict.hpp:222
void add_interactions_from(const Quadratic< IndexType, FloatType > &quadratic)
Add interactions and/or quadratic biases to a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:433
bool contains(const IndexType &v) const
Return true if the variable contains v.
Definition binary_quadratic_model_dict.hpp:232
BinaryQuadraticModel(const SparseMatrix &, const std::vector< IndexType > &, const FloatType &, const Vartype)
BinaryQuadraticModel constructor (with sparse matrix); this constructor is for developers.
Definition binary_quadratic_model_dict.hpp:161
static std::tuple< Linear< IndexType, FloatType >, Quadratic< IndexType, FloatType >, FloatType > binary_to_spin(const Linear< IndexType, FloatType > &linear, const Quadratic< IndexType, FloatType > &quadratic, const FloatType &offset)
Convert linear, quadratic and offset from binary to spin.
Definition binary_quadratic_model_dict.hpp:848
const Vartype & get_vartype() const
Get the vartype object.
Definition binary_quadratic_model_dict.hpp:292
Eigen::SparseMatrix< FloatType, Eigen::RowMajor > SparseMatrix
Definition binary_quadratic_model_dict.hpp:53
FloatType get_linear(IndexType label_i) const
Get the element of linear object.
Definition binary_quadratic_model_dict.hpp:247
void remove_offset()
Set the binary quadratic model's offset to zero.
Definition binary_quadratic_model_dict.hpp:528
void contract_variables(const IndexType &u, const IndexType &v)
Enforce u, v being the same variable in a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:691
void add_variables_from(const Linear< IndexType, FloatType > &linear)
Add variables and/or linear biases to a binary quadratic model.
Definition binary_quadratic_model_dict.hpp:361
Eigen::Matrix< FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > DenseMatrix
Definition binary_quadratic_model_dict.hpp:52
BinaryQuadraticModel< IndexType, FloatType, DataType > empty(Vartype vartype)
Create an empty BinaryQuadraticModel.
Definition binary_quadratic_model_dict.hpp:317
BinaryQuadraticModel change_vartype(const Vartype &vartype, bool inplace)
Create a binary quadratic model with the specified vartype.
Definition binary_quadratic_model_dict.hpp:772
static BinaryQuadraticModel from_qubo(const Quadratic< IndexType, FloatType > &Q, FloatType offset=0.0)
Create a binary quadratic model from a QUBO model.
Definition binary_quadratic_model_dict.hpp:941
const Quadratic< IndexType, FloatType > & get_quadratic() const
Get the quadratic object.
Definition binary_quadratic_model_dict.hpp:274
Class for dense binary quadratic model.
Definition binary_quadratic_model.hpp:148
size_t get_num_variables() const
get the number of variables
Definition binary_quadratic_model.hpp:1268
void add_variables_from(const Linear< IndexType, FloatType > &linear)
Add variables and/or linear biases to a binary quadratic model.
Definition binary_quadratic_model.hpp:1388
void add_interaction(const IndexType &u, const IndexType &v, const FloatType &bias)
Add an interaction and/or quadratic bias to a binary quadratic model.
Definition binary_quadratic_model.hpp:1401
Vartype get_vartype() const
Get the vartype object.
Definition binary_quadratic_model.hpp:1346
const std::vector< IndexType > & get_variables() const
Get variables.
Definition binary_quadratic_model.hpp:1355
FloatType m_offset
The energy offset associated with the model.
Definition binary_quadratic_model.hpp:207
void remove_interaction(const IndexType &u, const IndexType &v)
Remove interaction of variables u, v from a binary quadratic model.
Definition binary_quadratic_model.hpp:1445
void scale(const FloatType &scalar, const std::vector< IndexType > &ignored_variables={}, const std::vector< std::pair< IndexType, IndexType > > &ignored_interactions={}, const bool ignored_offset=false)
Multiply by the specified scalar all the biases and offset of a binary quadratic model.
Definition binary_quadratic_model.hpp:1486
void add_offset(const FloatType &offset)
Add specified value to the offset of a binary quadratic model.
Definition binary_quadratic_model.hpp:1467
void add_variable(const IndexType &v, const FloatType &bias)
Add variable v and/or its bias to a binary quadratic model.
Definition binary_quadratic_model.hpp:1377
Vartype m_vartype
The model's type.
Definition binary_quadratic_model.hpp:213
void remove_variable(const IndexType &v)
Remove variable v and all its interactions from a binary quadratic model.
Definition binary_quadratic_model.hpp:1424
void fix_variable(const IndexType &v, const int32_t &value)
Fix the value of a variable and remove it from a binary quadratic model.
Definition binary_quadratic_model.hpp:1564
FloatType energy(const Sample< IndexType > &sample) const
Determine the energy of the specified sample of a binary quadratic model.
Definition binary_quadratic_model.hpp:1686
void add_interactions_from(const Quadratic< IndexType, FloatType > &quadratic)
Add interactions and/or quadratic biases to a binary quadratic model.
Definition binary_quadratic_model.hpp:1413
void change_vartype(const Vartype &vartype)
‍**
Definition binary_quadratic_model.hpp:1649
void remove_interactions_from(const std::vector< std::pair< IndexType, IndexType > > &interactions)
Remove all specified interactions from the binary quadratic model.
Definition binary_quadratic_model.hpp:1456
void declare_BQM(py::module &m, const std::string &name)
Definition main.hpp:41
Definition binary_polynomial_model.hpp:139
bool check_vartype(const int32_t &var, const Vartype &vartype)
Check that the variable has appropriate value.
Definition vartypes.hpp:37
std::unordered_map< IndexType, int32_t > Sample
Type alias for sample, which represents the spin or binary configurations.
Definition binary_polynomial_model.hpp:162
void insert_or_assign(std::unordered_map< C_key, C_value, Hash > &um, const C_key &key, const C_value &val)
Insert or assign a element of unordered_map (for C++14 or C++11)
Definition utilities.hpp:34
std::unordered_map< IndexType, FloatType > Linear
Type alias for linear bias.
Definition binary_quadratic_model.hpp:111
std::unordered_map< std::pair< IndexType, IndexType >, FloatType, pair_hash > Quadratic
Type alias for quadratic bias.
Definition binary_quadratic_model.hpp:119
Vartype
Enum class for representing problem type.
Definition vartypes.hpp:24
Definition binary_quadratic_model_dict.hpp:41