Generates the vector representation of a Fock state.
Parameters: | N : int
n : int
offset : int (default 0)
|
---|---|
Returns: | state : qobj
|
Notes
A subtle incompatibility with the quantum optics toolbox: In QuTiP:
basis(N, 0) = ground state
but in the qotoolbox:
basis(N, 1) = ground state
Examples
>>> basis(5,2)
Quantum object: dims = [[5], [1]], shape = [5, 1], type = ket
Qobj data =
[[ 0.+0.j]
[ 0.+0.j]
[ 1.+0.j]
[ 0.+0.j]
[ 0.+0.j]]
Generates a coherent state with eigenvalue alpha.
Constructed using displacement operator on vacuum state.
Parameters: | N : int
alpha : float/complex
offset : int (default 0)
method : string {‘operator’, ‘analytic’}
|
---|---|
Returns: | state : qobj
|
Notes
Select method ‘operator’ (default) or ‘analytic’. With the ‘operator’ method, the coherent state is generated by displacing the vacuum state using the displacement operator defined in the truncated Hilbert space of size ‘N’. This method guarantees that the resulting state is normalized. With ‘analytic’ method the coherent state is generated using the analytical formula for the coherent state coefficients in the Fock basis. This method does not guarantee that the state is normalized if truncated to a small number of Fock states, but would in that case give more accurate coefficients.
Examples
>>> coherent(5,0.25j)
Quantum object: dims = [[5], [1]], shape = [5, 1], type = ket
Qobj data =
[[ 9.69233235e-01+0.j ]
[ 0.00000000e+00+0.24230831j]
[ -4.28344935e-02+0.j ]
[ 0.00000000e+00-0.00618204j]
[ 7.80904967e-04+0.j ]]
Density matrix representation of a coherent state.
Constructed via outer product of qutip.states.coherent
Parameters: | N : int
alpha : float/complex
offset : int (default 0)
method : string {‘operator’, ‘analytic’}
|
---|---|
Returns: | dm : qobj
|
Notes
Select method ‘operator’ (default) or ‘analytic’. With the ‘operator’ method, the coherent density matrix is generated by displacing the vacuum state using the displacement operator defined in the truncated Hilbert space of size ‘N’. This method guarantees that the resulting density matrix is normalized. With ‘analytic’ method the coherent density matrix is generated using the analytical formula for the coherent state coefficients in the Fock basis. This method does not guarantee that the state is normalized if truncated to a small number of Fock states, but would in that case give more accurate coefficients.
Examples
>>> coherent_dm(3,0.25j)
Quantum object: dims = [[3], [3]], shape = [3, 3], type = oper, isHerm = True
Qobj data =
[[ 0.93941695+0.j 0.00000000-0.23480733j -0.04216943+0.j ]
[ 0.00000000+0.23480733j 0.05869011+0.j 0.00000000-0.01054025j]
[-0.04216943+0.j 0.00000000+0.01054025j 0.00189294+0.j ]]
Bosonic Fock (number) state.
Same as qutip.states.basis.
Parameters: | N : int
n : int
|
---|---|
Returns: | Requested number state :math:`left|nright>`. : |
Examples
>>> fock(4,3)
Quantum object: dims = [[4], [1]], shape = [4, 1], type = ket
Qobj data =
[[ 0.+0.j]
[ 0.+0.j]
[ 0.+0.j]
[ 1.+0.j]]
Density matrix representation of a Fock state
Constructed via outer product of qutip.states.fock.
Parameters: | N : int
n : int
|
---|---|
Returns: | dm : qobj
|
Examples
>>> fock_dm(3,1)
Quantum object: dims = [[3], [3]], shape = [3, 3], type = oper, isHerm = True
Qobj data =
[[ 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 1.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j]]
Takes input ket or bra vector and returns density matrix formed by outer product.
Parameters: | Q : qobj
|
---|---|
Returns: | dm : qobj
|
Examples
>>> x=basis(3,2)
>>> ket2dm(x)
Quantum object: dims = [[3], [3]], shape = [3, 3], type = oper, isHerm = True
Qobj data =
[[ 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 1.+0.j]]
Basis states for a three level system (qutrit)
Returns: | qstates : array
|
---|
Density matrix for a thermal state of n particles
Parameters: | N : int
n : float
method : string {‘operator’, ‘analytic’}
|
---|---|
Returns: | dm : qobj
|
Notes
The ‘operator’ method (default) generates the thermal state using the truncated number operator num(N). This is the method that should be used in computations. The ‘analytic’ method uses the analytic coefficients derived in an infinite Hilbert space. The analytic form is not necessarily normalized, if truncated too aggressively.
Examples
>>> thermal_dm(5, 1)
Quantum object: dims = [[5], [5]], shape = [5, 5], type = oper, isHerm = True
Qobj data =
[[ 0.51612903 0. 0. 0. 0. ]
[ 0. 0.25806452 0. 0. 0. ]
[ 0. 0. 0.12903226 0. 0. ]
[ 0. 0. 0. 0.06451613 0. ]
[ 0. 0. 0. 0. 0.03225806]]
>>> thermal_dm(5, 1, 'analytic')
Quantum object: dims = [[5], [5]], shape = [5, 5], type = oper, isHerm = True
Qobj data =
[[ 0.5 0. 0. 0. 0. ]
[ 0. 0.25 0. 0. 0. ]
[ 0. 0. 0.125 0. 0. ]
[ 0. 0. 0. 0.0625 0. ]
[ 0. 0. 0. 0. 0.03125]]
An iterator that enumerate all the state number arrays (quantum numbers on the form [n1, n2, n3, ...]) for a system with dimensions given by dims.
Example:
>>> for state in state_number_enumerate([2,2]):
>>> print state
[ 0. 0.]
[ 0. 1.]
[ 1. 0.]
[ 1. 1.]
Parameters: | dims : list or array
state : list
idx : integer
|
---|---|
Returns: | state_number : list
|
Return the index of a quantum state corresponding to state, given a system with dimensions given by dims.
Example:
>>> state_number_index([2, 2, 2], [1, 1, 0])
6.0
Parameters: | dims : list or array
state : list
|
---|---|
Returns: | idx : list
|
Return a quantum number representation given a state index, for a system of composite structure defined by dims.
Example:
>>> state_index_number([2, 2, 2], 6)
[1, 1, 0]
Parameters: | dims : list or array
index : integer
|
---|---|
Returns: | state : list
|
Return a Qobj representation of a quantum state specified by the state array state.
Example:
>>> state_number_qobj([2, 2, 2], [1, 0, 1])
Quantum object: dims = [[2, 2, 2], [1, 1, 1]], shape = [8, 1], type = ket
Qobj data =
[[ 0.]
[ 0.]
[ 0.]
[ 0.]
[ 0.]
[ 1.]
[ 0.]
[ 0.]]
Parameters: | dims : list or array
state : list
|
---|---|
Returns: | state : qutip.Qobj.qobj
|
Basis vector for the mth phase of the Pegg-Barnett phase operator.
Parameters: | N : int
m : int
phi0 : float (default=0)
|
---|---|
Returns: | state : qobj
|
Notes
The Pegg-Barnett basis states form a complete set over the truncated Hilbert space.
Creation (raising) operator.
Parameters: | N : int
|
---|---|
Returns: | oper : qobj
offset : int (default 0)
|
Examples
>>> create(4)
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isHerm = False
Qobj data =
[[ 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j]
[ 1.00000000+0.j 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j]
[ 0.00000000+0.j 1.41421356+0.j 0.00000000+0.j 0.00000000+0.j]
[ 0.00000000+0.j 0.00000000+0.j 1.73205081+0.j 0.00000000+0.j]]
Destruction (lowering) operator.
Parameters: | N : int
offset : int (default 0)
|
---|---|
Returns: | oper : qobj
|
Examples
>>> destroy(4)
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isHerm = False
Qobj data =
[[ 0.00000000+0.j 1.00000000+0.j 0.00000000+0.j 0.00000000+0.j]
[ 0.00000000+0.j 0.00000000+0.j 1.41421356+0.j 0.00000000+0.j]
[ 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j 1.73205081+0.j]
[ 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j]]
Single-mode displacement operator.
Parameters: | N : int
alpha : float/complex
offset : int (default 0)
|
---|---|
Returns: | oper : qobj
|
Examples
>>> displace(4,0.25)
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isHerm = False
Qobj data =
[[ 0.96923323+0.j -0.24230859+0.j 0.04282883+0.j -0.00626025+0.j]
[ 0.24230859+0.j 0.90866411+0.j -0.33183303+0.j 0.07418172+0.j]
[ 0.04282883+0.j 0.33183303+0.j 0.84809499+0.j -0.41083747+0.j]
[ 0.00626025+0.j 0.07418172+0.j 0.41083747+0.j 0.90866411+0.j]]
Higher-order spin operators:
Parameters: | j : float
args : str
|
---|---|
Returns: | jmat : qobj/list
|
Notes
If no ‘args’ input, then returns array of [‘x’,’y’,’z’] operators.
Examples
>>> jmat(1)
[ Quantum object: dims = [[3], [3]], shape = [3, 3], type = oper, isHerm = True
Qobj data =
[[ 0. 0.70710678 0. ]
[ 0.70710678 0. 0.70710678]
[ 0. 0.70710678 0. ]]
Quantum object: dims = [[3], [3]], shape = [3, 3], type = oper, isHerm = True
Qobj data =
[[ 0.+0.j 0.+0.70710678j 0.+0.j ]
[ 0.-0.70710678j 0.+0.j 0.+0.70710678j]
[ 0.+0.j 0.-0.70710678j 0.+0.j ]]
Quantum object: dims = [[3], [3]], shape = [3, 3], type = oper, isHerm = True
Qobj data =
[[ 1. 0. 0.]
[ 0. 0. 0.]
[ 0. 0. -1.]]]
Quantum object for number operator.
Parameters: | N : int
offset : int (default 0)
|
---|---|
Returns: | oper: qobj :
|
Examples
>>> num(4)
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isHerm = True
Qobj data =
[[0 0 0 0]
[0 1 0 0]
[0 0 2 0]
[0 0 0 3]]
Identity operator
Parameters: | N : int
|
---|---|
Returns: | oper : qobj
|
Examples
>>> qeye(3)
Quantum object: dims = [[3], [3]], shape = [3, 3], type = oper, isHerm = True
Qobj data =
[[ 1. 0. 0.]
[ 0. 1. 0.]
[ 0. 0. 1.]]
Identity operator. Alternative name to qeye.
Parameters: | N : int
|
---|---|
Returns: | oper : qobj
|
Operators for a three level system (qutrit).
Returns: | opers: array :
|
---|
Annihilation operator for Pauli spins.
Examples
>>> sigmam()
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = False
Qobj data =
[[ 0. 0.]
[ 1. 0.]]
Creation operator for Pauli spins.
Examples
>>> sigmam()
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = False
Qobj data =
[[ 0. 1.]
[ 0. 0.]]
Pauli spin 1/2 sigma-x operator
Examples
>>> sigmax()
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = False
Qobj data =
[[ 0. 1.]
[ 1. 0.]]
Pauli spin 1/2 sigma-y operator.
Examples
>>> sigmay()
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = True
Qobj data =
[[ 0.+0.j 0.-1.j]
[ 0.+1.j 0.+0.j]]
Pauli spin 1/2 sigma-z operator.
Examples
>>> sigmaz()
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = True
Qobj data =
[[ 1. 0.]
[ 0. -1.]]
Single-mode Squeezing operator.
Parameters: | N : int
sp : float/complex
offset : int (default 0)
|
---|---|
Returns: | oper : qutip.qobj.Qobj
|
Examples
>>> squeeze(4,0.25)
Quantum object: dims = [[4], [4]], shape = [4, 4], type = oper, isHerm = False
Qobj data =
[[ 0.98441565+0.j 0.00000000+0.j 0.17585742+0.j 0.00000000+0.j]
[ 0.00000000+0.j 0.95349007+0.j 0.00000000+0.j 0.30142443+0.j]
[-0.17585742+0.j 0.00000000+0.j 0.98441565+0.j 0.00000000+0.j]
[ 0.00000000+0.j -0.30142443+0.j 0.00000000+0.j 0.95349007+0.j]]
Generalized squeezing operator.
Parameters: | a1 : qutip.qobj.Qobj
a2 : qutip.qobj.Qobj
z : float/complex
|
---|---|
Returns: | oper : qutip.qobj.Qobj
|
Single-mode Pegg-Barnett phase operator.
Parameters: | N : int
phi0 : float
|
---|---|
Returns: | oper : qobj
|
Notes
The Pegg-Barnett phase operator is Hermitian on a truncated Hilbert space.
This module is a collection of random state and operator generators. The sparsity of the ouput Qobj’s is controlled by varing the density parameter.
Creates a random NxN density matrix.
Parameters: | N : int
density : float
dims : list
|
---|---|
Returns: | oper : qobj
|
Notes
For small density matrices., choosing a low density will result in an error as no diagonal elements will be generated such that \(Tr(\rho)=1\).
Creates a random NxN sparse Hermitian quantum object.
Uses \(H=X+X^{+}\) where \(X\) is a randomly generated quantum operator with a given density.
Parameters: | N : int
density : float
dims : list
|
---|---|
Returns: | oper : qobj
|
Creates a random Nx1 sparse ket vector.
Parameters: | N : int
density : float
dims : list
|
---|---|
Returns: | oper : qobj
|
Creates a random NxN sparse unitary quantum object.
Uses \(\exp(-iH)\) where H is a randomly generated Hermitian operator.
Parameters: | N : int
density : float
dims : list
|
---|---|
Returns: | oper : qobj
|
This module provides functions that are useful for simulating the three level atom with QuTiP. A three level atom (qutrit) has three states, which are linked by dipole transitions so that 1 <-> 2 <-> 3. Depending on there relative energies they are in the ladder, lambda or vee configuration. The structure of the relevant operators is the same for any of the three configurations:
Ladder: Lambda: Vee:
|two> |three>
-------|three> ------- -------
| / \ |one> /
| / \ ------- /
| / \ \ /
-------|two> / \ \ /
| / \ \ /
| / \ \ /
| / -------- \ /
-------|one> ------- |three> -------
|one> |two>
Create a vector representation of a quantum operator given the matrix representation.
Create a matrix representation given a quantum operator in vector form.
Assembles the Liouvillian superoperator from a Hamiltonian and a list of collapse operators. Like liouvillian, but with an experimental implementation which avoids creating extra Qobj instances, which can be advantageous for large systems.
Parameters: | H : qobj
c_ops : array_like
|
---|---|
Returns: | L : qobj
|
Superoperator formed from post-multiplication by operator A
Parameters: | A : qobj
|
---|---|
Returns: | super : qobj
|
Superoperator formed from pre-multiplication by operator A.
Parameters: | A : qobj
|
---|---|
Returns: | super :qobj :
|
Lindblad dissipator (generalized) for a single pair of collapse operators (a, b), or for a single collapse operator (a) when b is not specified:
Parameters: | a : qobj
b : qobj (optional)
|
---|---|
Returns: | D : qobj
|
This module implements transformations between superoperator representations, including supermatrix, Kraus, Choi and Chi (process) matrix formalisms.
Converts a Qobj representing a quantum map to the Choi representation, such that the trace of the returned operator is equal to the dimension of the system.
Parameters: | q_oper : Qobj
|
---|---|
Returns: | choi : Qobj
|
Raises: | TypeError: if the given quantum object is not a map, or cannot be converted :
|
Converts a Qobj representing a quantum map to the supermatrix (Liouville) representation.
Parameters: | q_oper : Qobj
|
---|---|
Returns: | superop : Qobj
|
Raises: | TypeError: if the given quantum object is not a map, or cannot be converted :
|
Converts a Qobj representing a quantum map to a list of quantum objects, each representing an operator in the Kraus decomposition of the given map.
Parameters: | q_oper : Qobj
|
---|---|
Returns: | kraus_ops : list of Qobj
|
Raises: | TypeError: if the given quantum object is not a map, or cannot be :
|
Module for the creation of composite quantum objects via the tensor product.
Calculates the tensor product of input operators.
Parameters: | args : array_like
|
---|---|
Returns: | obj : qobj
|
Examples
>>> tensor([sigmax(), sigmax()])
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = True
Qobj data =
[[ 0.+0.j 0.+0.j 0.+0.j 1.+0.j]
[ 0.+0.j 0.+0.j 1.+0.j 0.+0.j]
[ 0.+0.j 1.+0.j 0.+0.j 0.+0.j]
[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j]]
Calculates the expectation value for operator(s) and state(s).
Parameters: | oper : qobj/array-like
state : qobj/array-like
|
---|---|
Returns: | expt : float/complex/array-like
|
Examples
>>> expect(num(4), basis(4, 3))
3
Variance of an operator for the given state vector or density matrix.
Parameters: | oper : qobj
state : qobj/list
|
---|---|
Returns: | var : float
|
Return the partial transpose of a Qobj instance rho, where mask is an array/list with length that equals the number of components of rho (that is, the length of rho.dims[0]), and the values in mask indicates whether or not the corresponding subsystem is to be transposed. The elements in mask can be boolean or integers 0 or 1, where True/1 indicates that the corresponding subsystem should be tranposed.
Parameters: | rho : qutip.qobj
mask : list / array
method : str
|
---|---|
Returns: | rho_pr: :class:`qutip.qobj` :
|
Calculate the concurrence entanglement measure for a two-qubit state.
Parameters: | state : qobj
|
---|---|
Returns: | concur : float
|
References
[R2] | http://en.wikipedia.org/wiki/Concurrence_(quantum_computing) |
Calculates the conditional entropy \(S(A|B)=S(A,B)-S(B)\) of a slected density matrix component.
Parameters: | rho : qobj
selB : int/list
base : {e,2}
sparse : {False,True}
|
---|---|
Returns: | ent_cond : float
|
Linear entropy of a density matrix.
Parameters: | rho : qobj
|
---|---|
Returns: | entropy : float
|
Examples
>>> rho=0.5*fock_dm(2,0)+0.5*fock_dm(2,1)
>>> entropy_linear(rho)
0.5
Calculates the mutual information S(A:B) between selection components of a system density matrix.
Parameters: | rho : qobj
selA : int/list
selB : int/list
base : {e,2}
sparse : {False,True}
|
---|---|
Returns: | ent_mut : float
|
Von-Neumann entropy of density matrix
Parameters: | rho : qobj
base : {e,2}
sparse : {False,True}
|
---|---|
Returns: | entropy : float
|
Examples
>>> rho=0.5*fock_dm(2,0)+0.5*fock_dm(2,1)
>>> entropy_vn(rho,2)
1.0
This module contains a collection of functions for calculating metrics (distance measures) between states and operators.
Calculates the fidelity (pseudo-metric) between two density matrices. See: Nielsen & Chuang, “Quantum Computation and Quantum Information”
Parameters: | A : qobj
B : qobj
|
---|---|
Returns: | fid : float
|
Examples
>>> x=fock_dm(5,3)
>>> y=coherent_dm(5,1)
>>> fidelity(x,y)
0.24104350624628332
Calculates the trace distance between two density matrices.. See: Nielsen & Chuang, “Quantum Computation and Quantum Information”
Parameters: | A : qobj
B : qobj
tol : float
sparse : {False, True}
|
---|---|
Returns: | tracedist : float
|
Examples
>>> x=fock_dm(5,3)
>>> y=coherent_dm(5,1)
>>> tracedist(x,y)
0.9705143161472971
This module contains a collection functions for calculating continuous variable quantities from fock-basis representation of the state of multi-mode fields.
Given a basis set of operators \(\{a\}_n\), calculate the correlation matrix:
Parameters: | basis : list of qutip.qobj.Qobj
rho : qutip.qobj.Qobj
|
---|---|
Returns: | corr_mat: *array* :
|
Given a basis set of operators \(\{a\}_n\), calculate the covariance matrix:
or, if of the optional argument symmetrized=False,
Parameters: | basis : list of qutip.qobj.Qobj
rho : qutip.qobj.Qobj
symmetrized : bool
|
---|---|
Returns: | corr_mat: *array* :
|
Calculate the correlation matrix for given field operators \(a_1\) and \(a_2\). If a density matrix is given the expectation values are calculated, otherwise a matrix with operators is returned.
Parameters: | a1 : qutip.qobj.Qobj
a2 : qutip.qobj.Qobj
rho : qutip.qobj.Qobj
|
---|---|
Returns: | cov_mat: *array* of complex numbers or :class:`qutip.qobj.Qobj` :
|
Calculate the quadrature correlation matrix with given field operators \(a_1\) and \(a_2\). If a density matrix is given the expectation values are calculated, otherwise a matrix with operators is returned.
Parameters: | a1 : qutip.qobj.Qobj
a2 : qutip.qobj.Qobj
rho : qutip.qobj.Qobj
|
---|---|
Returns: | corr_mat: *array* of complex numbers or :class:`qutip.qobj.Qobj` :
|
Calculate the Wigner covariance matrix \(V_{ij} = \frac{1}{2}(R_{ij} + R_{ji})\), given the quadrature correlation matrix \(R_{ij} = \langle R_{i} R_{j}\rangle - \langle R_{i}\rangle \langle R_{j}\rangle\), where \(R = (q_1, p_1, q_2, p_2)^T\) is the vector with quadrature operators for the two modes.
Alternatively, if R = None, and if annihilation operators a1 and a2 for the two modes are supplied instead, the quadrature correlation matrix is constructed from the annihilation operators before then the covariance matrix is calculated.
Parameters: | a1 : qutip.qobj.Qobj
a2 : qutip.qobj.Qobj
R : array
rho : qutip.qobj.Qobj
|
---|---|
Returns: | cov_mat: *array* :
|
Calculate the logarithmic negativity given the symmetrized covariance matrix, see qutip.continous_variables.covariance_matrix. Note that the two-mode field state that is described by V must be Gaussian for this function to applicable.
Parameters: | V : 2d array
|
---|---|
Returns: | N: *float*, the logarithmic negativity for the two-mode Gaussian state : that is described by the the Wigner covariance matrix V. : |
This module provides solvers for the unitary Schrodinger equation.
Schrodinger equation evolution of a state vector for a given Hamiltonian.
Evolve the state vector or density matrix (rho0) using a given Hamiltonian (H), by integrating the set of ordinary differential equations that define the system.
The output is either the state vector at arbitrary points in time (tlist), or the expectation values of the supplied operators (e_ops). If e_ops is a callback function, it is invoked for each time in tlist with time and the state as arguments, and the function does not use any return values.
Parameters: | H : qutip.qobj
rho0 : qutip.qobj
tlist : list / array
e_ops : list of qutip.qobj / callback function single
args : dictionary
options : qutip.Qdeoptions
|
---|---|
Returns: | output: :class:`qutip.solver` :
|
This module provides solvers for the Lindblad master equation and von Neumann equation.
Master equation evolution of a density matrix for a given Hamiltonian.
Evolve the state vector or density matrix (rho0) using a given Hamiltonian (H) and an [optional] set of collapse operators (c_op_list), by integrating the set of ordinary differential equations that define the system. In the absence of collapse operators the system is evolved according to the unitary evolution of the Hamiltonian.
The output is either the state vector at arbitrary points in time (tlist), or the expectation values of the supplied operators (e_ops). If e_ops is a callback function, it is invoked for each time in tlist with time and the state as arguments, and the function does not use any return values.
Time-dependent operators
For problems with time-dependent problems H and c_ops can be callback functions that takes two arguments, time and args, and returns the Hamiltonian or Liouvillian for the system at that point in time (callback format).
Alternatively, H and c_ops can be a specified in a nested-list format where each element in the list is a list of length 2, containing an operator (qutip.qobj) at the first element and where the second element is either a string (list string format), a callback function (list callback format) that evaluates to the time-dependent coefficient for the corresponding operator, or a numpy array (list array format) which specifies the value of the coefficient to the corresponding operator for each value of t in tlist.
Examples
H = [[H0, ‘sin(w*t)’], [H1, ‘sin(2*w*t)’]]
H = [[H0, sin(w*tlist)], [H1, sin(2*w*tlist)]]
H = [[H0, f0_t], [H1, f1_t]]
where f0_t and f1_t are python functions with signature f_t(t, args).
In the list string format and list callback format, the string expression and the callback function must evaluate to a real or complex number (coefficient for the corresponding operator).
In all cases of time-dependent operators, args is a dictionary of parameters that is used when evaluating operators. It is passed to the callback functions as second argument
Note
If an element in the list-specification of the Hamiltonian or the list of collapse operators are in super-operator for it will be added to the total Liouvillian of the problem with out further transformation. This allows for using mesolve for solving master equations that are not on standard Lindblad form.
Note
On using callback function: mesolve transforms all qutip.qobj objects to sparse matrices before handing the problem to the integrator function. In order for your callback function to work correctly, pass all qutip.qobj objects that are used in constructing the Hamiltonian via args. mesolve will check for qutip.qobj in args and handle the conversion to sparse matrices. All other qutip.qobj objects that are not passed via args will be passed on to the integrator in scipy which will raise an NotImplemented exception.
Parameters: | H : qutip.Qobj
rho0 : qutip.Qobj
tlist : list / array
c_ops : list of qutip.Qobj
e_ops : list of qutip.Qobj / callback function single
args : dictionary
options : qutip.Options
|
---|---|
Returns: | output: :class:`qutip.solver` :
|
Monte-Carlo evolution of a state vector \(|\psi \rangle\) for a given Hamiltonian and sets of collapse operators, and possibly, operators for calculating expectation values. Options for the underlying ODE solver are given by the Options class.
mcsolve supports time-dependent Hamiltonians and collapse operators using either Python functions of strings to represent time-dependent coefficients. Note that, the system Hamiltonian MUST have at least one constant term.
As an example of a time-dependent problem, consider a Hamiltonian with two terms H0 and H1, where H1 is time-dependent with coefficient sin(w*t), and collapse operators C0 and C1, where C1 is time-dependent with coeffcient exp(-a*t). Here, w and a are constant arguments with values W and A.
Using the Python function time-dependent format requires two Python functions, one for each collapse coefficient. Therefore, this problem could be expressed as:
def H1_coeff(t,args):
return sin(args['w']*t)
def C1_coeff(t,args):
return exp(-args['a']*t)
H=[H0,[H1,H1_coeff]]
c_op_list=[C0,[C1,C1_coeff]]
args={'a':A,'w':W}
or in String (Cython) format we could write:
H=[H0,[H1,'sin(w*t)']]
c_op_list=[C0,[C1,'exp(-a*t)']]
args={'a':A,'w':W}
Constant terms are preferably placed first in the Hamiltonian and collapse operator lists.
Parameters: | H : qobj
psi0 : qobj
tlist : array_like
ntraj : int
c_ops : array_like
e_ops : array_like
args : dict
options : Options
|
---|---|
Returns: | results : Result
|
Monte-Carlo wave function solver with fortran 90 backend. Usage is identical to qutip.mcsolve, for problems without explicit time-dependence, and with some optional input:
Parameters: | H : qobj
psi0 : qobj
tlist : array_like
ntraj : int
c_ops : array_like
e_ops : array_like
options : Options
sparse_dms : boolean
serial : boolean
ptrace_sel: list :
calc_entropy : boolean
|
---|---|
Returns: | results : Result
|
Evolution of a state vector or density matrix (rho0) for a given Hamiltonian (H) and set of collapse operators (c_op_list), by expressing the ODE as an exponential series. The output is either the state vector at arbitrary points in time (tlist), or the expectation values of the supplied operators (e_ops).
Parameters: | H : qobj/function_type
rho0 : qutip.qobj
tlist : list/array
c_op_list : list of qutip.qobj
e_ops : list of qutip.qobj
|
---|---|
Returns: | expt_array : array
.. note:: This solver does not support time-dependent Hamiltonians. : |
Creates an exponential series that describes the time evolution for the initial density matrix (or state vector) rho0, given the Liouvillian (or Hamiltonian) L.
Parameters: | L : qobj
rho0 : qobj
|
---|---|
Returns: | eseries : qutip.eseries
|
Solve the dynamics for the system using the Bloch-Redfeild master equation.
Note
This solver does not currently support time-dependent Hamiltonian or collapse operators.
Parameters: | H : qutip.qobj
rho0 / psi0: :class:`qutip.qobj` :
tlist : list / array
a_ops : list of qutip.qobj
e_ops : list of qutip.qobj / callback function
args : dictionary
options : qutip.Qdeoptions
|
---|---|
Returns: | output: :class:`qutip.solver` :
|
Calculate the Bloch-Redfield tensor for a system given a set of operators and corresponding spectral functions that describes the system’s coupling to its environment.
Parameters: | H : qutip.qobj
a_ops : list of qutip.qobj
spectra_cb : list of callback functions
use_secular : bool
|
---|---|
Returns: | R, kets: :class:`qutip.qobj`, list of :class:`qutip.qobj` :
|
Evolve the ODEs defined by Bloch-Redfield master equation. The Bloch-Redfield tensor can be calculated by the function bloch_redfield_tensor.
Parameters: | R : qutip.qobj
ekets : array of qutip.qobj
rho0 : qutip.qobj
tlist : list / array
e_ops : list of qutip.qobj / callback function
options : qutip.Qdeoptions
|
---|---|
Returns: | output: :class:`qutip.solver` :
|
Solve the dynamics for the system using the Floquet-Markov master equation.
Note
This solver currently does not support multiple collapse operators.
Parameters: | H : qutip.qobj
rho0 / psi0 : qutip.qobj
tlist : list / array
c_ops : list of qutip.qobj
e_ops : list of qutip.qobj / callback function
spectra_cb : list callback functions
T : float
args : dictionary
options : qutip.solver
k_max : int
|
---|---|
Returns: | output : qutip.solver
|
Calculate the initial Floquet modes Phi_alpha(0) for a driven system with period T.
Returns a list of qutip.qobj instances representing the Floquet modes and a list of corresponding quasienergies, sorted by increasing quasienergy in the interval [-pi/T, pi/T]. The optional parameter sort decides if the output is to be sorted in increasing quasienergies or not.
Parameters: | H : qutip.qobj
args : dictionary
T : float
U : qutip.qobj
|
---|---|
Returns: | output : list of kets, list of quasi energies
|
Calculate the Floquet modes at times tlist Phi_alpha(tlist) propagting the initial Floquet modes Phi_alpha(0)
Parameters: | f_modes_0 : list of qutip.qobj (kets)
f_energies : list
t : float
H : qutip.qobj
args : dictionary
T : float
|
---|---|
Returns: | output : list of kets
|
Pre-calculate the Floquet modes for a range of times spanning the floquet period. Can later be used as a table to look up the floquet modes for any time.
Parameters: | f_modes_0 : list of qutip.qobj (kets)
f_energies : list
tlist : array
H : qutip.qobj
T : float
args : dictionary
|
---|---|
Returns: | output : nested list
|
Lookup the floquet mode at time t in the pre-calculated table of floquet modes in the first period of the time-dependence.
Parameters: | f_modes_table_t : nested list of qutip.qobj (kets)
t : float
T : float
|
---|---|
Returns: | output : nested list
|
Evaluate the floquet states at time t given the initial Floquet modes.
Parameters: | f_modes_t : list of qutip.qobj (kets)
f_energies : array
t : float
H : qutip.qobj
T : float
args : dictionary
|
---|---|
Returns: | output : list
|
Evaluate the wavefunction for a time t using the Floquet state decompositon, given the initial Floquet modes.
Parameters: | f_modes_t : list of qutip.qobj (kets)
f_energies : array
f_coeff : array
t : float
H : qutip.qobj
T : float
args : dictionary
|
---|---|
Returns: | output : qutip.qobj
|
Decompose the wavefunction psi (typically an initial state) in terms of the Floquet states, \(\psi = \sum_\alpha c_\alpha \psi_\alpha(0)\).
Parameters: | f_states : list of qutip.qobj (kets)
f_energies : array
psi : qutip.qobj
|
---|---|
Returns: | output : array
|
Solve the Schrodinger equation using the Floquet formalism.
Parameters: | H : qutip.qobj.Qobj
psi0 : qutip.qobj
tlist : list / array
e_ops : list of qutip.qobj / callback function
T : float
args : dictionary
Tsteps : integer
|
---|---|
Returns: | output : qutip.solver.Result
|
This module contains experimental functions for solving stochastic schrodinger and master equations. The API should not be considered stable, and is subject to change when we work more on optimizing this module for performance and features.
Todo:
- parallelize
Solve stochastic master equation. Dispatch to specific solvers depending on the value of the solver keyword argument.
Parameters: | H : qutip.Qobj
rho0 : qutip.Qobj
times : list / array
c_ops : list of qutip.Qobj
sc_ops : list of qutip.Qobj
e_ops : list of qutip.Qobj / callback function single
kwargs : dictionary
|
---|---|
Returns: | output: :class:`qutip.solver.SolverResult` :
|
Solve stochastic Schrödinger equation. Dispatch to specific solvers depending on the value of the solver keyword argument.
Parameters: | H : qutip.Qobj
psi0 : qutip.Qobj
times : list / array
sc_ops : list of qutip.Qobj
e_ops : list of qutip.Qobj
kwargs : dictionary
|
---|---|
Returns: | output: :class:`qutip.solver.SolverResult` :
|
A stochastic (piecewse deterministic process) PDP solver for density matrix evolution.
Parameters: | H : qutip.Qobj
rho0 : qutip.Qobj
times : list / array
c_ops : list of qutip.Qobj
sc_ops : list of qutip.Qobj
e_ops : list of qutip.Qobj / callback function single
kwargs : dictionary
|
---|---|
Returns: | output: :class:`qutip.solver.SolverResult` :
|
A stochastic (piecewse deterministic process) PDP solver for wavefunction evolution. For most purposes, use qutip.mcsolve instead for quantum trajectory simulations.
Parameters: | H : qutip.Qobj
psi0 : qutip.Qobj
times : list / array
c_ops : list of qutip.Qobj
e_ops : list of qutip.Qobj / callback function single
kwargs : dictionary
|
---|---|
Returns: | output: :class:`qutip.solver.SolverResult` :
|
Calculate a two-operator two-time correlation function on the form \(\left<A(t+\tau)B(t)\right>\) or \(\left<A(t)B(t+\tau)\right>\) (if reverse=True), using the quantum regression theorem and the evolution solver indicated by the solver parameter.
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
tlist : list / array
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj
b_op : qutip.qobj
solver : str
|
---|---|
Returns: | corr_mat: *array* :
|
Calculate a two-operator two-time correlation function \(\left<A(\tau)B(0)\right>\) or \(\left<A(0)B(\tau)\right>\) (if reverse=True), using the quantum regression theorem and the evolution solver indicated by the solver parameter.
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj.Qobj
b_op : qutip.qobj.Qobj
reverse : bool
solver : str
|
---|---|
Returns: | corr_vec: *array* :
|
Calculate a two-operator two-time correlation function \(\left<A(\tau)B(0)\right>\) or \(\left<A(0)B(\tau)\right>\) (if reverse=True), using the quantum regression theorem and the evolution solver indicated by the solver parameter.
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj.Qobj
b_op : qutip.qobj.Qobj
reverse : bool
solver : str
|
---|---|
Returns: | corr_vec: *array* :
|
Calculate a two-operator two-time correlation function on the form \(\left<A(t+\tau)B(t)\right>\) or \(\left<A(t)B(t+\tau)\right>\) (if reverse=True), using the quantum regression theorem and the evolution solver indicated by the solver parameter.
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
tlist : list / array
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj.Qobj
b_op : qutip.qobj.Qobj
solver : str
reverse : bool
|
---|---|
Returns: | corr_mat: *array* :
|
Calculate the four-operator two-time correlation function on the from \(\left<A(0)B(\tau)C(\tau)D(0)\right>\) using the quantum regression theorem and the solver indicated by the ‘solver’ parameter.
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj.Qobj
b_op : qutip.qobj.Qobj
c_op : qutip.qobj.Qobj
d_op : qutip.qobj.Qobj
solver : str
|
---|---|
Returns: | corr_vec: *array* :
|
References
See, Gardiner, Quantum Noise, Section 5.2.1.
Calculate the four-operator two-time correlation function on the from \(\left<A(t)B(t+\tau)C(t+\tau)D(t)\right>\) using the quantum regression theorem and the solver indicated by the ‘solver’ parameter.
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
tlist : list / array
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj.Qobj
b_op : qutip.qobj.Qobj
c_op : qutip.qobj.Qobj
d_op : qutip.qobj.Qobj
solver : str
|
---|---|
Returns: | corr_mat: *array* :
|
References
See, Gardiner, Quantum Noise, Section 5.2.1.
Calculate the spectrum corresponding to a correlation function \(\left<A(\tau)B(0)\right>\), i.e., the Fourier transform of the correlation function:
Parameters: | H : qutip.qobj
wlist : list / array
c_ops : list of qutip.qobj
a_op : qutip.qobj
b_op : qutip.qobj
|
---|---|
Returns: | spectrum: *array* :
|
Calculate the spectrum corresponding to a correlation function \(\left<A(\tau)B(0)\right>\), i.e., the Fourier transform of the correlation function:
Parameters: | H : qutip.qobj
wlist : list / array
c_ops : list of qutip.qobj
a_op : qutip.qobj
b_op : qutip.qobj
|
---|---|
Returns: | s_vec: *array* :
|
Calculate the power spectrum corresponding to a two-time correlation function using FFT.
Parameters: | tlist : list / array
y : list / array
|
---|---|
Returns: | w, S : tuple
|
Calculate the first-order quantum coherence function:
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj.Qobj
solver : str
|
---|---|
Returns: | g1, G2: tuble of *array* :
|
Calculate the second-order quantum coherence function:
Parameters: | H : qutip.qobj.Qobj
rho0 : qutip.qobj.Qobj
taulist : list / array
c_ops : list of qutip.qobj.Qobj
a_op : qutip.qobj.Qobj
solver : str
|
---|---|
Returns: | g2, G2: tuble of *array* :
|
Module contains functions for solving for the steady state density matrix of open quantum systems defined by a Liouvillian or Hamiltonian and a list of collapse operators.
Calculates the steady state for quantum evolution subject to the supplied Hamiltonian or Liouvillian operator and (if given a Hamiltonian) a list of collapse operators.
If the user passes a Hamiltonian then it, along with the list of collapse operators, will be converted into a Liouvillian operator in Lindblad form.
Parameters: | A : qobj
c_op_list : list
method : str {‘direct’, ‘eigen’, ‘iterative-bicg’, ‘iterative-gmres’, ‘svd’, ‘power’}
sparse : bool, optional, default=True
use_rcm : bool, optional, default=True
use_wbm : bool, optional, default=False
weight : float, optional
use_umfpack : bool {False, True}
maxiter : int, optional, default=10000
tol : float, optional, default=1e-9
permc_spec : str, optional, default=’NATURAL’
use_precond : bool optional, default = True
M : {sparse matrix, dense matrix, LinearOperator}, optional
fill_factor : float, optional, default=10
drop_tol : float, optional, default=1e-3
diag_pivot_thresh : float, optional, default=None
ILU_MILU : str, optional, default=’smilu_2’
|
---|---|
Returns: | dm : qobj
|
Notes
The SVD method works only for dense operators (i.e. small systems).
Calculate the propagator U(t) for the density matrix or wave function such that \(\psi(t) = U(t)\psi(0)\) or \(\rho_{\mathrm vec}(t) = U(t) \rho_{\mathrm vec}(0)\) where \(\rho_{\mathrm vec}\) is the vector representation of the density matrix.
Parameters: | H : qobj or list
t : float or array-like
c_op_list : list
args : list/array/dictionary
options : qutip.Options
|
---|---|
Returns: | a : qobj
|
Find the steady state for successive applications of the propagator \(U\).
Parameters: | U : qobj
|
---|---|
Returns: | a : qobj
|
Generates the Cython functions needed for solving the dynamics of a given system using the mesolve function inside a parfor loop.
Parameters: | H : qobj
c_ops : list
args : dict
options : Options
name: str :
cleanup: bool :
|
---|
Notes
Using this function with any solver other than the mesolve function will result in an error.
Resets the string-format time-dependent Hamiltonian parameters.
Returns: | Nothing, just clears data from internal config module. : |
---|
Q-function of a given state vector or density matrix at points xvec + i * yvec.
Parameters: | state : qobj
xvec : array_like
yvec : array_like
g : float
|
---|---|
Returns: | Q : array
|
Wigner function for a state vector or density matrix at points xvec + i * yvec.
Parameters: | state : qobj
xvec : array_like
yvec : array_like
g : float
method : string {‘iterative’, ‘laguerre’, ‘fft’}
parfor : bool {False, True}
|
---|---|
Returns: | W : array
yvex : array
|
Notes
The ‘fft’ method accepts only an xvec input for the x-coordinate. The y-coordinates are calculated internally.
References
Ulf Leonhardt, Measuring the Quantum State of Light, (Cambridge University Press, 1997)
Functions for visualizing results of quantum dynamics simulations, visualizations of quantum states and processes.
Draws a Hinton diagram for visualizing a density matrix.
Parameters: | rho : qobj
xlabels : list of strings
ylabels : list of strings
title : string
ax : a matplotlib axes instance
|
---|---|
Returns: | fig, ax : tuple
|
Raises: | ValueError :
|
Draw a histogram for the matrix M, with the given x and y labels and title.
Parameters: | M : Matrix of Qobj
xlabels : list of strings
ylabels : list of strings
title : string
limits : list/array with two float numbers
ax : a matplotlib axes instance
|
---|---|
Returns: | fig, ax : tuple
|
Raises: | ValueError :
|
Draw a histogram for the amplitudes of matrix M, using the argument of each element for coloring the bars, with the given x and y labels and title.
Parameters: | M : Matrix of Qobj
xlabels : list of strings
ylabels : list of strings
title : string
limits : list/array with two float numbers
phase_limits : list/array with two float numbers
ax : a matplotlib axes instance
threshold: float (None) :
|
---|---|
Returns: | fig, ax : tuple
|
Raises: | ValueError :
|
Plot the energy level diagrams for a list of Hamiltonians. Include up to N energy levels. For each element in H_list, the energy levels diagram for the cummulative Hamiltonian sum(H_list[0:n]) is plotted, where n is the index of an element in H_list.
Parameters: | H_list : List of Qobj
|
---|---|
Returns: | fig, ax : tuple
|
Raises: | ValueError :
|
A custom colormap that emphasizes negative values by creating a nonlinear colormap.
Parameters: | W : array
levels : int
shift : float
invert : bool
|
---|---|
Returns: | Returns a Matplotlib colormap instance for use in plotting. : |
Notes
The ‘shift’ parameter allows you to vary where the colormap begins to highlight negative colors. This is beneficial in cases where there are small negative Wigner elements due to numerical round-off and/or truncation.
Plot the Fock distribution for a density matrix (or ket) that describes an oscillator mode.
Parameters: | rho : qutip.qobj.Qobj
fig : a matplotlib Figure instance
ax : a matplotlib axes instance
title : string
figsize : (width, height)
|
---|---|
Returns: | fig, ax : tuple
|
Plot the Fock distribution and the Wigner function for a density matrix (or ket) that describes an oscillator mode.
Parameters: | rho : qutip.qobj.Qobj
fig : a matplotlib Figure instance
axes : a list of two matplotlib axes instances
figsize : (width, height)
cmap : a matplotlib cmap instance
alpha_max : float
colorbar : bool
method : string {‘iterative’, ‘laguerre’, ‘fft’}
projection: string {‘2d’, ‘3d’} :
|
---|---|
Returns: | fig, ax : tuple
|
Plot the the Wigner function for a density matrix (or ket) that describes an oscillator mode.
Parameters: | rho : qutip.qobj.Qobj
fig : a matplotlib Figure instance
ax : a matplotlib axes instance
figsize : (width, height)
cmap : a matplotlib cmap instance
alpha_max : float
colorbar : bool
method : string {‘iterative’, ‘laguerre’, ‘fft’}
projection: string {‘2d’, ‘3d’} :
|
---|---|
Returns: | fig, ax : tuple
|
Plots a matrix of values on a sphere
Parameters: | theta : float
phi : float
values : array
fig : a matplotlib Figure instance
ax : a matplotlib axes instance
save : bool {False , True}
|
---|---|
Returns: | fig, ax : tuple
|
Plotting scheme related to Schmidt decomposition. Converts a state into a matrix (A_ij -> A_i^j), where rows are first particles and columns - last.
See also: plot_qubism with how=’before_after’ for a similar plot.
Parameters: | ket : Qobj
splitting : int
theme : ‘light’ (default) or ‘dark’
labels_iteration : int or pair of ints (default (3,2))
fig : a matplotlib figure instance
ax : a matplotlib axis instance
figsize : (width, height)
|
---|---|
Returns: | fig, ax : tuple
|
Qubism plot for pure states of many qudits. Works best for spin chains, especially with even number of particles of the same dimension. Allows to see entanglement between first 2*k particles and the rest.
Parameters: | ket : Qobj
theme : ‘light’ (default) or ‘dark’
how : ‘pairs’ (default), ‘pairs_skewed’ or ‘before_after’
grid_iteration : int (default 1)
legend_iteration : int (default 0) or ‘grid_iteration’ or ‘all’
fig : a matplotlib figure instance
ax : a matplotlib axis instance
figsize : (width, height)
|
---|---|
Returns: | fig, ax : tuple
|
Visualize the results (expectation values) for an evolution solver. results is assumed to be an instance of Result, or a list of Result instances.
Parameters: | results : (list of) qutip.solver.Result
ylabels : list of strings
title : string
show_legend : bool
fig : a matplotlib Figure instance
axes : a matplotlib axes instance
figsize : (width, height)
|
---|---|
Returns: | fig, ax : tuple
|
Calculates an angular wave function on a sphere. psi = orbital(theta,phi,ket1,ket2,...) calculates the angular wave function on a sphere at the mesh of points defined by theta and phi which is \(\sum_{lm} c_{lm} Y_{lm}(theta,phi)\) where \(C_{lm}\) are the coefficients specified by the list of kets. Each ket has 2l+1 components for some integer l.
Parameters: | theta : list/array
phi : list/array
args : list/array
|
---|---|
Returns: | ``array`` for angular wave function : |
Calculate the quantum process tomography chi matrix for a given (possibly nonunitary) transformation matrix U, which transforms a density matrix in vector form according to:
vec(rho) = U * vec(rho0)
or
rho = vec2mat(U * mat2vec(rho0))
U can be calculated for an open quantum system using the QuTiP propagator function.
Parameters: | U : Qobj
op_basis_list : list
|
---|---|
Returns: | chi : array
|
Visualize the quantum process tomography chi matrix. Plot the real and imaginary parts separately.
Parameters: | chi : array
lbls_list : list
title : string
fig : figure instance
axes : list of figure axis instance
|
---|---|
Returns: | fig, ax : tuple
|
Visualize the quantum process tomography chi matrix. Plot bars with height and color corresponding to the absolute value and phase, respectively.
Parameters: | chi : array
lbls_list : list
title : string
fig : figure instance
ax : figure axis instance
threshold: float (None) :
|
---|---|
Returns: | fig, ax : tuple
|
Single-qubit rotation for operator sigmax with angle phi.
Returns: | result : qobj
|
---|
Single-qubit rotation for operator sigmay with angle phi.
Returns: | result : qobj
|
---|
Single-qubit rotation for operator sigmaz with angle phi.
Returns: | result : qobj
|
---|
Single-qubit square root NOT gate.
Returns: | result : qobj
|
---|
Quantum object representing the SNOT (Hadamard) gate.
Returns: | snot_gate : qobj
|
---|
Examples
>>> snot()
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = True
Qobj data =
[[ 0.70710678+0.j 0.70710678+0.j]
[ 0.70710678+0.j -0.70710678+0.j]]
Returns quantum object representing the phase shift gate.
Parameters: | theta : float
|
---|---|
Returns: | phase_gate : qobj
|
Examples
>>> phasegate(pi/4)
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = False
Qobj data =
[[ 1.00000000+0.j 0.00000000+0.j ]
[ 0.00000000+0.j 0.70710678+0.70710678j]]
Returns quantum object representing the phase shift gate.
Parameters: | theta : float
N : integer
control : integer
target : integer
|
---|---|
Returns: | U : qobj
|
Quantum object representing the CNOT gate.
Returns: | cnot_gate : qobj
|
---|
Examples
>>> cnot()
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = True
Qobj data =
[[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 1.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 1.+0.j]
[ 0.+0.j 0.+0.j 1.+0.j 0.+0.j]]
Quantum object representing the CSIGN gate.
Returns: | csign_gate : qobj
|
---|
Examples
>>> csign()
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = True
Qobj data =
[[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 1.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 1.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j -1.+0.j]]
Quantum object representing the Berkeley gate.
Returns: | berkeley_gate : qobj
|
---|
Examples
>>> berkeley()
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = True
Qobj data =
[[ cos(pi/8).+0.j 0.+0.j 0.+0.j 0.+sin(pi/8).j]
[ 0.+0.j cos(3pi/8).+0.j 0.+sin(3pi/8).j 0.+0.j]
[ 0.+0.j 0.+sin(3pi/8).j cos(3pi/8).+0.j 0.+0.j]
[ 0.+sin(pi/8).j 0.+0.j 0.+0.j cos(pi/8).+0.j]]
Quantum object representing the SWAPalpha gate.
Returns: | swapalpha_gate : qobj
|
---|
Examples
>>> swapalpha(alpha)
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = True
Qobj data =
[[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.5*(1 + exp(j*pi*alpha) 0.5*(1 - exp(j*pi*alpha) 0.+0.j]
[ 0.+0.j 0.5*(1 - exp(j*pi*alpha) 0.5*(1 + exp(j*pi*alpha) 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 1.+0.j]]
Quantum object representing the SWAP gate.
Returns: | swap_gate : qobj
|
---|
Examples
>>> swap()
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = True
Qobj data =
[[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 1.+0.j 0.+0.j]
[ 0.+0.j 1.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 1.+0.j]]
Quantum object representing the iSWAP gate.
Returns: | iswap_gate : qobj
|
---|
Examples
>>> iswap()
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = False
Qobj data =
[[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+1.j 0.+0.j]
[ 0.+0.j 0.+1.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 1.+0.j]]
Quantum object representing the square root SWAP gate.
Returns: | sqrtswap_gate : qobj
|
---|
Quantum object representing the square root iSWAP gate.
Returns: | sqrtiswap_gate : qobj
|
---|
Examples
>>> sqrtiswap()
Quantum object: dims = [[2, 2], [2, 2]], shape = [4, 4], type = oper, isHerm = False
Qobj data =
[[ 1.00000000+0.j 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j]
[ 0.00000000+0.j 0.70710678+0.j 0.00000000-0.70710678j 0.00000000+0.j]
[ 0.00000000+0.j 0.00000000-0.70710678j 0.70710678+0.j 0.00000000+0.j]
[ 0.00000000+0.j 0.00000000+0.j 0.00000000+0.j 1.00000000+0.j]]
Quantum object representing the Fredkin gate.
Returns: | fredkin_gate : qobj
|
---|
Examples
>>> fredkin()
Quantum object: dims = [[2, 2, 2], [2, 2, 2]], shape = [8, 8], type = oper, isHerm = True
Qobj data =
[[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j]]
Quantum object representing the Toffoli gate.
Returns: | toff_gate : qobj
|
---|
Examples
>>> toffoli()
Quantum object: dims = [[2, 2, 2], [2, 2, 2]], shape = [8, 8], type = oper, isHerm = True
Qobj data =
[[ 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j 0.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j]
[ 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 0.+0.j 1.+0.j 0.+0.j]]
Single-qubit rotation for operator op with angle phi.
Returns: | result : qobj
|
---|
Create an N-qubit controlled gate from a single-qubit gate U with the given control and target qubits.
Parameters: | U : Qobj
N : integer
control : integer
target : integer
control_value : integer (1)
|
---|---|
Returns: | result : qobj
|
Returns quantum object representing the global phase shift gate.
Parameters: | theta : float
|
---|---|
Returns: | phase_gate : qobj
|
Examples
>>> phasegate(pi/4)
Quantum object: dims = [[2], [2]], shape = [2, 2], type = oper, isHerm = False
Qobj data =
[[ 0.70710678+0.70710678j 0.00000000+0.j]
[ 0.00000000+0.j 0.70710678+0.70710678j]]
Quantum object representing the N-qubit Hadamard gate.
Returns: | q : qobj
|
---|
Calculate the overall unitary matrix for a given list of unitary operations
Parameters: | U_list : list
left_to_right: Boolean :
|
---|---|
Returns: | U_overall: qobj :
|
Create a Qobj representing a one-qubit gate that act on a system with N qubits.
Parameters: | U : Qobj
N : integer
target : integer
|
---|---|
Returns: | gate : qobj
|
Create a Qobj representing a two-qubit gate that act on a system with N qubits.
Parameters: | U : Qobj
N : integer
control : integer
target : integer
targets : list
|
---|---|
Returns: | gate : qobj
|
Create a Qobj representing a three-qubit gate that act on a system with N qubits.
Parameters: | U : Qobj
N : integer
controls : list
target : integer
|
---|---|
Returns: | gate : qobj
|
Function to define initial state of the qubits.
Parameters: | N: Integer :
states: List :
|
---|---|
Returns: | qstates: Qobj :
|
Quantum Fourier Transform operator on N qubits.
Parameters: | N : int
|
---|---|
Returns: | QFT: qobj :
|
Quantum Fourier Transform operator on N qubits returning the individual steps as unitary matrices operating from left to right.
Parameters: | N: int :
swap: boolean :
|
---|---|
Returns: | U_step_list: list of qobj :
|
Quantum Fourier Transform operator on N qubits returning the gate sequence.
Parameters: | N: int :
swap: boolean :
|
---|---|
Returns: | qc: instance of QubitCircuit :
|
This module contains a collection of graph theory routines used mainly to reorder matrices for iterative steady state solvers.
Breadth-First-Search (BFS) of a graph in CSR or CSC matrix format starting from a given node (row). Takes Qobjs and CSR or CSC matrices as inputs.
This function requires a matrix with symmetric structure. Use A+trans(A) if original matrix is not symmetric or not sure.
Parameters: | A : csc_matrix, csr_matrix
start : int
|
---|---|
Returns: | order : array
levels : array
|
Returns the degree for the nodes (rows) of a symmetric graph in sparse CSR or CSC format, or a qobj.
Parameters: | A : qobj, csr_matrix, csc_matrix
|
---|---|
Returns: | degree : array
|
This module contains utility functions that are commonly needed in other qutip modules.
Return the number of photons in thermal equilibrium for an harmonic oscillator mode with frequency ‘w’, at the temperature described by ‘w_th’ where \(\omega_{\rm th} = k_BT/\hbar\).
Parameters: | w : float or array
w_th : float
|
---|---|
Returns: | n_avg : float or array
|
Return an array of numbers sampled over specified interval with additional elements added.
Returns num spaced array with elements from elems inserted if not already included in set.
Returned sample array is not evenly spaced if addtional elements are added.
Parameters: | start : int
stop : int
num : int, optional
elems : list/ndarray, optional
|
---|---|
Returns: | samples : ndadrray
|
Calculates the Clebsch-Gordon coefficient for coupling (j1,m1) and (j2,m2) to give (j3,m3).
Parameters: | j1 : float
j2 : float
j3 : float
m1 : float
m2 : float
m3 : float
|
---|---|
Returns: | cg_coeff : float
|
Convert an energy from unit orig to unit to.
Parameters: | value : float / array
orig : string
to : string
|
---|---|
Returns: | value_new_unit : float / array
|
Retrieves an array of data from the requested file.
Parameters: | filename : str
sep : str
|
---|---|
Returns: | data : array_like
|
Stores a matrix of data to a file to be read by an external program.
Parameters: | filename : str
data: array_like :
numtype : str {‘complex, ‘real’}
numformat : str {‘decimal’,’exp’}
sep : str
|
---|
Loads data file from file named ‘filename.qu’ in current directory.
Parameters: | name : str
|
---|---|
Returns: | qobject : instance / array_like
|
Saves given data to file named ‘filename.qu’ in current directory.
Parameters: | data : instance/array_like
filename : str
|
---|
This module contains utility functions for using QuTiP with IPython notebooks.
Call the function tast for each value in task_vec using a cluster of IPython engines. The function task should have the signature task(value, args) or task(value) if args=None.
The client and view are the IPython.parallel client and load-balanced view that will be used in the parfor execution. If these are None, new instances will be created.
Parameters: | task: a Python function :
task_vec: array / list :
args: list / dictionary :
client: IPython.parallel.Client :
view: a IPython.parallel.Client view :
show_scheduling: bool {False, True}, default False :
show_progressbar: bool {False, True}, default False :
|
---|---|
Returns: | result : list
|
Print an HTML-formatted table with version numbers for QuTiP and its dependencies. Use it in a IPython notebook to show which versions of different packages that were used to run the notebook. This should make it possible to reproduce the environment and the calculation later on.
Returns: | version_table: string :
|
---|
Executes a multi-variable function in parallel on the local machine.
Parallel execution of a for-loop over function func for multiple input arguments and keyword arguments.
Parameters: | func : function_type
The following keyword argument is reserved: : num_cpus : int
|
---|---|
Returns: | result : list
|
About box for qutip. Gives version numbers for QuTiP, NumPy, SciPy, Cython, and MatPlotLib.
Simulateous diagonalization of communting Hermitian matrices..
Parameters: | ops : list/array
|
---|---|
Returns: | eigs : tuple
|