Binary Quadratic Model Dict

Contents

Binary Quadratic Model Dict#

namespace cimod
struct Dict#
template<typename IndexType, typename FloatType>
class BinaryQuadraticModel<IndexType, FloatType, Dict>#

Class for legacy binary quadratic model (dict datastructure).

Public Types

using DenseMatrix = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>#
using SparseMatrix = Eigen::SparseMatrix<FloatType, Eigen::RowMajor>#
using Matrix = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>#

Eigen Matrix.

using json = nlohmann::json#

Public Functions

inline BinaryQuadraticModel(const Linear<IndexType, FloatType> &linear, const Quadratic<IndexType, FloatType> &quadratic, const FloatType &offset, const Vartype vartype)#

BinaryQuadraticModel constructor.

Parameters:
  • linear

  • quadratic

  • offset

  • vartype

inline BinaryQuadraticModel(const Linear<IndexType, FloatType> &linear, const Quadratic<IndexType, FloatType> &quadratic, const Vartype vartype)#

BinaryQuadraticModel constructor.

Parameters:
  • linear

  • quadratic

  • vartype

inline 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.

Parameters:
  • mat

  • labels_vec

  • offset

  • vartype

inline BinaryQuadraticModel(const Eigen::Ref<const DenseMatrix>&, const std::vector<IndexType>&, const Vartype, bool)#

BinaryQuadraticModel constructor (with matrix); This constructor is not be implemented.

Parameters:
  • mat

  • labels_vec

  • vartype

inline BinaryQuadraticModel(const SparseMatrix&, const std::vector<IndexType>&, const FloatType&, const Vartype)#

BinaryQuadraticModel constructor (with sparse matrix); this constructor is for developers.

Parameters:
  • mat

  • labels_vec

  • offset

  • vartype

inline BinaryQuadraticModel(const SparseMatrix&, const std::vector<IndexType>&, const Vartype)#

BinaryQuadraticModel constructor (with sparse matrix); this constructor is for developers.

Parameters:
  • mat

  • labels_vec

  • vartype

BinaryQuadraticModel(const BinaryQuadraticModel&) = default#

Copy constructor of BinaryQuadraticModel.

Parameters:

bqm

inline std::vector<IndexType> _generate_indices() const#

generate indices

Returns:

generated indices

inline size_t length() const#

Return the number of variables.

Deprecated:

use get_num_variables.

Returns:

The number of variables.

inline size_t get_num_variables() const#

Return the number of variables.

Deprecated:

use get_num_variables.

Returns:

The number of variables.

inline bool contains(const IndexType &v) const#

Return true if the variable contains v.

Parameters:

v

Returns:

Return true if the variable contains v.

inline FloatType get_linear(IndexType label_i) const#

Get the element of linear object.

Parameters:

label_i

Returns:

inline const Linear<IndexType, FloatType> &get_linear() const#

Get the linear object.

Returns:

A linear bias.

inline FloatType get_quadratic(IndexType label_i, IndexType label_j) const#

Get the element of quadratic object.

Returns:

A quadratic bias.

inline const Quadratic<IndexType, FloatType> &get_quadratic() const#

Get the quadratic object.

Returns:

A quadratic bias.

inline const FloatType &get_offset() const#

Get the offset.

Returns:

An offset.

inline const Vartype &get_vartype() const#

Get the vartype object.

Returns:

Type of the model.

inline const std::vector<IndexType> get_variables() const#

Get variables.

Returns:

variables

inline BinaryQuadraticModel<IndexType, FloatType, DataType> empty(Vartype vartype)#

Create an empty BinaryQuadraticModel.

inline void add_variable(const IndexType &v, const FloatType &bias)#

Add variable v and/or its bias to a binary quadratic model.

Parameters:
  • v

  • bias

inline void add_variables_from(const Linear<IndexType, FloatType> &linear)#

Add variables and/or linear biases to a binary quadratic model.

Parameters:

linear

inline 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.

Parameters:
  • arg_u

  • arg_v

  • bias

inline void add_interactions_from(const Quadratic<IndexType, FloatType> &quadratic)#

Add interactions and/or quadratic biases to a binary quadratic model.

Parameters:

quadratic

inline void remove_variable(const IndexType &v)#

Remove variable v and all its interactions from a binary quadratic model.

Parameters:

v

inline void remove_variables_from(const std::vector<IndexType> &variables)#

Remove specified variables and all of their interactions from a binary quadratic model.

Parameters:

variables

inline void remove_interaction(const IndexType &arg_u, const IndexType &arg_v)#

Remove interaction of variables u, v from a binary quadratic model.

Parameters:
  • arg_u

  • arg_v

inline void remove_interactions_from(const std::vector<std::pair<IndexType, IndexType>> &interactions)#

Remove all specified interactions from the binary quadratic model.

Parameters:

interactions

inline void add_offset(const FloatType &offset)#

Add specified value to the offset of a binary quadratic model.

Parameters:

offset

inline void remove_offset()#

Set the binary quadratic model’s offset to zero.

inline 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.

Parameters:
  • scalar

  • ignored_variables

  • ignored_interactions

  • ignored_offset

inline 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), and adjusts the offset appropriately.

Parameters:
  • bias_range

  • use_quadratic_range

  • quadratic_range

  • ignored_variables

  • ignored_interactions

  • ignored_offset

inline void fix_variable(const IndexType &v, const int32_t &value)#

Fix the value of a variable and remove it from a binary quadratic model.

Parameters:
  • v

  • value

inline 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.

Parameters:

fixed

inline void flip_variable(const IndexType &v)#

Flip variable v in a binary quadratic model.

Parameters:

v

inline void contract_variables(const IndexType &u, const IndexType &v)#

Enforce u, v being the same variable in a binary quadratic model.

Parameters:
  • u

  • v

inline void change_vartype(const Vartype &vartype)#

Create a binary quadratic model with the specified vartype. This function does not return any object.

Parameters:

vartype

inline BinaryQuadraticModel change_vartype(const Vartype &vartype, bool inplace)#

Create a binary quadratic model with the specified vartype.

Parameters:
  • vartype

  • inplace

Returns:

A new instance of the BinaryQuadraticModel class.

inline FloatType energy(const Sample<IndexType> &sample) const#

Determine the energy of the specified sample of a binary quadratic model.

Parameters:

sample

Returns:

An energy with respect to the sample.

inline std::vector<FloatType> energies(const std::vector<Sample<IndexType>> &samples_like) const#

Determine the energies of the given samples.

Parameters:

samples_like

Returns:

A vector including energies with respect to the samples.

inline std::tuple<Quadratic<IndexType, FloatType>, FloatType> to_qubo()#

Convert a binary quadratic model to QUBO format.

Returns:

A tuple including a quadratic bias and an offset.

inline std::tuple<Linear<IndexType, FloatType>, Quadratic<IndexType, FloatType>, FloatType> to_ising()#

Convert a binary quadratic model to Ising format.

Returns:

A tuple including a linear bias, a quadratic bias and an offset.

inline Matrix interaction_matrix(const std::vector<IndexType> &indices) const#

generate interaction matrix with given list of indices The generated matrix will be the following symmetric matrix:

(J~0,0J~0,1J~0,2J~1,0J~1,1J~1,2J~2,0J~2,1J~2,2) \begin{pmatrix} \tilde{J}_{0,0} & \tilde{J}_{0,1} & \tilde{J}_{0,2} & \cdots \\ \tilde{J}_{1,0} & \tilde{J}_{1,1} & \tilde{J}_{1,2} & \cdots \\ \tilde{J}_{2,0} & \tilde{J}_{2,1} & \tilde{J}_{2,2} & \cdots \\ \vdots & \vdots & \vdots & \ddots \\ \end{pmatrix}

where in the Ising case,

J~f(i),f(j)=Jij+Jji,J~f(i),f(i)=hi, \tilde{J}_{f(i),f(j)} = J_{ij} + J_{ji}, \\ \tilde{J}_{f(i),f(i)} = h_{i},
and the QUBO case,
J~f(i),f(j)=Qij+Qji,J~f(i),f(i)=Qii, \tilde{J}_{f(i),f(j)} = Q_{ij} + Q_{ji}, \\ \tilde{J}_{f(i),f(i)} = Q_{ii},
and the function ff denotes a mapping from index to the corresponding number specified by the argument indices. For instance, if indices is [‘a’, ‘b’, ‘c’], The following equations, f(a)=0,f(b)=1,andf(c)=2f(a) = 0, f(b)=1, \mathrm{and} f(c)=2 hold.

The original Hamiltonian can be rewritten with Jij~\tilde{J_{ij}} as

EIsing=iJ~f(i),f(i)si+i<jJ~f(i),f(j)sisj+δIsing, E_{\mathrm{Ising}} = \sum_{i} \tilde{J}_{f(i),f(i)} s_i + \sum_{i < j} \tilde{J}_{f(i), f(j)} s_i s_j + \delta_{\mathrm{Ising}},
and
EQUBO=iJ~f(i),f(i)xi+i<jJ~f(i),f(j)xixj+δQUBO. E_{\mathrm{QUBO}} = \sum_{i} \tilde{J}_{f(i), f(i)}x_i + \sum_{i < j} \tilde{J}_{f(i), f(j)} x_i x_j + \delta_{\mathrm{QUBO}}.
Deprecated:

use interaction_matrix()

Parameters:

indices

Returns:

corresponding interaction matrix (Eigen)

inline Matrix interaction_matrix() const#

generate interaction matrix with given list of indices The generated matrix will be the following symmetric matrix:

(J0,0J0,1J0,N1h00J1,1J1,N1h100JN1,N1hN10001) \begin{pmatrix} J_{0,0} & J_{0,1} & \cdots & J_{0,N-1} & h_{0}\\ 0 & J_{1,1} & \cdots & J_{1,N-1} & h_{1}\\ \vdots & \vdots & \vdots & \vdots & \vdots \\ 0 & 0 & \cdots & J_{N-1,N-1} & h_{N-1}\\ 0 & 0 & \cdots & 0 & 1 \\ \end{pmatrix}

inline json to_serializable() const#

Convert the binary quadratic model to a serializable object user_bytes is assume to be set to False.

Returns:

An object that can be serialized (nlohmann::json)

Public Static Functions

static inline 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. Does no checking of vartype. Copies all of the values into new objects.

Parameters:
  • linear

  • quadratic

  • offset

Returns:

A tuple including a linear bias, a quadratic bias and an offset.

static inline 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. Does no checking of vartype. Copies all of the values into new objects.

Parameters:
  • linear

  • quadratic

  • offset

Returns:

A tuple including a linear bias, a quadratic bias and an offset.

static inline BinaryQuadraticModel from_qubo(const Quadratic<IndexType, FloatType> &Q, FloatType offset = 0.0)#

Create a binary quadratic model from a QUBO model.

Parameters:
  • Q

  • offset

Returns:

Binary quadratic model with vartype set to .Vartype.BINARY.

static inline 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.

Parameters:
  • linear

  • quadratic

  • offset

Returns:

Binary quadratic model with vartype set to .Vartype.SPIN.

template<typename IndexType_serial = IndexType, typename FloatType_serial = FloatType>
static inline BinaryQuadraticModel<IndexType_serial, FloatType_serial, DataType> from_serializable(const json &input)#

Create a BinaryQuadraticModel instance from a serializable object.

Template Parameters:
  • IndexType_serial

  • FloatType_serial

Parameters:

input

Returns:

BinaryQuadraticModel<IndexType_serial, FloatType_serial, DataType>

Protected Attributes

Linear<IndexType, FloatType> m_linear#

Linear biases as a dictionary.

Quadratic<IndexType, FloatType> m_quadratic#

Quadratic biases as a dictionary.

FloatType m_offset#

The energy offset associated with the model.

Vartype m_vartype = Vartype::NONE#

The model’s type.

Private Types

using DataType = Dict#