quspin.operators.quantum_LinearOperator
- class quspin.operators.quantum_LinearOperator(*args, **kwargs)[source]
Bases:
LinearOperator
Applies a quantum operator directly onto a state, without constructing the operator matrix.
The quantum_LinearOperator class uses the basis.Op() function to calculate the matrix vector product on the fly, greatly reducing the amount of memory needed for a calculation at the cost of speed.
This object is useful for doing large scale Lanczos calculations using the eigsh method.
Notes
The class does NOT yet support time-dependent operators.
Examples
The following example shows how to construct and use quantum_LinearOperator objects.
1from quspin.basis import spin_basis_1d # Hilbert space spin basis 2import numpy as np # generic math functions 3 4# 5##### define model parameters ##### 6L = 4 # system size 7J = 1.0 # spin interaction 8g = 0.809 # transverse field 9h = 0.9045 # parallel field 10# 11##### construct basis 12basis = spin_basis_1d(L=L) 13# define PBC site-coupling lists for operators 14x_field = [[g, i] for i in range(L)] 15z_field = [[h, i] for i in range(L)] 16J_nn = [[J, i, (i + 1) % L] for i in range(L)] # PBC 17# static list 18static = [["zz", J_nn], ["z", z_field], ["x", x_field]] 19# 20##### construct Hamiltonian, NO dynamic list passed 21H = quantum_LinearOperator(static, basis=basis, dtype=np.float64) 22# 23##### apply operator ono state 24# compute domain wall initial state 25dw_str = "".join("1" for i in range(L // 2)) + "".join("0" for i in range(L - L // 2)) 26i_0 = basis.index(dw_str) # find index of product state in basis 27psi = np.zeros(basis.Ns) # allocate space for state 28psi[i_0] = 1.0 # set MB state to be the given product state 29# apply H on psi 30psi_new = H.dot(psi) 31# diagonalise operator 32E, V = H.eigsh(k=1, which="SA") 33# calculate operator squared 34H_squared = H.dot(H)
- __init__(static_list, N=None, basis=None, diagonal=None, check_symm=True, check_herm=True, check_pcon=True, dtype=<class 'numpy.complex128'>, copy=False, **basis_args)[source]
Intializes the quantum_LinearOperator object.
- Parameters:
- static_listlist
Contains list of objects to calculate the static part of a quantum_LinearOperator operator. Same as the static argument of the quantum_operator class. The format goes like:
>>> static_list=[[opstr_1,[indx_11,...,indx_1m]],matrix_2,...]
- Nint, optional
number of sites to create the default spin basis with.
- basis
basis
, optional basis object to construct quantum operator with.
- diagonalarray_like
array containing diagonal matrix elements precalculated by other means.
- dtype‘type’
Data type (e.g. numpy.float64) to construct the operator with.
- check_symmbool, optional
Enable/Disable symmetry check on static_list and dynamic_list.
- check_hermbool, optional
Enable/Disable hermiticity check on static_list and dynamic_list.
- check_pconbool, optional
Enable/Disable particle conservation check on static_list and dynamic_list.
- basis_argsdict
Optional additional arguments to pass to the basis class, if not already using a basis object to create the operator.
Methods
__init__
(static_list[, N, basis, diagonal, ...])Intializes the quantum_LinearOperator object.
adjoint
()Hermitian adjoint.
conj
()Conjugates quantum_LinearOperator operator.
Conjugates quantum_LinearOperator operator.
copy
()Returns a deep copy of quantum_LinearOperator object.
dot
(other[, out, a])Matrix-vector multiplication of quantum_LinearOperator operator, with state V.
eigsh
(**eigsh_args)Computes SOME eigenvalues and eigenvectors of hermitian quantum_LinearOperator operator using SPARSE hermitian methods.
expt_value
(V[, enforce_pure])Calculates expectation value of quantum_LinearOperator object in state V.
getH
([copy])Calculates hermitian conjugate of quantum_LinearOperator operator.
matmat
(X)Matrix-matrix multiplication.
matrix_ele
(Vl, Vr[, diagonal])Calculates matrix element of quantum_LinearOperator object between states Vl and Vr.
matvec
(x)Matrix-vector multiplication.
quant_fluct
(V[, enforce_pure, check, time])Calculates the quantum fluctuations (variance) of quantum_LinearOperator object in state V.
rmatmat
(X)Adjoint matrix-matrix multiplication.
rmatvec
(x)Adjoint matrix-vector multiplication.
set_diagonal
(diagonal[, copy])Sets the diagonal part of the quantum_LinearOperator.
transpose
([copy])Transposes quantum_LinearOperator operator.
Attributes
transposes and conjugates the operator matrix, \(H_{ij}\mapsto H_{ji}^*\).
number of states in the (symmetry-reduced) Hilbert space spanned by basis.
transposes the operator matrix, \(H_{ij}\mapsto H_{ji}\).
basis used to build the hamiltonian object.
static diagonal part of the linear operator.
data type of quantum_LinearOperator object.
shape of the quantum_LinearOperator object, always equal to (Ns,Ns).
number of dimensions, always equal to 2.
shape of linear operator.
operator list used to create this object.
- property H
transposes and conjugates the operator matrix, \(H_{ij}\mapsto H_{ji}^*\).
- Type:
- property Ns
number of states in the (symmetry-reduced) Hilbert space spanned by basis.
- Type:
int
- property T
transposes the operator matrix, \(H_{ij}\mapsto H_{ji}\).
- Type:
- adjoint()
Hermitian adjoint.
Returns the Hermitian adjoint of self, aka the Hermitian conjugate or Hermitian transpose. For a complex matrix, the Hermitian adjoint is equal to the conjugate transpose.
Can be abbreviated self.H instead of self.adjoint().
- Returns:
- A_HLinearOperator
Hermitian adjoint of self.
- property basis
basis used to build the hamiltonian object.
Defaults to None if operator has no basis (i.e. was created externally and passed as a precalculated array).
- Type:
- conj()[source]
Conjugates quantum_LinearOperator operator.
- Returns:
quantum_LinearOperator
\(H_{ij}\mapsto H_{ij}^*\)
Notes
This function does NOT transpose the operator.
Examples
>>> H_conj = H.conj()
- conjugate()[source]
Conjugates quantum_LinearOperator operator.
- Returns:
quantum_LinearOperator
\(H_{ij}\mapsto H_{ij}^*\)
Notes
This function does NOT transpose the operator.
Examples
>>> H_conj = H.conj()
- property diagonal
static diagonal part of the linear operator.
- Type:
numpy.ndarray
- dot(other, out=None, a=1.0)[source]
Matrix-vector multiplication of quantum_LinearOperator operator, with state V.
\[aH|V\rangle\]- Parameters:
- othernumpy.ndarray
Vector (quantums tate) to multiply the quantum_LinearOperator operator with.
- outarray_like, optional
specify the output array for the the result.
- ascalar, optional
scalar to multiply the final product with: \(B = aHA\).
- Returns:
- numpy.ndarray
Vector corresponding to the hamiltonian operator applied on the state V.
Examples
>>> B = H.dot(A,check=True)
corresponds to \(B = HA\).
- property dtype
data type of quantum_LinearOperator object.
- Type:
type
- eigsh(**eigsh_args)[source]
Computes SOME eigenvalues and eigenvectors of hermitian quantum_LinearOperator operator using SPARSE hermitian methods.
This function method solves for eigenvalues and eigenvectors, but can only solve for a few of them accurately. It calls scipy.sparse.linalg.eigsh, which is a wrapper for ARPACK.
- Parameters:
- eigsh_args
For all additional arguments see documentation of scipy.sparse.linalg.eigsh.
- Returns:
- tuple
Tuple containing the (eigenvalues, eigenvectors) of the quantum_LinearOperator operator.
Notes
Assumes the operator is hermitian! If the flat check_hermiticity = False is used, we advise the user to reassure themselves of the hermiticity properties before use.
Examples
>>> eigenvalues,eigenvectors = H.eigsh(**eigsh_args)
- expt_value(V, enforce_pure=False)[source]
Calculates expectation value of quantum_LinearOperator object in state V.
\[\langle V|H|V\rangle\]- Parameters:
- Vnumpy.ndarray
Depending on the shape, can be a single state or a collection of pure or mixed states [see enforce_pure argument of basis.ent_entropy].
- enforce_purebool, optional
Flag to enforce pure expectation value of V is a square matrix with multiple pure states in the columns.
- Returns:
- float
Expectation value of hamiltonian operator in state V.
Examples
>>> H_expt = H.expt_value(V,time=0,diagonal=False,check=True)
corresponds to \(H_{expt} = \langle V|H(t=0)|V\rangle\).
- getH(copy=False)[source]
Calculates hermitian conjugate of quantum_LinearOperator operator.
- Parameters:
- copybool, optional
Whether to return a deep copy of the original object. Default is copy = False.
- Returns:
quantum_LinearOperator
\(H_{ij}\mapsto H_{ji}^*\)
Examples
>>> H_herm = H.getH()
- property get_shape
shape of the quantum_LinearOperator object, always equal to (Ns,Ns).
- Type:
tuple
- matmat(X)
Matrix-matrix multiplication.
Performs the operation y=A*X where A is an MxN linear operator and X dense N*K matrix or ndarray.
- Parameters:
- X{matrix, ndarray}
An array with shape (N,K).
- Returns:
- Y{matrix, ndarray}
A matrix or ndarray with shape (M,K) depending on the type of the X argument.
Notes
This matmat wraps any user-specified matmat routine or overridden _matmat method to ensure that y has the correct type.
- matrix_ele(Vl, Vr, diagonal=False)[source]
Calculates matrix element of quantum_LinearOperator object between states Vl and Vr.
\[\langle V_l|H|V_r\rangle\]- Parameters:
- Vlnumpy.ndarray
Vector(s)/state(s) to multiple with on left side.
- Vlnumpy.ndarray
Vector(s)/state(s) to multiple with on right side.
- diagonalbool, optional
When set to True, returs only diagonal part of expectation value. Default is diagonal = False.
- Returns:
- float
Matrix element of quantum_operator quantum_operator between the states Vl and Vr.
Notes
Taking the conjugate or transpose of the state Vl is done automatically.
Examples
>>> H_lr = H.expt_value(Vl,Vr,pars=pars,diagonal=False,check=True)
corresponds to \(H_{lr} = \langle V_l|H(\lambda=0)|V_r\rangle\).
- matvec(x)
Matrix-vector multiplication.
Performs the operation y=A*x where A is an MxN linear operator and x is a column vector or 1-d array.
- Parameters:
- x{matrix, ndarray}
An array with shape (N,) or (N,1).
- Returns:
- y{matrix, ndarray}
A matrix or ndarray with shape (M,) or (M,1) depending on the type and shape of the x argument.
Notes
This matvec wraps the user-specified matvec routine or overridden _matvec method to ensure that y has the correct shape and type.
- property ndim
number of dimensions, always equal to 2.
- Type:
int
- quant_fluct(V, enforce_pure=False, check=True, time=0)[source]
Calculates the quantum fluctuations (variance) of quantum_LinearOperator object in state V.
\[\langle V|H^2|V\rangle - \langle V|H|V\rangle^2\]- Parameters:
- Vnumpy.ndarray
Depending on the shape, can be a single state or a collection of pure or mixed states [see enforce_pure].
- enforce_purebool, optional
Flag to enforce pure expectation value of V is a square matrix with multiple pure states in the columns.
- Returns:
- float
Quantum fluctuations of hamiltonian operator in state V.
Examples
>>> H_fluct = H.quant_fluct(V,time=0,check=True)
corresponds to \(\left(\Delta H\right)^2 = \langle V|H^2(t=\texttt{time})|V\rangle - \langle V|H(t=\texttt{time})|V\rangle^2\).
- rmatmat(X)
Adjoint matrix-matrix multiplication.
Performs the operation y = A^H * x where A is an MxN linear operator and x is a column vector or 1-d array, or 2-d array. The default implementation defers to the adjoint.
- Parameters:
- X{matrix, ndarray}
A matrix or 2D array.
- Returns:
- Y{matrix, ndarray}
A matrix or 2D array depending on the type of the input.
Notes
This rmatmat wraps the user-specified rmatmat routine.
- rmatvec(x)
Adjoint matrix-vector multiplication.
Performs the operation y = A^H * x where A is an MxN linear operator and x is a column vector or 1-d array.
- Parameters:
- x{matrix, ndarray}
An array with shape (M,) or (M,1).
- Returns:
- y{matrix, ndarray}
A matrix or ndarray with shape (N,) or (N,1) depending on the type and shape of the x argument.
Notes
This rmatvec wraps the user-specified rmatvec routine or overridden _rmatvec method to ensure that y has the correct shape and type.
- set_diagonal(diagonal, copy=True)[source]
Sets the diagonal part of the quantum_LinearOperator.
- Parameters:
- diagonal: array_like
array_like object containing the new diagonal.
- property shape
shape of linear operator.
- Type:
tuple
- property static_list
operator list used to create this object.
- Type:
list
- transpose(copy=False)[source]
Transposes quantum_LinearOperator operator.
- Returns:
quantum_LinearOperator
\(H_{ij}\mapsto H_{ji}\)
Notes
This function does NOT conjugate the operator.
Examples
>>> H_tran = H.transpose()