quspin.operators.quantum_operator
- class quspin.operators.quantum_operator(input_dict, N=None, basis=None, shape=None, copy=True, check_symm=True, check_herm=True, check_pcon=True, matrix_formats={}, dtype=<class 'numpy.complex128'>, **basis_args)[source]
Bases:
object
Constructs parameter-dependent (hermitian and nonhermitian) operators.
The quantum_operator class maps quantum operators to keys of a dictionary. When calling various methods of quantum_operator, it allows one to ‘dynamically’ specify the pre-factors of these operators.
Examples
It is often required to be able to handle a parameter-dependent Hamiltonian \(H(\lambda)=H_0 + \lambda H_1\), e.g.
\[H_0=\sum_j J_{zz}S^z_jS^z_{j+2} + h_xS^x_j, \qquad H_1=\sum_j S^z_j\]The following code snippet shows how to use the quantum_operator class to vary the parameter \(\lambda\) without having to re-build the Hamiltonian every time.
1from quspin.basis import spin_basis_1d # Hilbert spaces 2import numpy as np # general math functions 3from numpy.random import uniform 4 5# 6##### setting parameters for simulation 7# physical parameters 8L = 4 # system size 9J = 1.0 # interaction strength 10hx = np.sqrt(2) # transverse field strength 11# 12###### create the basis 13basis = spin_basis_1d(L, pblock=1, pauli=False) # up basis 14##### create model 15# define static (fixed) site-coupling lists 16J_list = [[J, i, (i + 2) % L] for i in range(L)] # nnn interaction PBC 17hx_list = [[hx, i] for i in range(L)] # onsite field 18# create static lists for H0 19operator_list_0 = [["zz", J_list], ["x", hx_list]] 20# define parametric lists for H1 (corresponding to operators the coupling of which will be changed) 21hz_list = [[1.0, i] for i in range(L)] # onsite field 22operator_list_1 = [["z", hz_list]] 23# 24###### create operator dictionary for quantum_operator class 25# add keys for TFI operators tring list and the longitudinal field 26operator_dict = dict(H0=operator_list_0, H1=operator_list_1) 27# 28###### setting up `quantum_operator` hamiltonian dictionary 29H = quantum_operator(operator_dict, basis=basis) 30# print Hamiltonian H = H0 + H1 31params_dict = dict( 32 H0=1.0, H1=1.0 33) # dict containing the couplings to operators in keys of operator_dict 34H_lmbda1 = H.tohamiltonian(params_dict) 35print(H_lmbda1) 36# change z-coupling strength: print Hamiltonian H = H0 + 2H1 37params_dict = dict( 38 H0=1.0, H1=2.0 39) # dict containing the couplings to operators in keys of operator_dict 40H_lmbda2 = H.tohamiltonian(params_dict) 41print(H_lmbda2)
- __init__(input_dict, N=None, basis=None, shape=None, copy=True, check_symm=True, check_herm=True, check_pcon=True, matrix_formats={}, dtype=<class 'numpy.complex128'>, **basis_args)[source]
Intializes the quantum_operator object (parameter dependent quantum quantum_operators).
- Parameters:
- input_dictdict
The values of this dictionary contain quantum_operator lists, in the same format as the static_list argument of the hamiltonian class.
The keys of this dictionary correspond to the parameter values, e.g. \(J_{zz},h_x\), and are used to specify the coupling strength during calls of the quantum_operator class methods.
>>> # use "Jzz" and "hx" keys to specify the zz and x coupling strengths, respectively >>> input_dict = { "Jzz": [["zz",Jzz_bonds]], "hx" : [["x" ,hx_site ]] }
- Nint, optional
Number of lattice sites for the hamiltonian object.
- dtype‘type’
Data type (e.g. numpy.float64) to construct the quantum_operator with.
- shapetuple, optional
Shape to create the hamiltonian object with. Default is shape = None.
- copy: bool, optional
If set to True, this option creates a copy of the input array.
- 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.
- matrix_formats: dict, optional
Dictionary of key,value pairs which, given a key associated with an operator in input_dict, the value of this key specifies the sparse matrix format {“csr”,”csc”,”dia”,”dense”}.
- kw_argsdict
Optional additional arguments to pass to the basis class, if not already using a basis object to create the quantum_operator.
Methods
__init__
(input_dict[, N, basis, shape, ...])Intializes the quantum_operator object (parameter dependent quantum quantum_operators).
aslinearoperator
([pars])Returns copy of a quantum_operator object for parameters pars as a scipy.sparse.linalg.Linearquantum_operator.
astype
(dtype[, copy, casting])Changes data type of quantum_operator object.
conj
()Conjugates quantum_operator quantum_operator.
Conjugates quantum_operator quantum_operator.
copy
()Returns a deep copy of quantum_operator object.
diagonal
([pars])Returns diagonal of quantum_operator quantum_operator for parameters pars.
dot
(V[, pars, check, out, overwrite_out, a])Matrix-vector multiplication of quantum_operator quantum_operator for parameters pars, with state V.
eigh
([pars])Computes COMPLETE eigensystem of hermitian quantum_operator quantum_operator using DENSE hermitian methods.
eigsh
([pars])Computes SOME eigenvalues and eigenvectors of hermitian quantum_operator quantum_operator using SPARSE hermitian methods.
eigvalsh
([pars])Computes ALL eigenvalues of hermitian quantum_operator quantum_operator using DENSE hermitian methods.
expt_value
(V[, pars, check, enforce_pure])Calculates expectation value of of quantum_operator object for parameters pars, in state V.
getH
([copy])Calculates hermitian conjugate of quantum_operator quantum_operator.
get_operators
(key)matmat
(X)Matrix-matrix multiplication.
matrix_ele
(Vl, Vr[, pars, diagonal, check])Calculates matrix element of quantum_operator object for parameters pars in states Vl and Vr.
matvec
(x)Matrix-vector multiplication.
quant_fluct
(V[, pars, check, enforce_pure])Calculates the quantum fluctuations (variance) of quantum_operator object for parameters pars, in state V.
rdot
(V[, pars, check, out, overwrite_out, a])Vector-matrix multiplication of quantum_operator quantum_operator for parameters pars, with state V.
rmatvec
(x)Adjoint matrix-vector multiplication.
toarray
([pars, out])Returns copy of a quantum_operator object for parameters pars as a dense array.
tocsc
([pars])Returns copy of a quantum_operator object for parameters pars as a scipy.sparse.csc_matrix.
tocsr
([pars])Returns copy of a quantum_operator object for parameters pars as a scipy.sparse.csr_matrix.
todense
([pars, out])Returns copy of a quantum_operator object for parameters pars as a dense array.
tohamiltonian
([pars, copy])Returns copy of a quantum_operator object for parameters pars as a hamiltonian object.
trace
([pars])Calculates trace of quantum_operator quantum_operator for parameters pars.
transpose
([copy])Transposes quantum_operator quantum_operator.
update_matrix_formats
(matrix_formats)Change the internal structure of the matrices in-place.
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.
data type of quantum_operator object.
shape of the quantum_operator object, always equal to (Ns,Ns).
True if quantum_operator contains a dense matrix as a componnent of either the static or dynamic list.
number of dimensions, always equal to 2.
shape of the quantum_operator object, always equal to (Ns,Ns).
- 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:
- aslinearoperator(pars={})[source]
Returns copy of a quantum_operator object for parameters pars as a scipy.sparse.linalg.Linearquantum_operator.
Casts the quantum_operator object as a scipy.sparse.linalg.Linearquantum_operator object.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- Returns:
scipy.sparse.linalg.Linearquantum_operator
Examples
>>> H_aslinop=H.aslinearquantum_operator(pars=pars)
- astype(dtype, copy=False, casting='unsafe')[source]
Changes data type of quantum_operator object.
- Parameters:
- dtype‘type’
The data type (e.g. numpy.float64) to cast the Hamiltonian with.
- Returns
- `quantum_operator`
quantum_operator with altered data type.
Examples
>>> H_cpx=H.astype(np.complex128)
- property basis
basis used to build the hamiltonian object. Defaults to None if quantum_operator has no basis (i.e. was created externally and passed as a precalculated array).
- Type:
- conj()[source]
Conjugates quantum_operator quantum_operator.
- Returns:
quantum_operator
\(H_{ij}\mapsto H_{ij}^*\)
Notes
This function does NOT transpose the quantum_operator.
Examples
>>> H_conj = H.conj()
- conjugate()[source]
Conjugates quantum_operator quantum_operator.
- Returns:
quantum_operator
\(H_{ij}\mapsto H_{ij}^*\)
Notes
This function does NOT transpose the quantum_operator.
Examples
>>> H_conj = H.conj()
- diagonal(pars={})[source]
Returns diagonal of quantum_operator quantum_operator for parameters pars.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- Returns:
- numpy.ndarray
array containing the diagonal part of the operator \(diag_j = H_{jj}(\lambda)\).
Examples
>>> H_diagonal = H.diagonal(pars=pars)
- dot(V, pars={}, check=True, out=None, overwrite_out=True, a=1.0)[source]
Matrix-vector multiplication of quantum_operator quantum_operator for parameters pars, with state V.
\[aH(\lambda)|V\rangle\]- Parameters:
- Vnumpy.ndarray
Vector (quantums tate) to multiply the quantum_operator quantum_operator with.
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- checkbool, optional
Whether or not to do checks for shape compatibility.
- outarray_like, optional
specify the output array for the the result. This is not supported if V is a sparse matrix.
- overwrite_outbool, optional
flag used to toggle between two different ways to treat out. If set to True all values in out will be overwritten with the result of the dot product. If False the result of the dot product will be added to the values of out.
- ascalar, optional
scalar to multiply the final product with: \(B = aHV\).
- Returns:
- numpy.ndarray
Vector corresponding to the quantum_operator quantum_operator applied on the state V.
Examples
>>> B = H.dot(A,pars=pars,check=True)
corresponds to \(B = HA\).
- property dtype
data type of quantum_operator object.
- Type:
type
- eigh(pars={}, **eigh_args)[source]
Computes COMPLETE eigensystem of hermitian quantum_operator quantum_operator using DENSE hermitian methods.
This function method solves for all eigenvalues and eigenvectors. It calls numpy.linalg.eigh, and uses wrapped LAPACK functions which are contained in the module py_lapack.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- eigh_args
For all additional arguments see documentation of numpy.linalg.eigh.
- Returns:
- tuple
Tuple containing the (eigenvalues, eigenvectors) of the quantum_operator quantum_operator.
Notes
Assumes the quantum_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.eigh(pars=pars,**eigh_args)
- eigsh(pars={}, **eigsh_args)[source]
Computes SOME eigenvalues and eigenvectors of hermitian quantum_operator quantum_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:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- eigsh_args
For all additional arguments see documentation of scipy.sparse.linalg.eigsh.
- Returns:
- tuple
Tuple containing the (eigenvalues, eigenvectors) of the quantum_operator quantum_operator.
Notes
Assumes the quantum_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(pars=pars,**eigsh_args)
- eigvalsh(pars={}, **eigvalsh_args)[source]
Computes ALL eigenvalues of hermitian quantum_operator quantum_operator using DENSE hermitian methods.
This function method solves for all eigenvalues. It calls numpy.linalg.eigvalsh, and uses wrapped LAPACK functions which are contained in the module py_lapack.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- eigvalsh_args
For all additional arguments see documentation of numpy.linalg.eigvalsh.
- Returns:
- numpy.ndarray
Eigenvalues of the quantum_operator quantum_operator.
Notes
Assumes the quantum_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 = H.eigvalsh(pars=pars,**eigvalsh_args)
- expt_value(V, pars={}, check=True, enforce_pure=False)[source]
Calculates expectation value of of quantum_operator object for parameters pars, in state V.
\[\langle V|H(\lambda)|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].
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- enforce_purebool, optional
Flag to enforce pure expectation value of V is a square matrix with multiple pure states in the columns.
- checkbool, optional
- 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_operator quantum_operator.
- Parameters:
- copybool, optional
Whether to return a deep copy of the original object. Default is copy = False.
- Returns:
quantum_operator
\(H_{ij}\mapsto H_{ij}^*\)
Examples
>>> H_herm = H.getH()
- property get_shape
shape of the quantum_operator object, always equal to (Ns,Ns).
- Type:
tuple
- property is_dense
True if quantum_operator contains a dense matrix as a componnent of either the static or dynamic list.
- Type:
bool
- matmat(X)[source]
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, pars={}, diagonal=False, check=True)[source]
Calculates matrix element of quantum_operator object for parameters pars in states Vl and Vr.
\[\langle V_l|H(\lambda)|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.
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- diagonalbool, optional
When set to True, returs only diagonal part of expectation value. Default is diagonal = False.
- checkbool,
- 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)[source]
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, pars={}, check=True, enforce_pure=False)[source]
Calculates the quantum fluctuations (variance) of quantum_operator object for parameters pars, in state V.
\[\langle V|H(\lambda)^2|V\rangle - \langle V|H(\lambda)|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].
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- enforce_purebool, optional
Flag to enforce pure expectation value of V is a square matrix with multiple pure states in the columns.
- checkbool, optional
- Returns:
- float
Quantum fluctuations of hamiltonian operator in state V.
Examples
>>> H_fluct = H.quant_fluct(V,time=0,diagonal=False,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\).
- rdot(V, pars={}, check=False, out=None, overwrite_out=True, a=1.0)[source]
Vector-matrix multiplication of quantum_operator quantum_operator for parameters pars, with state V.
\[a\langle V]H(\lambda)\]- Parameters:
- Vnumpy.ndarray
Vector (quantums tate) to multiply the quantum_operator quantum_operator with.
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- checkbool, optional
Whether or not to do checks for shape compatibility.
- outarray_like, optional
specify the output array for the the result. This is not supported if V is a sparse matrix.
- overwrite_outbool, optional
flag used to toggle between two different ways to treat out. If set to True all values in out will be overwritten with the result. If False the result of the dot product will be added to the values of out.
- ascalar, optional
scalar to multiply the final product with: \(B = aVH\).
- Returns:
- numpy.ndarray
Vector corresponding to the quantum_operator quantum_operator applied on the state V.
Examples
>>> B = H.dot(A,pars=pars,check=True)
corresponds to \(B = AH\).
- rmatvec(x)[source]
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.
- property shape
shape of the quantum_operator object, always equal to (Ns,Ns).
- Type:
tuple
- toarray(pars={}, out=None)[source]
Returns copy of a quantum_operator object for parameters pars as a dense array.
This function can overflow memory if not used carefully!
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- outnumpy.ndarray
Array to fill in with the output.
- Returns:
- numpy.ndarray
Dense array.
Examples
>>> H_dense=H.toarray(pars=pars)
- tocsc(pars={})[source]
Returns copy of a quantum_operator object for parameters pars as a scipy.sparse.csc_matrix.
Casts the quantum_operator object as a scipy.sparse.csc_matrix object.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- Returns:
scipy.sparse.csc_matrix
Examples
>>> H_csc=H.tocsc(pars=pars)
- tocsr(pars={})[source]
Returns copy of a quantum_operator object for parameters pars as a scipy.sparse.csr_matrix.
Casts the quantum_operator object as a scipy.sparse.csr_matrix object.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- Returns:
scipy.sparse.csr_matrix
Examples
>>> H_csr=H.tocsr(pars=pars)
- todense(pars={}, out=None)[source]
Returns copy of a quantum_operator object for parameters pars as a dense array.
This function can overflow memory if not used carefully!
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- outnumpy.ndarray
Array to fill in with the output.
- Returns:
- obj
Depending of size of array, can be either one of
numpy.ndarray.
numpy.matrix.
Notes
If the array dimension is too large, scipy may choose to cast the quantum_operator quantum_operator as a numpy.matrix instead of a numpy.ndarray. In such a case, one can use the quantum_operator.toarray() method.
Examples
>>> H_dense=H.todense(pars=pars)
- tohamiltonian(pars={}, copy=True)[source]
Returns copy of a quantum_operator object for parameters pars as a hamiltonian object.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- copybool, optional
Explicitly copy matrices when constructing the hamiltonian object, default is True.
- Returns:
Examples
>>> H_ham=H.tohamiltonian(pars=pars)
- trace(pars={})[source]
Calculates trace of quantum_operator quantum_operator for parameters pars.
- Parameters:
- parsdict, optional
Dictionary with same keys as input_dict and coupling strengths as values. Any missing keys are assumed to be set to unity.
- Returns:
- float
Trace of quantum_operator \(\sum_{j=1}^{Ns} H_{jj}(\lambda)\).
Examples
>>> H_tr = H.trace(pars=pars)
- transpose(copy=False)[source]
Transposes quantum_operator quantum_operator.
- Returns:
quantum_operator
\(H_{ij}\mapsto H_{ji}\)
Notes
This function does NOT conjugate the quantum_operator.
Examples
>>> H_tran = H.transpose()
- update_matrix_formats(matrix_formats)[source]
Change the internal structure of the matrices in-place.
- Parameters:
- matrix_formats: dict, optional
Dictionary of key,value pairs which, given a key associated with an operator in input_dict, the value of this key specifies the sparse matrix format {“csr”,”csc”,”dia”,”dense”}.
Examples
Given O which has two operators defined by strings: ‘Hx’ for transverse field, and ‘Hising’ is the Ising part. The Ising part must be diagonal therefore it is useful to cast it to a DIA matrix format, while the transverse field is not diagonal so it is most efficient to use CSR matrix format. This can be accomplished by the following:
>>> O.update_matrix_formats(dict(Hx="csr",Hising="dia"))