Classes

Qobj

class Qobj(inpt=None, dims=, [[][]]shape=, []type=None, isherm=None, fast=False, superrep=None)

A class for representing quantum objects, such as quantum operators and states.

The Qobj class is the QuTiP representation of quantum operators and state vectors. This class also implements math operations +,-,* between Qobj instances (and / by a C-number), as well as a collection of common operator/state operations. The Qobj constructor optionally takes a dimension list and/or shape list as arguments.

Parameters:

inpt : array_like

Data for vector/matrix representation of the quantum object.

dims : list

Dimensions of object used for tensor products.

shape : list

Shape of underlying data structure (matrix shape).

fast : bool

Flag for fast qobj creation when running ode solvers. This parameter is used internally only.

Attributes

data array_like Sparse matrix characterizing the quantum object.
dims list List of dimensions keeping track of the tensor structure.
superrep str Representation used if type is ‘super’. One of ‘super’ (Liouville form) or ‘choi’ (Choi matrix with tr = dimension).

Methods

conj() Conjugate of quantum object.
dag() Adjoint (dagger) of quantum object.
eigenenergies(sparse=False, sort=’low’, eigvals=0, tol=0, maxiter=100000) Returns eigenenergies (eigenvalues) of a quantum object.
eigenstates(sparse=False, sort=’low’, eigvals=0, tol=0, maxiter=100000) Returns eigenenergies and eigenstates of quantum object.
expm() Matrix exponential of quantum object.
full() Returns dense array of quantum object data attribute.
groundstate(sparse=False,tol=0,maxiter=100000) Returns eigenvalue and eigenket for the groundstate of a quantum object.
matrix_element(bra, ket) Returns the matrix element of operator between bra and ket vectors.
norm(norm=’tr’, sparse=False, tol=0, maxiter=100000) Returns norm of a ket or an operator.
permute(order) Returns composite qobj with indices reordered.
ptrace(sel) Returns quantum object for selected dimensions after performing partial trace.
sqrtm() Matrix square root of quantum object.
tidyup(atol=1e-12) Removes small elements from quantum object.
tr() Trace of quantum object.
trans() Transpose of quantum object.
transform(inpt, inverse=False) Performs a basis transformation defined by inpt matrix.
unit(norm=’tr’, sparse=False, tol=0, maxiter=100000) Returns normalized quantum object.
checkherm()

Check if the quantum object is hermitian.

Returns:

isherm: bool :

Returns the new value of isherm property.

conj()

Conjugate operator of quantum object.

dag()

Adjoint operator of quantum object.

diag()

Diagonal elements of quantum object.

Returns:

diags: array :

Returns array of real values if operators is Hermitian, otherwise complex values are returned.

eigenenergies(sparse=False, sort='low', eigvals=0, tol=0, maxiter=100000)

Eigenenergies of a quantum object.

Eigenenergies (eigenvalues) are defined for operators or superoperators only.

Parameters:

sparse : bool

Use sparse Eigensolver

sort : str

Sort eigenvalues ‘low’ to high, or ‘high’ to low.

eigvals : int

Number of requested eigenvalues. Default is all eigenvalues.

tol : float

Tolerance used by sparse Eigensolver (0=machine precision). The sparse solver may not converge if the tolerance is set too low.

maxiter : int

Maximum number of iterations performed by sparse solver (if used).

Returns:

eigvals: array :

Array of eigenvalues for operator.

Notes

The sparse eigensolver is much slower than the dense version. Use sparse only if memory requirements demand it.

eigenstates(sparse=False, sort='low', eigvals=0, tol=0, maxiter=100000)

Eigenstates and eigenenergies.

Eigenstates and eigenenergies are defined for operators and superoperators only.

Parameters:

sparse : bool

Use sparse Eigensolver

sort : str

Sort eigenvalues (and vectors) ‘low’ to high, or ‘high’ to low.

eigvals : int

Number of requested eigenvalues. Default is all eigenvalues.

tol : float

Tolerance used by sparse Eigensolver (0 = machine precision). The sparse solver may not converge if the tolerance is set too low.

maxiter : int

Maximum number of iterations performed by sparse solver (if used).

Returns:

eigvals : array

Array of eigenvalues for operator.

eigvecs : array

Array of quantum operators representing the oprator eigenkets. Order of eigenkets is determined by order of eigenvalues.

Notes

The sparse eigensolver is much slower than the dense version. Use sparse only if memory requirements demand it.

eliminate_states(states_inds, normalize=False)

Creates a new quantum object with states in state_inds eliminated.

Parameters:

states_inds : list of integer

The states that should be removed.

normalize : True / False

Weather or not the new Qobj instance should be normalized (default is False). For Qobjs that represents density matrices or state vectors normalized should probably be set to True, but for Qobjs that represents operators in for example an Hamiltonian, normalize should be False.

Returns:

q : qutip.Qobj

A new instance of qutip.Qobj that contains only the states corresponding to indices that are not in state_inds.

.. note:: :

Experimental.

static evaluate(qobj_list, t, args)

Evaluate a time-dependent quantum object in list format. For example,

qobj_list = [H0, [H1, func_t]]

is evaluated to

Qobj(t) = H0 + H1 * func_t(t, args)

and

qobj_list = [H0, [H1, ‘sin(w * t)’]]

is evaluated to

Qobj(t) = H0 + H1 * sin(args[‘w’] * t)
Parameters:

qobj_list : list

A nested list of Qobj instances and corresponding time-dependent coefficients.

t : float

The time for which to evaluate the time-dependent Qobj instance.

args : dictionary

A dictionary with parameter values required to evaluate the time-dependent Qobj intance.

Returns:

output : Qobj

A Qobj instance that represents the value of qobj_list at time t.

expm(method=None)

Matrix exponential of quantum operator.

Input operator must be square.

Parameters:

method : str {‘dense’, ‘sparse’, ‘scipy-dense’, ‘scipy-sparse’}

Use set method to use to calculate the matrix exponentiation. The available choices includes ‘dense’ and ‘sparse’ for using QuTiP’s implementation of expm using dense and sparse matrices, respectively, and ‘scipy-dense’ and ‘scipy-sparse’ for using the scipy.linalg.expm (dense) and scipy.sparse.linalg.expm (sparse). If no method is explicitly given a heuristic will be used to try and automatically select the most appropriate solver.

Returns:

oper : qobj

Exponentiated quantum operator.

Raises:

TypeError :

Quantum operator is not square.

extract_states(states_inds, normalize=False)

Qobj with states in state_inds only.

Parameters:

states_inds : list of integer

The states that should be kept.

normalize : True / False

Weather or not the new Qobj instance should be normalized (default is False). For Qobjs that represents density matrices or state vectors normalized should probably be set to True, but for Qobjs that represents operators in for example an Hamiltonian, normalize should be False.

Returns:

q : qutip.Qobj

A new instance of qutip.Qobj that contains only the states corresponding to the indices in state_inds.

.. note:: :

Experimental.

full(squeeze=False)

Dense array from quantum object.

Returns:

data : array

Array of complex data from quantum objects data attribute.

groundstate(sparse=False, tol=0, maxiter=100000)

Ground state Eigenvalue and Eigenvector.

Defined for quantum operators or superoperators only.

Parameters:

sparse : bool

Use sparse Eigensolver

tol : float

Tolerance used by sparse Eigensolver (0 = machine precision). The sparse solver may not converge if the tolerance is set too low.

maxiter : int

Maximum number of iterations performed by sparse solver (if used).

Returns:

eigval : float

Eigenvalue for the ground state of quantum operator.

eigvec : qobj

Eigenket for the ground state of quantum operator.

Notes

The sparse eigensolver is much slower than the dense version. Use sparse only if memory requirements demand it.

matrix_element(bra, ket)

Calculates a matrix element.

Gives the matrix element for the quantum object sandwiched between a bra and ket vector.

Parameters:

bra : qobj

Quantum object of type ‘bra’.

ket : qobj

Quantum object of type ‘ket’.

Returns:

elem : complex

Complex valued matrix element.

Raises:

TypeError :

Can only calculate matrix elements between a bra and ket quantum object.

norm(norm=None, sparse=False, tol=0, maxiter=100000)

Norm of a quantum object.

Default norm is L2-norm for kets and trace-norm for operators. Other ket and operator norms may be specified using the norm and argument.

Parameters:

norm : str

Which norm to use for ket/bra vectors: L2 ‘l2’, max norm ‘max’, or for operators: trace ‘tr’, Frobius ‘fro’, one ‘one’, or max ‘max’.

sparse : bool

Use sparse eigenvalue solver for trace norm. Other norms are not affected by this parameter.

tol : float

Tolerance for sparse solver (if used) for trace norm. The sparse solver may not converge if the tolerance is set too low.

maxiter : int

Maximum number of iterations performed by sparse solver (if used) for trace norm.

Returns:

norm : float

The requested norm of the operator or state quantum object.

Notes

The sparse eigensolver is much slower than the dense version. Use sparse only if memory requirements demand it.

overlap(state)

Overlap between two state vectors.

Gives the overlap (scalar product) for the quantum object and state state vector.

Parameters:

state : qobj

Quantum object for a state vector of type ‘ket’ or ‘bra’.

Returns:

overlap : complex

Complex valued overlap.

Raises:

TypeError :

Can only calculate overlap between a bra and ket quantum objects.

permute(order)

Permutes a composite quantum object.

Parameters:

order : list/array

List specifying new tensor order.

Returns:

P : qobj

Permuted quantum object.

i :

ptrace(sel)

Partial trace of the quantum object.

Parameters:

sel : int/list

An int or list of components to keep after partial trace.

Returns:

oper: qobj :

Quantum object representing partial trace with selected components remaining.

Notes

This function is identical to the qutip.qobj.ptrace function that has been deprecated.

sqrtm(sparse=False, tol=0, maxiter=100000)

Sqrt of a quantum operator.

Operator must be square.

Parameters:

sparse : bool

Use sparse eigenvalue/vector solver.

tol : float

Tolerance used by sparse solver (0 = machine precision).

maxiter : int

Maximum number of iterations used by sparse solver.

Returns:

oper: qobj :

Matrix square root of operator.

Raises:

TypeError :

Quantum object is not square.

Notes

The sparse eigensolver is much slower than the dense version. Use sparse only if memory requirements demand it.

tidyup(atol=None)

Removes small elements from the quantum object.

Parameters:

atol : float

Absolute tolerance used by tidyup. Default is set via qutip global settings parameters.

Returns:

oper: qobj :

Quantum object with small elements removed.

tr()

Trace of a quantum object.

Returns:

trace: float :

Returns real if operator is Hermitian, returns complex otherwise.

trans()

Transposed operator.

Returns:

oper : qobj

Transpose of input operator.

transform(inpt, inverse=False)

Basis transform defined by input array.

Input array can be a matrix defining the transformation, or a list of kets that defines the new basis.

Parameters:

inpt : array_like

A matrix or list of kets defining the transformation.

inverse : bool

Whether to return inverse transformation.

Returns:

oper : qobj

Operator in new basis.

Notes

This function is still in development.

unit(norm=None, sparse=False, tol=0, maxiter=100000)

Operator or state normalized to unity.

Uses norm from Qobj.norm().

Parameters:

norm : str

Requested norm for states / operators.

sparse : bool

Use sparse eigensolver for trace norm. Does not affect other norms.

tol : float

Tolerance used by sparse eigensolver.

maxiter: int :

Number of maximum iterations performed by sparse eigensolver.

Returns:

oper : qobj

Normalized quantum object.

eseries

class eseries(q=array(, []dtype=float64), s=array(, []dtype=float64))

Class representation of an exponential-series expansion of time-dependent quantum objects.

Attributes

ampl ndarray Array of amplitudes for exponential series.
rates ndarray Array of rates for exponential series.
dims list Dimensions of exponential series components
shape list Shape corresponding to exponential series components

Methods

value(tlist) Evaluate an exponential series at the times listed in tlist
spec(wlist) Evaluate the spectrum of an exponential series at frequencies in wlist.
tidyup() Returns a tidier version of the exponential series
spec(wlist)

Evaluate the spectrum of an exponential series at frequencies in wlist.

Parameters:

wlist : array_like

Array/list of frequenies.

Returns:

val_list : ndarray

Values of exponential series at frequencies in wlist.

tidyup(*args)

Returns a tidier version of exponential series.

value(tlist)

Evaluates an exponential series at the times listed in tlist.

Parameters:

tlist : ndarray

Times at which to evaluate exponential series.

Returns:

val_list : ndarray

Values of exponential at times in tlist.

Bloch sphere

class Bloch(fig=None, axes=None, view=None, figsize=None, background=False)

Class for plotting data on the Bloch sphere. Valid data can be either points, vectors, or qobj objects.

Attributes

axes instance {None} User supplied Matplotlib axes for Bloch sphere animation.
fig instance {None} User supplied Matplotlib Figure instance for plotting Bloch sphere.
font_color str {‘black’} Color of font used for Bloch sphere labels.
font_size int {20} Size of font used for Bloch sphere labels.
frame_alpha float {0.1} Sets transparency of Bloch sphere frame.
frame_color str {‘gray’} Color of sphere wireframe.
frame_width int {1} Width of wireframe.
point_color list {[“b”,”r”,”g”,”#CC6600”]} List of colors for Bloch sphere point markers to cycle through. i.e. By default, points 0 and 4 will both be blue (‘b’).
point_marker list {[“o”,”s”,”d”,”^”]} List of point marker shapes to cycle through.
point_size list {[25,32,35,45]} List of point marker sizes. Note, not all point markers look the same size when plotted!
sphere_alpha float {0.2} Transparency of Bloch sphere itself.
sphere_color str {‘#FFDDDD’} Color of Bloch sphere.
figsize list {[7,7]} Figure size of Bloch sphere plot. Best to have both numbers the same; otherwise you will have a Bloch sphere that looks like a football.
vector_color list {[“g”,”#CC6600”,”b”,”r”]} List of vector colors to cycle through.
vector_width int {5} Width of displayed vectors.
vector_style str {‘-|>’, ‘simple’, ‘fancy’, ‘’} Vector arrowhead style (from matplotlib’s arrow style).
vector_mutation int {20} Width of vectors arrowhead.
view list {[-60,30]} Azimuthal and Elevation viewing angles.
xlabel list {[“$x$”,”“]} List of strings corresponding to +x and -x axes labels, respectively.
xlpos list {[1.1,-1.1]} Positions of +x and -x labels respectively.
ylabel list {[“$y$”,”“]} List of strings corresponding to +y and -y axes labels, respectively.
ylpos list {[1.2,-1.2]} Positions of +y and -y labels respectively.
zlabel list {[r’$left|0right>$’,r’$left|1right>$’]} List of strings corresponding to +z and -z axes labels, respectively.
zlpos list {[1.2,-1.2]} Positions of +z and -z labels respectively.

Methods

add_annotation(state_or_vector, text, **kwargs)

Add a text or LaTeX annotation to Bloch sphere, parametrized by a qubit state or a vector.

Parameters:

state_or_vector : Qobj/array/list/tuple

Position for the annotaion. Qobj of a qubit or a vector of 3 elements.

text : str/unicode

Annotation text. You can use LaTeX, but remember to use raw string e.g. r”$langle x rangle$” or escape backslashes e.g. “$\langle x \rangle$”.

**kwargs : :

Options as for mplot3d.axes3d.text, including: fontsize, color, horizontalalignment, verticalalignment.

add_points(points, meth='s')

Add a list of data points to bloch sphere.

Parameters:

points : array/list

Collection of data points.

meth : str {‘s’, ‘m’, ‘l’}

Type of points to plot, use ‘m’ for multicolored, ‘l’ for points connected with a line.

add_states(state, kind='vector')

Add a state vector Qobj to Bloch sphere.

Parameters:

state : qobj

Input state vector.

kind : str {‘vector’,’point’}

Type of object to plot.

add_vectors(vectors)

Add a list of vectors to Bloch sphere.

Parameters:

vectors : array/list

Array with vectors of unit length or smaller.

clear()

Resets Bloch sphere data sets to empty.

make_sphere()

Plots Bloch sphere and data sets.

render(fig=None, axes=None)

Render the Bloch sphere and its data sets in on given figure and axes.

save(name=None, format='png', dirc=None)

Saves Bloch sphere to file of type format in directory dirc.

Parameters:

name : str

Name of saved image. Must include path and format as well. i.e. ‘/Users/Paul/Desktop/bloch.png’ This overrides the ‘format’ and ‘dirc’ arguments.

format : str

Format of output image.

dirc : str

Directory for output images. Defaults to current working directory.

Returns:

File containing plot of Bloch sphere. :

set_label_convention(convention)

Set x, y and z labels according to one of conventions.

Parameters:

convention : string

One of the following: - “original” - “xyz” - “sx sy sz” - “01” - “polarization jones” - “polarization jones letters”

show()

Display Bloch sphere and corresponding data sets.

vector_mutation = None

Sets the width of the vectors arrowhead

vector_style = None

Style of Bloch vectors, default = ‘-|>’ (or ‘simple’)

vector_width = None

Width of Bloch vectors, default = 5

class Bloch3d(fig=None)

Class for plotting data on a 3D Bloch sphere using mayavi. Valid data can be either points, vectors, or qobj objects corresponding to state vectors or density matrices. for a two-state system (or subsystem).

Notes

The use of mayavi for 3D rendering of the Bloch sphere comes with a few limitations: I) You can not embed a Bloch3d figure into a matplotlib window. II) The use of LaTex is not supported by the mayavi rendering engine. Therefore all labels must be defined using standard text. Of course you can post-process the generated figures later to add LaTeX using other software if needed.

Attributes

fig instance {None} User supplied Matplotlib Figure instance for plotting Bloch sphere.
font_color str {‘black’} Color of font used for Bloch sphere labels.
font_scale float {0.08} Scale for font used for Bloch sphere labels.
frame bool {True} Draw frame for Bloch sphere
frame_alpha float {0.05} Sets transparency of Bloch sphere frame.
frame_color str {‘gray’} Color of sphere wireframe.
frame_num int {8} Number of frame elements to draw.
frame_radius floats {0.005} Width of wireframe.
point_color list {[‘r’, ‘g’, ‘b’, ‘y’]} List of colors for Bloch sphere point markers to cycle through. i.e. By default, points 0 and 4 will both be blue (‘r’).
point_mode string {‘sphere’,’cone’,’cube’,’cylinder’,’point’} Point marker shapes.
point_size float {0.075} Size of points on Bloch sphere.
sphere_alpha float {0.1} Transparency of Bloch sphere itself.
sphere_color str {‘#808080’} Color of Bloch sphere.
size list {[500,500]} Size of Bloch sphere plot in pixels. Best to have both numbers the same otherwise you will have a Bloch sphere that looks like a football.
vector_color list {[‘r’, ‘g’, ‘b’, ‘y’]} List of vector colors to cycle through.
vector_width int {3} Width of displayed vectors.
view list {[45,65]} Azimuthal and Elevation viewing angles.
xlabel list {[‘|x>’, ‘’]} List of strings corresponding to +x and -x axes labels, respectively.
xlpos list {[1.07,-1.07]} Positions of +x and -x labels respectively.
ylabel list {[‘|y>’, ‘’]} List of strings corresponding to +y and -y axes labels, respectively.
ylpos list {[1.07,-1.07]} Positions of +y and -y labels respectively.
zlabel list {[‘|0>’, ‘|1>’]} List of strings corresponding to +z and -z axes labels, respectively.
zlpos list {[1.07,-1.07]} Positions of +z and -z labels respectively.

Methods

add_points(points, meth='s')

Add a list of data points to bloch sphere.

Parameters:

points : array/list

Collection of data points.

meth : str {‘s’,’m’}

Type of points to plot, use ‘m’ for multicolored.

add_states(state, kind='vector')

Add a state vector Qobj to Bloch sphere.

Parameters:

state : qobj

Input state vector.

kind : str {‘vector’,’point’}

Type of object to plot.

add_vectors(vectors)

Add a list of vectors to Bloch sphere.

Parameters:

vectors : array/list

Array with vectors of unit length or smaller.

clear()

Resets the Bloch sphere data sets to empty.

make_sphere()

Plots Bloch sphere and data sets.

plot_points()

Plots points on the Bloch sphere.

plot_vectors()

Plots vectors on the Bloch sphere.

save(name=None, format='png', dirc=None)

Saves Bloch sphere to file of type format in directory dirc.

Parameters:

name : str

Name of saved image. Must include path and format as well. i.e. ‘/Users/Paul/Desktop/bloch.png’ This overrides the ‘format’ and ‘dirc’ arguments.

format : str

Format of output image. Default is ‘png’.

dirc : str

Directory for output images. Defaults to current working directory.

Returns:

File containing plot of Bloch sphere. :

show()

Display the Bloch sphere and corresponding data sets.

Solver Options and Results

class Options(atol=1e-08, rtol=1e-06, method='adams', order=12, nsteps=1000, first_step=0, max_step=0, min_step=0, average_expect=True, average_states=False, tidy=True, num_cpus=0, norm_tol=0.001, norm_steps=5, rhs_reuse=False, rhs_filename=None, ntraj=500, gui=False, rhs_with_state=False, store_final_state=False, store_states=False, seeds=None, steady_state_average=False)

Class of options for evolution solvers such as qutip.mesolve and qutip.mcsolve. Options can be specified either as arguments to the constructor:

opts = Options(order=10, ...)

or by changing the class attributes after creation:

opts = Options()
opts.order = 10

Returns options class to be used as options in evolution solvers.

Attributes

atol float {1e-8} Absolute tolerance.
rtol float {1e-6} Relative tolerance.
method str {‘adams’,’bdf’} Integration method.
order int {12} Order of integrator (<=12 ‘adams’, <=5 ‘bdf’)
nsteps int {2500} Max. number of internal steps/call.
first_step float {0} Size of initial step (0 = automatic).
min_step float {0} Minimum step size (0 = automatic).
max_step float {0} Maximum step size (0 = automatic)
tidy bool {True,False} Tidyup Hamiltonian and initial state by removing small terms.
num_cpus int Number of cpus used by mcsolver (default = # of cpus).
norm_tol float Tolerance used when finding wavefunction norm in mcsolve.
norm_steps int Max. number of steps used to find wavefunction norm to within norm_tol in mcsolve.
average_states bool {True, False} Avg. expectation values in mcsolver.
ntraj int {500} Number of trajectories in stochastic solvers.
rhs_reuse bool {False,True} Reuse Hamiltonian data.
rhs_with_state bool {False,True} Whether or not to include the state in the Hamiltonian function callback signature.
rhs_filename str Name for compiled Cython file.
store_final_state bool {False, True} Whether or not to store the final state of the evolution in the result class.
store_states bool {False, True} Whether or not to store the state vectors or density matrices in the result class, even if expectation values operators are given. If no expectation are provided, then states are stored by default and this option has no effect.
class Result

Class for storing simulation results from any of the dynamics solvers.

Attributes

solver str Which solver was used [e.g., ‘mesolve’, ‘mcsolve’, ‘brmesolve’, ...]
times list/array Times at which simulation data was collected.
expect list/array Expectation values (if requested) for simulation.
states array State of the simulation (density matrix or ket) evaluated at times.
num_expect int Number of expectation value operators in simulation.
num_collapse int Number of collapse operators in simualation.
ntraj int/list Number of monte-carlo trajectories (if using mcsolve). List indicates that averaging of expectation values was done over a subset of total number of trajectories.
col_times list Times at which state collpase occurred. Only for Monte Carlo solver.
col_which list Which collapse operator was responsible for each collapse in col_times. Only for Monte Carlo solver.
class StochasticSolverOptions(H=None, state0=None, times=None, c_ops=, []sc_ops=, []e_ops=, []m_ops=None, args=None, ntraj=1, nsubsteps=1, d1=None, d2=None, d2_len=1, dW_factors=None, rhs=None, generate_A_ops=None, generate_noise=None, homogeneous=True, solver=None, method=None, distribution='normal', store_measurement=False, noise=None, normalize=True, options=None, progress_bar=None)

Class of options for stochastic solvers such as qutip.stochastic.ssesolve, qutip.stochastic.smesolve, etc. Options can be specified either as arguments to the constructor:

sso = StochasticSolverOptions(nsubsteps=100, ...)

or by changing the class attributes after creation:

sso = StochasticSolverOptions()
sso.nsubsteps = 1000

The stochastic solvers qutip.stochastic.ssesolve, qutip.stochastic.smesolve, qutip.stochastic.ssepdpsolve and qutip.stochastic.smepdpsolve all take the same keyword arguments as the constructor of these class, and internally they use these arguments to construct an instance of this class, so it is rarely needed to explicitly create an instance of this class.

Attributes

H qutip.Qobj System Hamiltonian.
state0 qutip.Qobj Initial state vector (ket) or density matrix.
times list / array List of times for \(t\). Must be uniformly spaced.
c_ops list of qutip.Qobj List of deterministic collapse operators.
sc_ops list of qutip.Qobj List of stochastic collapse operators. Each stochastic collapse operator will give a deterministic and stochastic contribution to the equation of motion according to how the d1 and d2 functions are defined.
e_ops list of qutip.Qobj Single operator or list of operators for which to evaluate expectation values.
m_ops list of qutip.Qobj List of operators representing the measurement operators. The expected format is a nested list with one measurement operator for each stochastic increament, for each stochastic collapse operator.
args dict / list List of dictionary of additional problem-specific parameters.
ntraj int Number of trajectors.
nsubsteps int Number of sub steps between each time-spep given in times.
d1 function Function for calculating the operator-valued coefficient to the deterministic increment dt.
d2 function Function for calculating the operator-valued coefficient to the stochastic increment(s) dW_n, where n is in [0, d2_len[.
d2_len int (default 1) The number of stochastic increments in the process.
dW_factors array Array of length d2_len, containing scaling factors for each measurement operator in m_ops.
rhs function Function for calculating the deterministic and stochastic contributions to the right-hand side of the stochastic differential equation. This only needs to be specified when implementing a custom SDE solver.
generate_A_ops function Function that generates a list of pre-computed operators or super- operators. These precomputed operators are used in some d1 and d2 functions.
generate_noise function Function for generate an array of pre-computed noise signal.
homogeneous bool (True) Wheter or not the stochastic process is homogenous. Inhomogenous processes are only supported for poisson distributions.
solver string Name of the solver method to use for solving the stochastic equations. Valid values are: ‘euler-maruyama’, ‘fast-euler-maruyama’, ‘milstein’, ‘fast-milstein’, ‘platen’.
method string (‘homodyne’, ‘heterodyne’, ‘photocurrent’) The name of the type of measurement process that give rise to the stochastic equation to solve. Specifying a method with this keyword argument is a short-hand notation for using pre-defined d1 and d2 functions for the corresponding stochastic processes.
distribution string (‘normal’, ‘poission’) The name of the distribution used for the stochastic increments.
store_measurements bool (default False) Whether or not to store the measurement results in the qutip.solver.SolverResult instance returned by the solver.
noise array Vector specifying the noise.
normalize bool (default True) Whether or not to normalize the wave function during the evolution.
options qutip.solver.Options Generic solver options.
progress_bar qutip.ui.BaseProgressBar Optional progress bar class instance.

Distribution functions

class Distribution(data=None, xvecs=, []xlabels=[])

A class for representation spatial distribution functions.

The Distribution class can be used to prepresent spatial distribution functions of arbitray dimension (although only 1D and 2D distributions are used so far).

It is indented as a base class for specific distribution function, and provide implementation of basic functions that are shared among all Distribution functions, such as visualization, calculating marginal distributions, etc.

Parameters:

data : array_like

Data for the distribution. The dimensions must match the lengths of the coordinate arrays in xvecs.

xvecs : list

List of arrays that spans the space for each coordinate.

xlabels : list

List of labels for each coordinate.

Methods

marginal(dim=0)

Calculate the marginal distribution function along the dimension dim. Return a new Distribution instance describing this reduced- dimensionality distribution.

Parameters:

dim : int

The dimension (coordinate index) along which to obtain the marginal distribution.

Returns:

d : Distributions

A new instances of Distribution that describes the marginal distribution.

project(dim=0)

Calculate the projection (max value) distribution function along the dimension dim. Return a new Distribution instance describing this reduced-dimensionality distribution.

Parameters:

dim : int

The dimension (coordinate index) along which to obtain the projected distribution.

Returns:

d : Distributions

A new instances of Distribution that describes the projection.

visualize(fig=None, ax=None, figsize=(8, 6), colorbar=True, cmap=None, style='colormap', show_xlabel=True, show_ylabel=True)

Visualize the data of the distribution in 1D or 2D, depending on the dimensionality of the underlaying distribution.

Parameters:

fig : matplotlib Figure instance
If given, use this figure instance for the visualization,
ax : matplotlib Axes instance
If given, render the visualization using this axis instance.
figsize : tuple
Size of the new Figure instance, if one needs to be created.
colorbar: Bool
Whether or not the colorbar (in 2D visualization) should be used.
cmap: matplotlib colormap instance
If given, use this colormap for 2D visualizations.
style : string
Type of visualization: ‘colormap’ (default) or ‘surface’.
Returns:

fig, ax : tuple

A tuple of matplotlib figure and axes instances.

class WignerDistribution(rho=None, extent=[[-5, 5], [-5, 5]], steps=250)

Methods

class QDistribution(rho=None, extent=[[-5, 5], [-5, 5]], steps=250)

Methods

class TwoModeQuadratureCorrelation(state=None, theta1=0.0, theta2=0.0, extent=[[-5, 5], [-5, 5]], steps=250)

Methods

update(state)

calculate probability distribution for quadrature measurement outcomes given a two-mode wavefunction or density matrix

update_psi(psi)

calculate probability distribution for quadrature measurement outcomes given a two-mode wavefunction

update_rho(rho)

calculate probability distribution for quadrature measurement outcomes given a two-mode density matrix

class HarmonicOscillatorWaveFunction(psi=None, omega=1.0, extent=[-5, 5], steps=250)

Methods

update(psi)

Calculate the wavefunction for the given state of an harmonic oscillator

class HarmonicOscillatorProbabilityFunction(rho=None, omega=1.0, extent=[-5, 5], steps=250)

Methods

update(rho)

Calculate the probability function for the given state of an harmonic oscillator (as density matrix)

Quantum information processing

class Gate(name, targets=None, controls=None, arg_value=None, arg_label=None)

Representation of a quantum gate, with its required parametrs, and target and control qubits.

class QubitCircuit(N, reverse_states=True)

Representation of a quantum program/algorithm, maintaining a sequence of gates.

Methods

add_1q_gate(name, start=0, end=None, qubits=None, arg_value=None, arg_label=None)

Adds a single qubit gate with specified parameters on a variable number of qubits in the circuit. By default, it applies the given gate to all the qubits in the register.

Parameters:

name: String :

Gate name.

start: Integer :

Starting location of qubits.

end: Integer :

Last qubit for the gate.

qubits: List :

Specific qubits for applying gates.

arg_value: Float :

Argument value(phi).

arg_label: String :

Label for gate representation.

add_circuit(qc, start=0)

Adds a block of a qubit circuit to the main circuit. Globalphase gates are not added.

Parameters:

qc: QubitCircuit :

The circuit block to be added to the main circuit.

start: Integer :

The qubit on which the first gate is applied.

add_gate(name, targets=None, controls=None, arg_value=None, arg_label=None)

Adds a gate with specified parameters to the circuit.

Parameters:

name: String :

Gate name.

targets: List :

Gate targets.

controls: List :

Gate controls.

arg_value: Float :

Argument value(phi).

arg_label: String :

Label for gate representation.

adjacent_gates()

Method to resolve two qubit gates with non-adjacent control/s or target/s in terms of gates with adjacent interactions.

Returns:

qc: QubitCircuit :

Returns QubitCircuit of resolved gates for the qubit circuit in the desired basis.

propagators()

Propagator matrix calculator for N qubits returning the individual steps as unitary matrices operating from left to right.

Returns:

U_list: list :

Returns list of unitary matrices for the qubit circuit.

remove_gate(index=None, name=None, remove='first')

Removes a gate with from a specific index or the first, last or all instances of a particular gate.

Parameters:

index: Integer :

Location of gate to be removed.

name: String :

Gate name to be removed.

remove: String :

If first or all gate are to be removed.

resolve_gates(basis=['CNOT', 'RX', 'RY', 'RZ'])

Unitary matrix calculator for N qubits returning the individual steps as unitary matrices operating from left to right in the specified basis.

Parameters:

basis: list. :

Basis of the resolved circuit.

Returns:

qc: QubitCircuit :

Returns QubitCircuit of resolved gates for the qubit circuit in the desired basis.

reverse_circuit()

Reverses an entire circuit of unitary gates.

Returns:

qc: QubitCircuit :

Returns QubitCircuit of resolved gates for the qubit circuit in the desired basis.

class CircuitProcessor(N, correct_global_phase)

Base class for representation of the physical implementation of a quantum program/algorithm on a specified qubit system.

Methods

adjacent_gates(qc, setup)

Function to take a quantum circuit/algorithm and convert it into the optimal form/basis for the desired physical system.

Parameters:

qc: QubitCircuit :

Takes the quantum circuit to be implemented.

setup: String :

Takes the nature of the spin chain; linear or circular.

Returns:

qc: QubitCircuit :

The resolved circuit representation.

get_ops_and_u()

Returns the Hamiltonian operators and corresponding values by stacking them together.

get_ops_labels()

Returns the Hamiltonian operators and corresponding labels by stacking them together.

load_circuit(qc)

Translates an abstract quantum circuit to its corresponding Hamiltonian for a specific model.

Parameters:

qc: QubitCircuit :

Takes the quantum circuit to be implemented.

optimize_circuit(qc)

Function to take a quantum circuit/algorithm and convert it into the optimal form/basis for the desired physical system.

Parameters:

qc: QubitCircuit :

Takes the quantum circuit to be implemented.

Returns:

qc: QubitCircuit :

The optimal circuit representation.

plot_pulses()

Maps the physical interaction between the circuit components for the desired physical system.

Returns:

fig, ax: Figure :

Maps the physical interaction between the circuit components.

pulse_matrix()

Generates the pulse matrix for the desired physical system.

Returns:

t, u, labels: :

Returns the total time and label for every operation.

run(qc=None)

Generates the propagator matrix by running the Hamiltonian for the appropriate time duration for the desired physical system.

Parameters:

qc: QubitCircuit :

Takes the quantum circuit to be implemented.

Returns:

U_list: list :

The propagator matrix obtained from the physical implementation.

run_state(qc=None, states=None)

Generates the propagator matrix by running the Hamiltonian for the appropriate time duration for the desired physical system with the given initial state of the qubit register.

Parameters:

qc: QubitCircuit :

Takes the quantum circuit to be implemented.

states: Qobj :

Initial state of the qubits in the register.

Returns:

U_list: list :

The propagator matrix obtained from the physical implementation.

class SpinChain(N, correct_global_phase=True, sx=None, sz=None, sxsy=None)

Representation of the physical implementation of a quantum program/algorithm on a spin chain qubit system.

Methods

adjacent_gates(qc, setup='linear')

Method to resolve 2 qubit gates with non-adjacent control/s or target/s in terms of gates with adjacent interactions for linear/circular spin chain system.

Parameters:

qc: QubitCircuit :

The circular spin chain circuit to be resolved

setup: Boolean :

Linear of Circular spin chain setup

Returns:

qc: QubitCircuit :

Returns QubitCircuit of resolved gates for the qubit circuit in the desired basis.

class LinearSpinChain(N, correct_global_phase=True, sx=None, sz=None, sxsy=None)

Representation of the physical implementation of a quantum program/algorithm on a spin chain qubit system arranged in a linear formation. It is a sub-class of SpinChain.

Methods

class CircularSpinChain(N, correct_global_phase=True, sx=None, sz=None, sxsy=None)

Representation of the physical implementation of a quantum program/algorithm on a spin chain qubit system arranged in a circular formation. It is a sub-class of SpinChain.

Methods

class DispersivecQED(N, correct_global_phase=True, Nres=None, deltamax=None, epsmax=None, w0=None, eps=None, delta=None, g=None)

Representation of the physical implementation of a quantum program/algorithm on a dispersive cavity-QED system.

Methods

dispersive_gate_correction(qc1, rwa=True)

Method to resolve ISWAP and SQRTISWAP gates in a cQED system by adding single qubit gates to get the correct output matrix.

Parameters:

qc: Qobj :

The circular spin chain circuit to be resolved

rwa: Boolean :

Specify if RWA is used or not.

Returns:

qc: QubitCircuit :

Returns QubitCircuit of resolved gates for the qubit circuit in the desired basis.