netket.operator.Ising#

class netket.operator.Ising[source]#

Bases: IsingBase

The Transverse-Field Ising Hamiltonian \(-h\sum_i \sigma_i^{(x)} +J\sum_{\langle i,j\rangle} \sigma_i^{(z)}\sigma_j^{(z)}\).

This implementation is considerably faster than the Ising hamiltonian constructed by summing LocalOperator s.

Inheritance
Inheritance diagram of netket.operator.Ising
__init__(hilbert, graph, h, J, dtype)[source]#

Constructs the Ising Operator from an hilbert space and a graph specifying the connectivity.

Parameters:
  • hilbert (AbstractHilbert) – Hilbert space the operator acts on.

  • h (float) – The strength of the transverse field.

  • J (float) – The strength of the coupling. Default is 1.0.

  • dtype (Any | None) – The dtype of the matrix elements.

  • graph (AbstractGraph | ndarray | Array)

Examples

Constructs an Ising operator for a 1D system.

>>> import netket as nk
>>> g = nk.graph.Hypercube(length=20, n_dim=1, pbc=True)
>>> hi = nk.hilbert.Spin(s=0.5, N=g.n_nodes)
>>> op = nk.operator.Ising(h=1.321, hilbert=hi, J=0.5, graph=g)
>>> print(op)
Ising(J=0.5, h=1.321; dim=20)
Attributes
H#

Returns the Conjugate-Transposed operator

J#

The magnitude of the hopping

T#

Returns the transposed operator

dtype#

The dtype of the matrix elements.

edges#

The (N_conns, 2) matrix of edges on which the interaction term is non-zero.

h#

The magnitude of the transverse field

hilbert#

The hilbert space associated to this observable.

is_hermitian#

A boolean stating whether this hamiltonian is hermitian.

max_conn_size#

The maximum number of non zero ⟨x|O|x’⟩ for every x.

Methods
__call__(v)[source]#

Call self as a function.

Return type:

ndarray

Parameters:

v (ndarray)

apply(v)[source]#
Return type:

ndarray

Parameters:

v (ndarray)

collect()[source]#

Returns a guaranteed concrete instance of an operator.

As some operations on operators return lazy wrappers (such as transpose, hermitian conjugate…), this is used to obtain a guaranteed non-lazy operator.

Return type:

AbstractOperator

conj(*, concrete=False)[source]#
Return type:

AbstractOperator

conjugate(*, concrete=True)[source]#

Returns the complex-conjugate of this operator.

Parameters:

concrete – if True returns a concrete operator and not a lazy wrapper

Returns:

if concrete is not True, self or a lazy wrapper; the complex-conjugated operator otherwise

copy(*, dtype=None)[source]#
Parameters:

dtype (Any | None)

get_conn(x)[source]#

Finds the connected elements of the Operator. Starting from a given quantum number x, it finds all other quantum numbers x’ such that the matrix element \(O(x,x')\) is different from zero. In general there will be several different connected states x’ satisfying this condition, and they are denoted here \(x'(k)\), for \(k=0,1...N_{\mathrm{connected}}\).

Parameters:

x (ndarray) – An array of shape (hilbert.size, ) containing the quantum numbers x.

Returns:

The connected states x’ of shape (N_connected,hilbert.size) array: An array containing the matrix elements \(O(x,x')\) associated to each x’.

Return type:

matrix

Raises:

ValueError – If the given quantum number is not compatible with the hilbert space.

get_conn_flattened(x, sections, pad=False)[source]#

Finds the connected elements of the Operator. Starting from a given quantum number x, it finds all other quantum numbers x’ such that the matrix element \(O(x,x')\) is different from zero. In general there will be several different connected states x’ satisfying this condition, and they are denoted here \(x'(k)\), for \(k=0,1...N_{\mathrm{connected}}\).

This is a batched version, where x is a matrix of shape (batch_size,hilbert.size).

Parameters:
  • x (matrix) – A matrix of shape (batch_size,hilbert.size) containing the batch of quantum numbers x.

  • sections (array) – An array of size (batch_size) useful to unflatten the output of this function. See numpy.split for the meaning of sections.

  • pad (bool) – no effect here

Returns:

The connected states x’, flattened together in a single matrix. array: An array containing the matrix elements \(O(x,x')\) associated to each x’.

Return type:

matrix

get_conn_padded(x)[source]#

Finds the connected elements of the Operator.

Starting from a batch of quantum numbers \(x={x_1, ... x_n}\) of size \(B \times M\) where \(B\) size of the batch and \(M\) size of the hilbert space, finds all states \(y_i^1, ..., y_i^K\) connected to every \(x_i\).

Returns a matrix of size \(B \times K_{max} \times M\) where \(K_{max}\) is the maximum number of connections for every \(y_i\).

Parameters:

x (ndarray) – A N-tensor of shape \((...,hilbert.size)\) containing the batch/batches of quantum numbers \(x\).

Returns:

The connected states x’, in a N+1-tensor and an N-tensor containing the matrix elements \(O(x,x')\) associated to each x’ for every batch.

Return type:

(x_primes, mels)

n_conn(x, out=None)[source]#

Return the number of states connected to x.

Parameters:
  • x (matrix) – A matrix of shape (batch_size,hilbert.size) containing the batch of quantum numbers x.

  • out (array) – If None an output array is allocated.

Returns:

The number of connected states x’ for each x[i].

Return type:

array

to_dense()[source]#

Returns the dense matrix representation of the operator. Note that, in general, the size of the matrix is exponential in the number of quantum numbers, and this operation should thus only be performed for low-dimensional Hilbert spaces or sufficiently sparse operators.

This method requires an indexable Hilbert space.

Return type:

ndarray

Returns:

The dense matrix representation of the operator as a Numpy array.

to_jax_operator()[source]#

Returns the jax-compatible version of this operator, which is an instance of netket.operator.IsingJax.

Return type:

IsingJax

to_linear_operator()[source]#
to_local_operator()[source]#
to_qobj()[source]#

Convert the operator to a qutip’s Qobj.

Returns:

A qutip.Qobj object.

to_sparse()[source]#

Returns the sparse matrix representation of the operator. Note that, in general, the size of the matrix is exponential in the number of quantum numbers, and this operation should thus only be performed for low-dimensional Hilbert spaces or sufficiently sparse operators.

This method requires an indexable Hilbert space.

Return type:

csr_matrix

Returns:

The sparse matrix representation of the operator.

transpose(*, concrete=False)[source]#

Returns the transpose of this operator.

Parameters:

concrete – if True returns a concrete operator and not a lazy wrapper

Return type:

AbstractOperator

Returns:

if concrete is not True, self or a lazy wrapper; the transposed operator otherwise