quspin.operators.hamiltonian
- class quspin.operators.hamiltonian(static_list, dynamic_list, N=None, basis=None, shape=None, dtype=<class 'numpy.complex128'>, static_fmt=None, dynamic_fmt=None, copy=True, check_symm=True, check_herm=True, check_pcon=True, **basis_kwargs)[source]
Bases:
object
Constructs time-dependent (hermitian and nonhermitian) operators.
The hamiltonian class wraps most of the functionalty of the QuSpin package. This object allows the user to construct lattice Hamiltonians and operators, solve the time-dependent Schroedinger equation, do full/Lanczos diagonalization, etc.
The user can create both static and time-dependent, hermitian and non-hermitian operators for any particle type (boson, spin, fermion) specified by the basis constructor.
Notes
One can instantiate the class either by parsing a set of symmetries, or an instance of basis. Note that instantiation with a basis will automatically ignore all symmetry inputs.
Examples
Here is an example how to employ a basis object to construct the periodically driven XXZ Hamiltonian
\[H(t) = \sum_{j=0}^{L-1} \left( JS^z_{j+1}S^z_j + hS^z_j + g\cos(\Omega t)S^x_j \right)\]in the zero-momentum sector (kblock=0) of positive parity (pblock=1). We use periodic boundary conditions.
The code snippet below initiates the class, and is required to run the example codes for the function methods.
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 # longitudinal field 10##### define periodic drive ##### 11Omega = 4.5 # drive frequency 12 13 14def drive(t, Omega): 15 return np.cos(Omega * t) 16 17 18drive_args = [Omega] 19# 20##### construct basis in the 0-total momentum and +1-parity sector 21basis = spin_basis_1d(L=L, a=1, kblock=0, pblock=1) 22# define PBC site-coupling lists for operators 23x_field = [[g, i] for i in range(L)] 24z_field = [[h, i] for i in range(L)] 25J_nn = [[J, i, (i + 1) % L] for i in range(L)] # PBC 26# static and dynamic lists 27static = [["zz", J_nn], ["z", z_field]] 28dynamic = [["x", x_field, drive, drive_args]] 29###### construct Hamiltonian 30H = hamiltonian(static, dynamic, static_fmt="dia", dtype=np.float64, basis=basis) 31print(H.toarray())
- __init__(static_list, dynamic_list, N=None, basis=None, shape=None, dtype=<class 'numpy.complex128'>, static_fmt=None, dynamic_fmt=None, copy=True, check_symm=True, check_herm=True, check_pcon=True, **basis_kwargs)[source]
Intializes the hamtilonian object (any quantum operator).
- Parameters:
- static_listlist
Contains list of objects to calculate the static part of a hamiltonian operator. The format goes like:
>>> static_list=[[opstr_1,[indx_11,...,indx_1m]],matrix_2,...]
- dynamic_listlist
Contains list of objects to calculate the dynamic (time-dependent) part of a hamiltonian operator. The format goes like:
>>> dynamic_list=[[opstr_1,[indx_11,...,indx_1n],fun_1,fun_1_args],[matrix_2,fun_2,fun_2_args],...]
fun: function object which multiplies the matrix or operator given in the same list.
fun_args: tuple of the extra arguments which go into the function to evaluate it like:
>>> f_val = fun(t,*fun_args)
If the operator is time-INdependent, one must pass an empty list: dynamic_list = [].
- Nint, optional
Number of lattice sites for the hamiltonian object.
- dtypenumpy.datatype, optional
Data type (e.g. numpy.float64) to construct the operator with.
- static_fmtstr {“csr”,”csc”,”dia”,”dense”}, optional
Specifies format of static part of Hamiltonian.
- dynamic_fmt: str {“csr”,”csc”,”dia”,”dense”} or dict, keys: (func,func_args), values: str {“csr”,”csc”,”dia”,”dense”}
Specifies the format of the dynamic parts of the hamiltonian. To specify a particular dynamic part of the hamiltonian use a tuple (func,func_args) which matches a function+argument pair used in the construction of the hamiltonian as a key in the dictionary.
- 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.
- basis_kwargsdict
Optional additional arguments to pass to the basis class, if not already using a basis object to create the operator.
Methods
__init__
(static_list, dynamic_list[, N, ...])Intializes the hamtilonian object (any quantum operator).
as_dense_format
([copy])Casts hamiltonian operator to DENSE format.
as_sparse_format
([static_fmt, dynamic_fmt, copy])Casts hamiltonian operator to SPARSE format(s).
aslinearoperator
([time])Returns copy of a hamiltonian object at time time as a scipy.sparse.linalg.LinearOperator.
astype
(dtype[, copy, casting])Changes data type of hamiltonian object.
updates attribute _.is_dense.
conj
()Same functionality as
conjugate()
.Conjugates hamiltonian operator.
copy
()Returns a copy of hamiltonian object.
diagonal
([time])Calculates diagonal of hamiltonian operator at time time.
dot
(V[, time, check, out, overwrite_out, a])Matrix-vector multiplication of hamiltonian operator at time time, with state V.
eigh
([time])Computes COMPLETE eigensystem of hermitian hamiltonian operator using DENSE hermitian methods.
eigsh
([time])Computes SOME eigenvalues and eigenvectors of hermitian hamiltonian operator using SPARSE hermitian methods.
eigvalsh
([time])Computes ALL eigenvalues of hermitian hamiltonian operator using DENSE hermitian methods.
evolve
(v0, t0, times[, eom, solver_name, ...])Implements (imaginary) time evolution generated by the hamiltonian object.
expt_value
(V[, time, check, enforce_pure])Calculates expectation value of hamiltonian operator at time time, in state V.
getH
([copy])Calculates hermitian conjugate of hamiltonian operator.
matrix_ele
(Vl, Vr[, time, diagonal, check])Calculates matrix element of hamiltonian operator at time time in states Vl and Vr.
project_to
(proj)Projects/Transforms hamiltonian operator with projector/operator proj.
quant_fluct
(V[, time, check, enforce_pure])Calculates the quantum fluctuations (variance) of hamiltonian operator at time time, in state V.
rdot
(V[, time, check, out, overwrite_out, a])Vector-Matrix multiplication of hamiltonian operator at time time, with state V.
rotate_by
(other[, generator])Rotates/Transforms hamiltonian object by an operator other.
toarray
([time, order, out])Returns copy of a hamiltonian object at time time as a dense array.
tocsc
([time])Returns copy of a hamiltonian object at time time as a scipy.sparse.csc_matrix.
tocsr
([time])Returns copy of a hamiltonian object at time time as a scipy.sparse.csr_matrix.
todense
([time, order, out])Returns copy of a hamiltonian object at time time as a dense array.
trace
([time])Calculates trace of hamiltonian operator at time time.
transpose
([copy])Transposes hamiltonian operator.
update_matrix_formats
(static_fmt, dynamic_fmt)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 hamiltonian object.
contains dynamic parts of the operator as dict(func=Hdyn).
shape of the hamiltonian object, always equal to (Ns,Ns).
checks sparsity of operator matrix.
Total bytes consumed by the elements of the hamiltonian array.
number of dimensions, always equal to 2.
shape of the hamiltonian object, always equal to (Ns,Ns).
static part of the operator.
- 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:
- as_dense_format(copy=False)[source]
Casts hamiltonian operator to DENSE format.
- Parameters:
- copybool,optional
Whether to return a deep copy of the original object. Default is copy = False.
- Returns:
- obj
Either one of the following:
Shallow copy, if copy = False.
Deep copy, if copy = True.
Examples
>>> H_dense=H.as_dense_format()
- as_sparse_format(static_fmt='csr', dynamic_fmt={}, copy=False)[source]
Casts hamiltonian operator to SPARSE format(s).
- Parameters:
- static_fmtstr {“csr”,”csc”,”dia”,”dense”}
Specifies format of static part of Hamiltonian.
- dynamic_fmt: str {“csr”,”csc”,”dia”,”dense”} or dict, keys: (func,func_args), values: str {“csr”,”csc”,”dia”,”dense”}
Specifies the format of the dynamic parts of the hamiltonian. To specify a particular dynamic part of the hamiltonian use a tuple (func,func_args) which matches a function+argument pair used in the construction of the hamiltonian as a key in the dictionary.
- copybool,optional
Whether to return a deep copy of the original object. Default is copy = False.
- Returns:
- obj
Either one of the following:
whenever possible do not copy data, if copy = False.
explicitly copy all possible data, if copy = True.
Examples
>>> H_dia=H.as_sparse_format(static_fmt="csr",dynamic_fmt={(func,func_args):"dia"})
- aslinearoperator(time=0.0)[source]
Returns copy of a hamiltonian object at time time as a scipy.sparse.linalg.LinearOperator.
Casts the hamiltonian object as a scipy.sparse.linalg.LinearOperator object.
- Parameters:
- timefloat, optional
Time to evalute the time-dependent part of the operator at (if existent). Default is time = 0.0.
- Returns:
scipy.sparse.linalg.LinearOperator
Examples
>>> H_aslinop=H.aslinearoperator(time=time)
- astype(dtype, copy=False, casting='unsafe')[source]
Changes data type of hamiltonian object.
- Parameters:
- dtype‘type’
The data type (e.g. numpy.float64) to cast the Hamiltonian with.
- Returns:
- hamiltonian
Operator with altered data type.
Examples
- hamiltonian
Operator with altered data type.
>>> H_cpx=H.astype(np.complex128)
- 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]
Same functionality as
conjugate()
.
- conjugate()[source]
Conjugates hamiltonian operator.
- Returns:
hamiltonian
\(H_{ij}\mapsto H_{ij}^*\)
Notes
This function does NOT transpose the operator.
Examples
>>> H_conj = H.conjugate()
- diagonal(time=0)[source]
Calculates diagonal of hamiltonian operator at time time.
- Parameters:
- timefloat, optional
Time to evalute the time-dependent part of the operator at (if existent). Default is time = 0.0.
- Returns:
- numpy.ndarray
Diagonal part of operator \(H(t=\texttt{time})\).
Examples
>>> H_diag = H.diagonal(time=0.0)
- dot(V, time=0, check=True, out=None, overwrite_out=True, a=1.0)[source]
Matrix-vector multiplication of hamiltonian operator at time time, with state V.
\[aH(t=\texttt{time})|V\rangle\]- Parameters:
- V{numpy.ndarray, scipy.spmatrix}
Array containing the quantums state to multiply the hamiltonian operator with.
- timeobj, optional
Can be either one of the following:
- float: time to evalute the time-dependent part of the operator at (if operator has time dependence).
Default is time = 0.
- (N,) array_like: if V.shape[-1] == N, the hamiltonian operator is evaluated at the i-th time
and dotted into V[…,i] to get the i-th slice of the output array. Here V must be either 2- or 3-d array, where 2-d would be for pure states and 3-d would be for mixed states.
- 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 or if times is an array.
- 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 = aHV\).
- Returns:
- numpy.ndarray
Vector corresponding to the hamiltonian operator applied on the state V.
Notes
this function does the matrix multiplication with the state(s) and Hamiltonian as is, see Example 17 (Lidblad dynamics / Optical Bloch Equations)
for right-multiplication of quantum operators, see function rdot().
Examples
>>> B = H.dot(A,time=0,check=True)
corresponds to \(B = HA\).
- property dtype
data type of hamiltonian object.
- Type:
type
- property dynamic
contains dynamic parts of the operator as dict(func=Hdyn).
The key func is the memory address of the time-dependent function which can be called as func(time). The function arguments are hard-coded, and are not passed. The value Hdyn is the sparse matrix to which the drive couples.
- Type:
dict
- eigh(time=0, **eigh_args)[source]
Computes COMPLETE eigensystem of hermitian hamiltonian 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:
- timefloat
Time to evalute the hamiltonian operator at (if time dependent). Default is time = 0.0.
- eigh_args
For all additional arguments see documentation of numpy.linalg.eigh.
- Returns:
- tuple
Tuple containing the (eigenvalues, eigenvectors) of the hamiltonian 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.eigh(time=time,**eigh_args)
- eigsh(time=0.0, **eigsh_args)[source]
Computes SOME eigenvalues and eigenvectors of hermitian hamiltonian 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:
- timefloat
Time to evalute the hamiltonian operator at (if time dependent). Default is time = 0.0.
- eigsh_args
For all additional arguments see documentation of scipy.sparse.linalg.eigsh.
- Returns:
- tuple
Tuple containing the (eigenvalues, eigenvectors) of the hamiltonian 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(time=time,**eigsh_args)
- eigvalsh(time=0, **eigvalsh_args)[source]
Computes ALL eigenvalues of hermitian hamiltonian 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:
- timefloat
Time to evalute the hamiltonian operator at (if time dependent). Default is time = 0.0.
- eigvalsh_args
For all additional arguments see documentation of numpy.linalg.eigvalsh.
- Returns:
- numpy.ndarray
Eigenvalues of the hamiltonian 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 = H.eigvalsh(time=time,**eigvalsh_args)
- evolve(v0, t0, times, eom='SE', solver_name='dop853', stack_state=False, verbose=False, iterate=False, imag_time=False, **solver_args)[source]
Implements (imaginary) time evolution generated by the hamiltonian object.
The functions handles evolution generated by both time-dependent and time-independent Hamiltonians.
- Currently the following three built-in routines are supported (see parameter eom):
real-time Schroedinger equation: \(\partial_t|v(t)\rangle=-iH(t)|v(t)\rangle\).
imaginary-time Schroedinger equation: \(\partial_t|v(t)\rangle=-H(t)|v(t)\rangle\).
Liouvillian dynamics: \(\partial_t\rho(t)=-i[H,\rho(t)]\).
- Parameters:
- v0numpy.ndarray
Initial state \(|v(t)\rangle\) or density matrix (pure and mixed) \(\rho(t)\).
- t0float
Initial time.
- timesnumpy.ndarray
Vector of times to compute the time-evolved state at.
- eomstr, optional
Specifies the ODE type. Can be either one of
“SE”, real and imaginary-time Schroedinger equation.
“LvNE”, real-time Liouville equation.
Default is “eom = SE” (Schroedinger evolution).
- iteratebool, optional
If set to True, creates a generator object for the time-evolved the state. Default is False.
- solver_namestr, optional
Scipy solver integrator name. Default is dop853.
See scipy integrator (solver) for other options.
- solver_argsdict, optional
Dictionary with additional scipy integrator (solver).
- stack_statebool, optional
Flag to determine if f is real or complex-valued. Default is False (i.e. complex-valued).
- imag_timebool, optional
Must be set to True when f defines imaginary-time evolution, in order to normalise the state at each time in times. Default is False.
- verbosebool, optional
If set to True, prints normalisation of state at teach time in times.
- Returns:
- obj
Can be either one of the following:
numpy.ndarray containing evolved state against time.
generator object for time-evolved state (requires iterate = True).
Note that for Liouvillian dynamics the output is a square complex numpy.ndarray.
Notes
Supports evolution of multiple states simulataneously (eom=”SE”) and evolution of mixed and pure density matrices (`eom=”LvNE”). For a user-defined custom ODE solver which can handle non-linear equations, check out the `measurements.evolve() routine, which has a similar functionality but allows for a complete freedom over the differential equation to be solved.
Examples
>>> v_t = H.evolve(v0,t0,times,eom="SE",solver_name="dop853",verbose=False,iterate=False,imag_time=False,**solver_args)
- expt_value(V, time=0, check=True, enforce_pure=False)[source]
Calculates expectation value of hamiltonian operator at time time, in state V.
\[\langle V|H(t=\texttt{time})|V\rangle,\qquad \mathrm{tr}(V H(t=\texttt{time}))\]- 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].
- timeobj, optional
Can be either one of the following:
- float: time to evalute the time-dependent part of the operator at (if existent).
Default is time = 0.
- (N,) array_like: if V.shape[-1] == N, the hamiltonian operator is evaluated at the i-th time
and the expectation value is calculated with respect to V[…,i]. Here V must be either 2- or 3-d array, where 2-d would be for pure states and 3-d would be for mixed states.
- 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 hamiltonian operator.
- Parameters:
- copybool, optional
Whether to return a deep copy of the original object. Default is copy = False.
- Returns:
hamiltonian
\(H_{ij}\mapsto H_{ij}^*\)
Examples
>>> H_herm = H.getH()
- property get_shape
shape of the hamiltonian object, always equal to (Ns,Ns).
- Type:
tuple
- property is_dense
checks sparsity of operator matrix.
True if the operator contains a dense matrix as a component of either the static or dynamic lists.
- Type:
bool
- matrix_ele(Vl, Vr, time=0, diagonal=False, check=True)[source]
Calculates matrix element of hamiltonian operator at time time in states Vl and Vr.
\[\langle V_l|H(t=\texttt{time})|V_r\rangle\]- Parameters:
- Vl{numpy.ndarray, scipy.spmatrix}
Vector(s)/state(s) to multiple with on left side.
- Vl{numpy.ndarray, scipy.spmatrix}
Vector(s)/state(s) to multiple with on right side.
- timeobj, optional
Can be either one of the following:
- float: time to evalute the time-dependent part of the operator at (if existent).
Default is time = 0.
- (N,) array_like: if V.shape[1] == N, the hamiltonian operator is evaluated at the i-th time
and the fluctuations are calculated with respect to V[:,i]. Here V must be a 2-d array containing pure states in the columns of the array.
- diagonalbool, optional
When set to True, returs only diagonal part of expectation value. Default is diagonal = False.
- checkbool,
- Returns:
- float
Matrix element of hamiltonian 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,time=0,diagonal=False,check=True)
corresponds to \(H_{lr} = \langle V_l|H(t=0)|V_r\rangle\).
- property nbytes
Total bytes consumed by the elements of the hamiltonian array.
- Type:
float
- property ndim
number of dimensions, always equal to 2.
- Type:
int
- project_to(proj)[source]
Projects/Transforms hamiltonian operator with projector/operator proj.
Let us call the projector/transformation \(V\). Then, the function computes
\[V^\dagger H V\]- Parameters:
- projobj
Can be either one of the following:
hamiltonian object
exp_op object
numpy.ndarray
scipy.sparse array
The shape of proj need not be square, but has to comply with the matrix multiplication requirements in the definition above.
- Returns:
- obj
Projected/Transformed hamiltonian operator. The output object type depends on the object type of proj.
Notes
The proj argument can be a square array, in which case the function just transforms the hailtonian operator \(H\). Or it can be a projector which then projects \(H\) onto a smaller Hilbert space.
Projectors onto bases with symmetries other than H.basis can be conveniently obtain using the basis.get_proj() method of the basis constructor class.
Examples
>>> H_new = H.project_to(V)
correponds to \(V^\dagger H V\).
- quant_fluct(V, time=0, check=True, enforce_pure=False)[source]
Calculates the quantum fluctuations (variance) of hamiltonian operator at time time, in state V.
\[\langle V|H^2(t=\texttt{time})|V\rangle - \langle V|H(t=\texttt{time})|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].
- timeobj, optional
Can be either one of the following:
- float: time to evalute the time-dependent part of the operator at (if existent).
Default is time = 0.
- (N,) array_like: if V.shape[-1] == N, the hamiltonian operator is evaluated at the i-th time
and the fluctuations are calculated with respect to V[…,i]. Here V must be either 2- or 3-d array, where 2-d would be for pure states and 3-d would be for mixed states.
- 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, time=0, check=True, out=None, overwrite_out=True, a=1.0)[source]
Vector-Matrix multiplication of hamiltonian operator at time time, with state V.
\[a\langle V|H(t=\texttt{time})\]- Parameters:
- Vnumpy.ndarray
Vector (quantum state) to multiply the hamiltonian operator with on the left.
- timeobj, optional
Can be either one of the following:
- float: time to evalute the time-dependent part of the operator at (if existent).
Default is time = 0.
- (N,) array_like: if V.shape[-1] == N, the hamiltonian operator is evaluated at the i-th time
and the mattrix multiplication on the right is calculated with respect to V[…,i]. Here V must be either 2- or 3-d array, where 2-d would be for pure states and 3-d would be for mixed states.
- 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 or if times is an array.
- 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 hamiltonian operator applied on the state V.
Notes
this function does the matrix multiplication with the state(s) and Hamiltonian as is, see Example 17 (Lidblad dynamics / Optical Bloch Equations).
Examples
>>> B = H.rdot(A,time=0,check=True)
corresponds to \(B = AH\).
- rotate_by(other, generator=False, **exp_op_kwargs)[source]
Rotates/Transforms hamiltonian object by an operator other.
Let us denote the transformation by \(V\). With generator=False, other corresponds to the transformation \(V\), and this function implements
\[V^\dagger H V\]while for generator=True, other corresponds to a generator \(K\), and the function implements
\[\exp(a^*K^\dagger) H \exp(a K)\]- Parameters:
- otherobj
Can be either one of the following:
hamiltonian object
exp_op object
numpy.ndarray
scipy.sparse array
- generatorbool, optional
If set to True, this flag renders other a generator, and implements the calculation of
\[\exp(a^*K^\dagger) H \exp(a K)\]If set to False, the function implements
\[V^\dagger H V\]Default is generator = False.
- All other optional arguments are the same as for the `exp_op` class.
- Returns:
- obj
Transformed hamiltonian operator. The output object type depends on the object type of other.
Notes
If generator = False, this function calls project_to.
Examples
>>> H_new = H.rotate_by(V,generator=False)
corresponds to \(V^\dagger H V\).
>>> H_new = H.rotate_by(K,generator=True,**exp_op_kwargs)
corresponds to \(\exp(K^\dagger) H \exp(K)\).
- property shape
shape of the hamiltonian object, always equal to (Ns,Ns).
- Type:
tuple
- property static
static part of the operator.
- Type:
scipy.sparse.csr
- toarray(time=0, order=None, out=None)[source]
Returns copy of a hamiltonian object at time time as a dense array.
This function can overflow memory if not used carefully!
- Parameters:
- timefloat, optional
Time to evalute the time-dependent part of the operator at (if existent). Default is time = 0.0.
- orderstr, optional
Whether to store multi-dimensional data in C (rom-major) or Fortran (molumn-major) order in memory.
Default is order = None, indicating the NumPy default of C-ordered.
Cannot be specified in conjunction with the out argument.
- outnumpy.ndarray
Array to fill in with the output.
- Returns:
- numpy.ndarray
Dense array.
Examples
>>> H_dense=H.toarray(time=time)
- tocsc(time=0)[source]
Returns copy of a hamiltonian object at time time as a scipy.sparse.csc_matrix.
Casts the hamiltonian object as a scipy.sparse.csc_matrix object.
- Parameters:
- timefloat, optional
Time to evalute the time-dependent part of the operator at (if existent). Default is time = 0.0.
- Returns:
scipy.sparse.csc_matrix
Examples
>>> H_csc=H.tocsc(time=time)
- tocsr(time=0)[source]
Returns copy of a hamiltonian object at time time as a scipy.sparse.csr_matrix.
Casts the hamiltonian object as a scipy.sparse.csr_matrix object.
- Parameters:
- timefloat, optional
Time to evalute the time-dependent part of the operator at (if existent). Default is time = 0.0.
- Returns:
scipy.sparse.csr_matrix
Examples
>>> H_csr=H.tocsr(time=time)
- todense(time=0, order=None, out=None)[source]
Returns copy of a hamiltonian object at time time as a dense array.
This function can overflow memory if not used carefully!
- Parameters:
- timefloat, optional
Time to evalute the time-dependent part of the operator at (if existent). Default is time = 0.0.
- orderstr, optional
Whether to store multi-dimensional data in C (rom-major) or Fortran (molumn-major) order in memory.
Default is order = None, indicating the NumPy default of C-ordered.
Cannot be specified in conjunction with the out argument.
- 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 hamiltonian operator as a numpy.matrix instead of a numpy.ndarray. In such a case, one can use the hamiltonian.toarray() method.
Examples
>>> H_dense=H.todense(time=time)
- trace(time=0)[source]
Calculates trace of hamiltonian operator at time time.
- Parameters:
- timefloat, optional
Time to evalute the time-dependent part of the operator at (if existent). Default is time = 0.0.
- Returns:
- float
Trace of operator \(\sum_{j=1}^{Ns} H_{jj}(t=\texttt{time})\).
Examples
>>> H_tr = H.tr(time=0.0)
- transpose(copy=False)[source]
Transposes hamiltonian operator.
- Returns:
hamiltonian
\(H_{ij}\mapsto H_{ji}\)
Notes
This function does NOT conjugate the operator.
Examples
>>> H_tran = H.transpose()
- update_matrix_formats(static_fmt, dynamic_fmt)[source]
Change the internal structure of the matrices in-place.
- Parameters:
- static_fmtstr {“csr”,”csc”,”dia”,”dense”}
Specifies format of static part of Hamiltonian.
- dynamic_fmt: str {“csr”,”csc”,”dia”,”dense”} or dict, keys: (func,func_args), values: str {“csr”,”csc”,”dia”,”dense”}
Specifies the format of the dynamic parts of the hamiltonian. To specify a particular dynamic part of the hamiltonian use a tuple (func,func_args) which matches a function+argument pair used in the construction of the hamiltonian as a key in the dictionary.
- copybool,optional
Whether to return a deep copy of the original object. Default is copy = False.
Examples
make the dynamic part of the hamiltonian object to be DIA matrix format and have the static part be CSR matrix format:
>>> H.update_matrix_formats(static_fmt="csr",dynamic_fmt={(func,func_args):"dia"})