quimb.gen.operators¶
Functions for generating quantum operators.
Attributes¶
Functions¶
|
Generate a general spin-operator. |
|
Generates the pauli operators for dimension 2 or 3. |
|
The Hadamard gate. |
|
The generalized qubit phase-gate, which adds phase |
|
The T-gate (pi/8 gate). |
|
The S-gate (phase gate). |
|
The single qubit rotation gate. |
|
Arbitrary unitary single qubit gate. |
|
Rx(pi / 2). |
|
Ry(pi / 2). |
|
Rz(pi / 2). |
|
R[X + Y](pi / 2). |
|
The SWAP operator acting on subsystems of dimension dim. |
|
The 'fermionic simulation' gate: |
|
The 'fermionic simulation' gate, with: |
|
|
|
Build an n-qubit controlled gate. The control qubits are the |
|
Construct a controlled pauli gate for two qubits. |
|
The controlled-not gate. |
|
The controlled-X gate. |
|
The controlled-Y gate. |
|
The controlled-Z gate. |
|
The double controlled X gate, or Toffoli gate. |
|
The double controlled Y gate. |
|
The double controlled Z gate. |
|
The controlled swap or Fredkin gate. The control qubit is the first |
Wrap a function to perform some generic postprocessing and take the |
|
|
Constructs the nearest neighbour 1d heisenberg spin-1/2 hamiltonian. |
|
Generate the quantum transverse field ising model hamiltonian. This is a |
|
Generate the quantum transverse field XY model hamiltonian. This is a |
|
Generate the XXZ-model hamiltonian. This is a |
|
Generate the j1-j2 hamiltonian, i.e. next nearest neighbour |
|
|
|
Constructs a heisenberg hamiltonian with isotropic coupling and |
|
Construct the 2D spin-1/2 heisenberg model hamiltonian: |
|
Generate all the unique permutations of sequence |
|
Construct the projector onto spin-z subpspaces. |
|
The creation operator acting on an n-level system. |
|
The annihilation operator acting on an n-level system. |
|
The number operator acting on an n-level system. |
|
Generate the spinless fermion hopping hamiltonian. |
Module Contents¶
- quimb.gen.operators.spin_operator(label, S=1 / 2, **kwargs)[source]¶
Generate a general spin-operator.
- Parameters:
label (str) –
The type of operator, can be one of six options:
{'x', 'X'}
, x-spin operator.{'y', 'Y'}
, y-spin operator.{'z', 'Z'}
, z-spin operator.{'+', 'p'}
, Raising operator.{'-', 'm'}
, Lowering operator.{'i', 'I'}
, identity operator.
S (float, optional) – The spin of particle to act on, default to spin-1/2.
kwargs – Passed to
quimbify()
.
- Returns:
S – The spin operator.
- Return type:
immutable operator
See also
Examples
>>> spin_operator('x') qarray([[0. +0.j, 0.5+0.j], [0.5+0.j, 0. +0.j]])
>>> qu.spin_operator('+', S=1) qarray([[0. +0.j, 1.41421356+0.j, 0. +0.j], [0. +0.j, 0. +0.j, 1.41421356+0.j], [0. +0.j, 0. +0.j, 0. +0.j]])
>>> qu.spin_operator('Y', sparse=True) <2x2 sparse matrix of type '<class 'numpy.complex128'>' with 2 stored elements in Compressed Sparse Row format>
- quimb.gen.operators.pauli(xyz, dim=2, **kwargs)[source]¶
Generates the pauli operators for dimension 2 or 3.
- Parameters:
- Returns:
P – The pauli operator.
- Return type:
immutable operator
See also
- quimb.gen.operators.phase_gate(phi, dtype=complex, sparse=False)[source]¶
The generalized qubit phase-gate, which adds phase
phi
to the|1>
state.
- quimb.gen.operators.rotation(phi, xyz='Z', dtype=complex, sparse=False)[source]¶
The single qubit rotation gate.
- quimb.gen.operators.Rx¶
- quimb.gen.operators.Ry¶
- quimb.gen.operators.Rz¶
- quimb.gen.operators.U_gate(theta, phi, lamda, dtype=complex, sparse=False)[source]¶
Arbitrary unitary single qubit gate.
\[\begin{split}U_3(\theta, \phi, \lambda) = \begin{bmatrix} \cos(\theta / 2) & - e^{i \lambda} \sin(\theta / 2) \\ e^{i \phi} \sin(\theta / 2) & e^{i(\lambda + \phi)}\cos(\theta / 2) \end{bmatrix}\end{split}\]
- quimb.gen.operators.Xsqrt(**qu_opts)[source]¶
Rx(pi / 2).
\[\begin{split}X^{\frac{1}{2}} = \frac{1}{\sqrt{2}} \begin{bmatrix} 1 & - i \\ - i & 1 \end{bmatrix}\end{split}\]
- quimb.gen.operators.Ysqrt(**qu_opts)[source]¶
Ry(pi / 2).
\[\begin{split}Y^{\frac{1}{2}} = \frac{1}{\sqrt{2}} \begin{bmatrix} 1 & - 1 \\ 1 & 1 \end{bmatrix}\end{split}\]
- quimb.gen.operators.Zsqrt(**qu_opts)[source]¶
Rz(pi / 2).
\[\begin{split}Z^{\frac{1}{2}} = \frac{1}{\sqrt{2}} \begin{bmatrix} 1 - i & 0 \\ 0 & 1 + i \end{bmatrix}\end{split}\]
- quimb.gen.operators.Wsqrt(**qu_opts)[source]¶
R[X + Y](pi / 2).
\[\begin{split}W^{\frac{1}{2}} = \frac{1}{\sqrt{2}} \begin{bmatrix} 1 & -\sqrt{i} \\ \sqrt{-i} & 1 \end{bmatrix}\end{split}\]
- quimb.gen.operators.swap(dim=2, dtype=complex, **kwargs)[source]¶
The SWAP operator acting on subsystems of dimension dim.
- quimb.gen.operators.fsim(theta, phi, dtype=complex, **kwargs)[source]¶
The ‘fermionic simulation’ gate:
\[\begin{split}\mathrm{fsim}(\theta, \phi) = \begin{bmatrix} 1 & 0 & 0 & 0\\ 0 & \cos(\theta) & -i sin(\theta) & 0\\ 0 & -i sin(\theta) & \cos(\theta) & 0\\ 0 & 0 & 0 & \exp(-i \phi) \end{bmatrix}\end{split}\]Note that
theta
andphi
should be specified in radians and the sign convention with this gate varies. Here for example,fsim(- pi / 2, 0) == iswap()
.
- quimb.gen.operators.fsimg(theta, zeta, chi, gamma, phi, dtype=complex, **kwargs)[source]¶
The ‘fermionic simulation’ gate, with:
\(\theta\) is the iSWAP angle
\(\phi\) is the controlled-phase angle
\(\zeta, \chi, \gamma\) are single-qubit phase angles.
\[\begin{split}\mathrm{fsimg}(\theta, \zeta, \chi, \gamma, \phi) = \begin{bmatrix} 1 & 0 & 0 & 0\\ 0 & \exp(-i(\gamma +\zeta )) \cos(\theta) & -i \exp(-i(\gamma - \chi )) sin(\theta) & 0\\ 0 & -i \exp(-i(\gamma + \chi )) sin(\theta) & \exp(-i(\gamma - \zeta )) \cos(\theta) & 0\\ 0 & 0 & 0 & \exp(-i (\phi +2 \gamma)) \end{bmatrix}\end{split}\]See Equation 18 of https://arxiv.org/abs/2010.07965. Note that
theta
,phi
,zeta
,chi
,gamma
should be specified in radians and the sign convention with this gate varies. Here for example,fsimg(- pi / 2, 0, 0, 0,0) == iswap()
.
- quimb.gen.operators.ncontrolled_gate(ncontrol, gate, dtype=complex, sparse=False)[source]¶
Build an n-qubit controlled gate. The control qubits are the first
ncontrol
qubits.
- quimb.gen.operators.controlled(s, dtype=complex, sparse=False)[source]¶
Construct a controlled pauli gate for two qubits.
- quimb.gen.operators.ccX(dtype=complex, sparse=False)[source]¶
The double controlled X gate, or Toffoli gate.
- quimb.gen.operators.controlled_swap(dtype=complex, sparse=False)[source]¶
The controlled swap or Fredkin gate. The control qubit is the first qubit, if in state |1> a swap is performed on the last two qubits.
- quimb.gen.operators.hamiltonian_builder(fn)[source]¶
Wrap a function to perform some generic postprocessing and take the kwargs
stype
andsparse
. This assumes the core function always builds the hamiltonian in sparse form. The wrapper then:Checks if the operator is real and, if so, discards imaginary part if no explicity dtype was given
Converts the operator to dense or the correct sparse form
Makes the operator immutable so it can be safely cached
- quimb.gen.operators.ham_heis(n, j=1.0, b=0.0, cyclic=False, parallel=False, nthreads=None, ownership=None)[source]¶
Constructs the nearest neighbour 1d heisenberg spin-1/2 hamiltonian.
- Parameters:
n (int) – Number of spins.
j (float or tuple(float, float, float), optional) – Coupling constant(s), with convention that positive = antiferromagnetic. Can supply scalar for isotropic coupling or vector
(jx, jy, jz)
.b (float or tuple(float, float, float), optional) – Magnetic field, defaults to z-direction only if tuple not given.
cyclic (bool, optional) – Whether to couple the first and last spins.
sparse (bool, optional) – Whether to return the hamiltonian in sparse form.
stype (str, optional) – What format of sparse operator to return if
sparse
.parallel (bool, optional) – Whether to build the operator in parallel. By default will do this for n > 16.
nthreads (int optional) – How mny threads to use in parallel to build the operator.
ownership ((int, int), optional) – If given, which range of rows to generate.
kwargs – Supplied to
quimbify()
.
- Returns:
H – The Hamiltonian.
- Return type:
immutable operator
- quimb.gen.operators.ham_ising(n, jz=1.0, bx=1.0, **ham_opts)[source]¶
Generate the quantum transverse field ising model hamiltonian. This is a simple alias for
ham_heis()
with Z-interactions and an X-field.
- quimb.gen.operators.ham_XY(n, jxy, bz, **ham_opts)[source]¶
Generate the quantum transverse field XY model hamiltonian. This is a simple alias for
ham_heis()
with X- and Y-interactions and a Z-field.
- quimb.gen.operators.ham_XXZ(n, delta, jxy=1.0, **ham_opts)[source]¶
Generate the XXZ-model hamiltonian. This is a simple alias for
ham_heis()
with matched X- and Y-interactions anddelta
Z coupling.
- quimb.gen.operators.ham_j1j2(n, j1=1.0, j2=0.5, bz=0.0, cyclic=False, ownership=None)[source]¶
Generate the j1-j2 hamiltonian, i.e. next nearest neighbour interactions.
- Parameters:
n (int) – Number of spins.
j1 (float, optional) – Nearest neighbour coupling strength.
j2 (float, optional) – Next nearest neighbour coupling strength.
bz (float, optional) – B-field strength in z-direction.
cyclic (bool, optional) – Cyclic boundary conditions.
sparse (bool, optional) – Return hamiltonian as sparse-csr operator.
ownership ((int, int), optional) – If given, which range of rows to generate.
kwargs – Supplied to
quimbify()
.
- Returns:
H – The Hamiltonian.
- Return type:
immutable operator
- quimb.gen.operators.ham_mbl(n, dh, j=1.0, bz=0.0, cyclic=False, seed=None, dh_dist='s', dh_dim=1, beta=None, ownership=None)[source]¶
Constructs a heisenberg hamiltonian with isotropic coupling and random fields acting on each spin - the many-body localized (MBL) spin hamiltonian.
- Parameters:
n (int) – Number of spins.
dh (float or (float, float, float)) – Strength of random fields (stdev of gaussian distribution), can be scalar (isotropic noise) or 3-vector for (x, y, z) directions.
j (float or (float, float, float), optional) – Coupling strength, can be scalar (isotropic) or 3-vector.
bz (float, optional) – Global magnetic field (in z-direction).
cyclic (bool, optional) – Whether to use periodic boundary conditions.
seed (int, optional) – Number to seed random number generator with.
dh_dist ({'g', 's', 'qr'}, optional) –
Type of random distribution for the noise:
”s”: square, with bounds
(-dh, dh)
”g”: gaussian, with standard deviation
dh
”qp”: quasi periodic, with amplitude
dh
and ‘wavenumber’beta
so that the field at sitei
isdh * cos(2 * pi * beta * i + delta)
withdelta
a random offset between(0, 2 * pi)
, possibly seeded byseed
.
dh_dim ({1, 2, 3} or str, optional) – The number of dimensions the noise acts in, or string specifier like
'yz'
.beta (float, optional) – The wave number if
dh_dist='qr'
, defaults to the golden ratio``(5**0.5 - 1) / 2``.sparse (bool, optional) – Whether to construct the hamiltonian in sparse form.
stype ({'csr', 'csc', 'coo'}, optional) – The sparse format.
ownership ((int, int), optional) – If given, which range of rows to generate.
kwargs – Supplied to
quimbify()
.
- Returns:
H – The MBL hamiltonian for spin-1/2.
- Return type:
operator
See also
MPO_ham_mbl
- quimb.gen.operators.ham_heis_2D(n, m, j=1.0, bz=0.0, cyclic=False, parallel=False, ownership=None)[source]¶
Construct the 2D spin-1/2 heisenberg model hamiltonian:
\[\hat{H} = \sum_{<i, j>} J_X S^X_i S^X_j + J_Y S^Y_i S^Y_j + J_Z S^Z_i S^Z_j\]where the sum runs over pairs \(<i,j>\) on a 2D square lattice.
- Parameters:
n (int) – The number of rows.
m (int) – The number of columns.
j (float or (float, float, float), optional) – The coupling strength(s). Isotropic if scalar else if vector
(Jx, Jy, Jz) = j
.bz (float, optional) – The z direction magnetic field.
cyclic (bool, optional) – Whether to use periodic boundary conditions.
sparse (bool, optional) – Whether to construct the hamiltonian in sparse form.
stype ({'csr', 'csc', 'coo'}, optional) – The sparse format.
parallel (bool, optional) – Construct the hamiltonian in parallel. Faster but might use more memory.
ownership ((int, int), optional) – If given, which range of rows to generate.
kwargs – Supplied to
quimbify()
.
- Returns:
H – The hamiltonian.
- Return type:
operator
- quimb.gen.operators.uniq_perms(xs)[source]¶
Generate all the unique permutations of sequence
xs
.Examples
>>> list(uniq_perms('0011')) [('0', '0', '1', '1'), ('0', '1', '0', '1'), ('0', '1', '1', '0'), ('1', '0', '0', '1'), ('1', '0', '1', '0'), ('1', '1', '0', '0')]
- quimb.gen.operators.zspin_projector(n, sz=0, stype='csr', dtype=float)[source]¶
Construct the projector onto spin-z subpspaces.
- Parameters:
- Returns:
prj – The (non-square) projector onto the specified subspace(s). The subspace size
D
is given byn choose (n / 2 + s)
for eachs
specified insz
.- Return type:
immutable sparse operator, shape (2**n, D)
Examples
>>> zspin_projector(n=2, sz=0).toarray() array([[0., 0.], [1., 0.], [0., 1.], [0., 0.]]
Project a 9-spin Heisenberg-Hamiltonian into its spin-1/2 subspace:
>>> H = ham_heis(9, sparse=True) >>> H.shape (512, 512)
>>> P = zspin_projector(n=9, sz=1 / 2) >>> H0 = P.T @ H @ P >>> H0.shape (126, 126)
- quimb.gen.operators.create(n=2, **qu_opts)[source]¶
The creation operator acting on an n-level system.
- quimb.gen.operators.destroy(n=2, **qu_opts)[source]¶
The annihilation operator acting on an n-level system.
- quimb.gen.operators.ham_hubbard_hardcore(n, t=0.5, V=1.0, mu=1.0, cyclic=False, parallel=False, ownership=None)[source]¶
Generate the spinless fermion hopping hamiltonian.
- Parameters:
n (int) – The number of sites.
t (float, optional) – The hopping energy.
V (float, optional) – The interaction energy.
mu (float, optional) – The chemical potential - defaults to half-filling.
cyclic (bool, optional) – Whether to use periodic boundary conditions.
parallel (bool, optional) – Construct the hamiltonian in parallel. Faster but might use more memory.
ownership ((int, int), optional) – If given, which range of rows to generate.
kwargs – Supplied to
quimbify()
.
- Returns:
H – The hamiltonian.
- Return type:
operator