netket.operator.LocalLiouvillian#

class netket.operator.LocalLiouvillian[source]#

Bases: AbstractSuperOperator

LocalLiouvillian super-operator, acting on the DoubledHilbert (tensor product) space ℋ⊗ℋ.

Internally it uses netket.operator.LocalOperator everywhere.

The Liouvillian is defined according to the definition:

\[\mathcal{L} = -i \left[ \hat{H}, \hat{\rho}\right] + \sum_i \left[ \hat{L}_i\hat{\rho}\hat{L}_i^\dagger - \left\{ \hat{L}_i^\dagger\hat{L}_i, \hat{\rho} \right\} \right]\]

which generates the dynamics according to the equation

\[\frac{d\hat{\rho}}{dt} = \mathcal{L}\hat{\rho}\]

Internally, it stores the non-hermitian hamiltonian

\[\hat{H}_{nh} = \hat{H} - \sum_i \frac{i}{2}\hat{L}_i^\dagger\hat{L}_i\]

That is then composed with the jump operators in the inner kernel with the formula:

\[\mathcal{L} = -i \hat{H}_{nh}\hat{\rho} +i\hat{\rho}\hat{H}_{nh}^\dagger + \sum_i \hat{L}_i\hat{\rho}\hat{L}_i^\dagger\]
Inheritance
Inheritance diagram of netket.operator.LocalLiouvillian
Attributes
H#

Returns the Conjugate-Transposed operator

T#

Returns the transposed operator

dtype#
hamiltonian#

The hamiltonian of this Liouvillian

hamiltonian_nh#

The non hermitian Local Operator part of the Liouvillian

hilbert#

The hilbert space associated to this observable.

hilbert_physical#

The physical hilbert space on which this super-operator acts.

is_hermitian#
jump_operators#

The list of local operators in this Liouvillian

max_conn_size#

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

Methods
__call__(v)#

Call self as a function.

Return type:

ndarray

Parameters:

v (ndarray)

add_jump_operator(op)[source]#
add_jump_operators(ops)[source]#
apply(v)#
Return type:

ndarray

Parameters:

v (ndarray)

collect()#

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)#
Return type:

AbstractOperator

conjugate(*, concrete=False)#

Returns the complex-conjugate 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 complex-conjugated operator otherwise

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 – 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 – A matrix of shape (batch_size, hilbert.size) containing the batch of quantum numbers x.

  • sections – An array of sections for the flattened x’. See numpy.split for the meaning of sections.

Returns:

The connected states x’, flattened together in

a single matrix. An array containing the matrix elements \(O(x,x')\) associated to each x’.

Return type:

(matrix, array)

get_conn_padded(x)#

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)#

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()#

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_linear_operator(*, sparse=True, append_trace=False)[source]#

Returns a lazy scipy linear_operator representation of the Lindblad Super-Operator.

The returned operator behaves like the M**2 x M**2 matrix obtained with to_dense/sparse, and accepts vectorised density matrices as input.

Parameters:
  • sparse (bool) – If True internally uses sparse matrices for the hamiltonian and jump operators, dense otherwise (default=True)

  • append_trace (bool) – If True (default=False) the resulting operator has size M**2 + 1, and the last element of the input vector is the trace of the input density matrix. This is useful when implementing iterative methods.

Return type:

LinearOperator

Returns:

A linear operator taking as input vectorised density matrices and returning the product L*rho

to_qobj()[source]#

Convert the operator to a qutip’s liouvillian Qobj.

Returns:

A qutip.liouvillian object.

to_sparse()#

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)#

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