openjij
Framework for the Ising model and QUBO.
Loading...
Searching...
No Matches
openjij::graph::Dense< FloatType > Class Template Reference

two-body all-to-all interactions The Hamiltonian is like More...

#include <dense.hpp>

Inheritance diagram for openjij::graph::Dense< FloatType >:
Collaboration diagram for openjij::graph::Dense< FloatType >:

Public Types

using Interactions = Eigen::Matrix< FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >
 interaction type (Eigen) The stored matrix has the following triangular form:
 
using value_type = FloatType
 float type
 

Public Member Functions

 Dense (std::size_t num_spins)
 Dense constructor.
 
 Dense (const json &j)
 Dense constructor (from nlohmann::json)
 
void set_interaction_matrix (const Interactions &interaction)
 set interaction matrix from Eigen Matrix.
 
 Dense (const Dense< FloatType > &)=default
 Dense copy constructor.
 
 Dense (Dense< FloatType > &&)=default
 Dense move constructor.
 
FloatType calc_energy (const Spins &spins) const
 calculate total energy
 
FloatType calc_energy (const Eigen::Matrix< FloatType, Eigen::Dynamic, 1, Eigen::ColMajor > &spins) const
 
FloatType energy (const Spins &spins) const
 calculate total energy
 
FloatType energy (const Eigen::Matrix< FloatType, Eigen::Dynamic, 1, Eigen::ColMajor > &spins) const
 
FloatTypeJ (Index i, Index j)
 access J_{ij}
 
const FloatTypeJ (Index i, Index j) const
 access J_{ij}
 
FloatTypeh (Index i)
 access h_{i} (local field)
 
const FloatTypeh (Index i) const
 access h_{i} (local field)
 
const Interactions get_interactions () const
 get interactions (Eigen Matrix)
 
- Public Member Functions inherited from openjij::graph::Graph
 Graph (std::size_t num_spins)
 Graph constructor.
 
template<typename RandomNumberEngine >
const Spins gen_spin (RandomNumberEngine &random_numder_engine) const
 generate spins randomly.
 
template<typename RandomNumberEngine >
const Binaries gen_binary (RandomNumberEngine &random_numder_engine) const
 generate spins randomly.
 
std::size_t get_num_spins () const noexcept
 get number of spins
 
std::size_t size () const noexcept
 get number of spins
 

Private Attributes

Interactions _J
 interactions
 

Detailed Description

template<typename FloatType>
class openjij::graph::Dense< FloatType >

two-body all-to-all interactions The Hamiltonian is like

\[ H = \sum_{i<j}J_{ij} \sigma_i \sigma_j + \sum_{i}h_{i} \sigma_i \]

Template Parameters
FloatTypefloat type of Sparse class (double or float)

Member Typedef Documentation

◆ Interactions

template<typename FloatType >
using openjij::graph::Dense< FloatType >::Interactions = Eigen::Matrix<FloatType, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>

interaction type (Eigen) The stored matrix has the following triangular form:

\[ \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} \]

◆ value_type

float type

Constructor & Destructor Documentation

◆ Dense() [1/4]

template<typename FloatType >
openjij::graph::Dense< FloatType >::Dense ( std::size_t  num_spins)
inlineexplicit

Dense constructor.

Parameters
num_spinsthe number of spins

References openjij::graph::Dense< FloatType >::_J.

◆ Dense() [2/4]

template<typename FloatType >
openjij::graph::Dense< FloatType >::Dense ( const json j)
inline

Dense constructor (from nlohmann::json)

Parameters
jJSON object this object must be a serialized object of dense BQM.

References openjij::graph::json_parse().

Here is the call graph for this function:

◆ Dense() [3/4]

template<typename FloatType >
openjij::graph::Dense< FloatType >::Dense ( const Dense< FloatType > &  )
default

Dense copy constructor.

◆ Dense() [4/4]

template<typename FloatType >
openjij::graph::Dense< FloatType >::Dense ( Dense< FloatType > &&  )
default

Dense move constructor.

Member Function Documentation

◆ calc_energy() [1/2]

template<typename FloatType >
FloatType openjij::graph::Dense< FloatType >::calc_energy ( const Eigen::Matrix< FloatType, Eigen::Dynamic, 1, Eigen::ColMajor > &  spins) const
inline

References openjij::graph::Dense< FloatType >::energy(), and openjij::graph::json_parse().

Here is the call graph for this function:

◆ calc_energy() [2/2]

template<typename FloatType >
FloatType openjij::graph::Dense< FloatType >::calc_energy ( const Spins spins) const
inline

calculate total energy

Parameters
spins
Deprecated:
use energy(spins)
Returns
corresponding energy

References openjij::graph::Dense< FloatType >::energy().

Referenced by openjij::declare_Dense().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ energy() [1/2]

template<typename FloatType >
FloatType openjij::graph::Dense< FloatType >::energy ( const Eigen::Matrix< FloatType, Eigen::Dynamic, 1, Eigen::ColMajor > &  spins) const
inline

◆ energy() [2/2]

template<typename FloatType >
FloatType openjij::graph::Dense< FloatType >::energy ( const Spins spins) const
inline

calculate total energy

Parameters
spins
Returns
corresponding energy

References openjij::graph::Dense< FloatType >::_J, openjij::graph::Graph::get_num_spins(), and openjij::graph::json_parse().

Referenced by openjij::graph::Dense< FloatType >::calc_energy(), openjij::graph::Dense< FloatType >::calc_energy(), and openjij::graph::Dense< FloatType >::energy().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_interactions()

template<typename FloatType >
const Interactions openjij::graph::Dense< FloatType >::get_interactions ( ) const
inline

get interactions (Eigen Matrix)

The returned matrix has the following symmetric form:

\[ \begin{pmatrix} J_{0,0} & J_{0,1} & \cdots & J_{0,N-1} & h_{0}\\ J_{0,1} & J_{1,1} & \cdots & J_{1,N-1} & h_{1}\\ \vdots & \vdots & \vdots & \vdots & \vdots \\ J_{0,N-1} & J_{N-1,1} & \cdots & J_{N-1,N-1} & h_{N-1}\\ h_{0} & h_{1} & \cdots & h_{N-1} & 1 \\ \end{pmatrix} \]

Returns
Eigen Matrix

References openjij::graph::json_parse().

Here is the call graph for this function:

◆ h() [1/2]

access h_{i} (local field)

Parameters
iIndex i
Returns
h_{i}

References openjij::graph::Graph::get_num_spins(), openjij::graph::Dense< FloatType >::J(), and openjij::graph::json_parse().

Referenced by openjij::declare_Dense(), and openjij::utility::gen_matrix_from_graph().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ h() [2/2]

template<typename FloatType >
const FloatType & openjij::graph::Dense< FloatType >::h ( Index  i) const
inline

access h_{i} (local field)

Parameters
iIndex i
Returns
h_{i}

References openjij::graph::Graph::get_num_spins(), openjij::graph::Dense< FloatType >::J(), and openjij::graph::json_parse().

Here is the call graph for this function:

◆ J() [1/2]

template<typename FloatType >
FloatType & openjij::graph::Dense< FloatType >::J ( Index  i,
Index  j 
)
inline

access J_{ij}

Parameters
iIndex i
jIndex j
Returns
J_{ij}

References openjij::graph::Dense< FloatType >::_J, openjij::graph::Graph::get_num_spins(), and openjij::graph::json_parse().

Referenced by openjij::declare_Dense(), openjij::utility::gen_matrix_from_graph(), openjij::graph::Dense< FloatType >::h(), and openjij::graph::Dense< FloatType >::h().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ J() [2/2]

template<typename FloatType >
const FloatType & openjij::graph::Dense< FloatType >::J ( Index  i,
Index  j 
) const
inline

access J_{ij}

Parameters
iIndex i
jIndex j
Returns
J_{ij}

References openjij::graph::Dense< FloatType >::_J, openjij::graph::Graph::get_num_spins(), and openjij::graph::json_parse().

Here is the call graph for this function:

◆ set_interaction_matrix()

template<typename FloatType >
void openjij::graph::Dense< FloatType >::set_interaction_matrix ( const Interactions interaction)
inline

set interaction matrix from Eigen Matrix.

Parameters
interactionEigen matrix

References openjij::graph::Dense< FloatType >::_J, openjij::graph::Graph::get_num_spins(), and openjij::graph::json_parse().

Here is the call graph for this function:

Member Data Documentation

◆ _J


The documentation for this class was generated from the following file: