Introduction
Struqture is a Rust (struqture) and Python (struqture-py) library by HQS Quantum Simulations to represent quantum mechanical operators, Hamiltonians and open quantum systems. The library supports building spin objects, fermionic objects, bosonic objects and mixed system objects that contain arbitrary many spin, fermionic and bosonic subsystems.
Struqture has been developed to create and exchange definitions of operators, Hamiltonians and open systems. A special focus is the use as input to quantum computing simulation software.
To best support this use case, struqture has a number of design goals:
- Support for arbitrary spin, bosonic, fermionic and mixed systems
- Full serialisation support to json and other formats
- Preventing construction of unphysical objects by using well defined types for all objects in struqture
- Support of symbolic values in operators, Hamiltonians and open systems
Following these design goals, we prioritize using distinctive types to construct objects over a less verbose syntax. Similarly the support of symbolic expression leads to a trade-off in speed compared to an implementation using only floating point values. The symbolic expression support is achieved by using CalculatorComplex and CalculatorFloat values instead of complex and float values (respectively), which are imported from qoqo_calculator. Struqture is designed to also support the construction and (de)serialisation of large systems but for the use in numeric algorithms we recommend transforming Operators and Hamiltonians into a sparse matrix form.
This documentation is split into two parts. The first part covers the basic usage for spins, bosons, fermions and mixed systems. The second part covers the shared design patterns between spins, bosons, fermions and mixed systems. A real-world example is also included in.
Note: the package will be faster in Rust than Python, as Rust is a compiled language. This should only make a big difference, however, if you are performing hundreds of multiplication operations and a large amount of getter/setter calls.
Installation
Python
You can install struqture_py
from PyPi. For x86 Linux, Windows and macOS systems pre-built wheels are available.
On other platforms a local Rust toolchain is required to compile the Python source distribution.
pip install struqture-py
Rust
You can use struqture in your Rust project by adding
struqture = { version = "1.0.1" }
to your Cargo.toml file.
API Documentation
This user documentation is intended to give a high level overview of the design and usage of struqture. For a full list of the available data types and functions see the API-Documentaions of struqture and struqture-py.
Physical Types
In this part of the user documentation we show the basic usage of operators, Hamiltonians and open systems for spins, bosons, fermions and mixed systems. Stuqture is designed to use the same patterns to construct objects across all physical types. There is a large overlap in the user documentation between all physical types.
Spins
Building blocks
All spin objects in struqture are expressed based on products of either Pauli operators (X, Y, Z) or operators suited to express decoherence (X, iY, Z). The products are built by setting the operators acting on separate spins.
PauliProducts
PauliProducts are combinations of SingleSpinOperators on specific qubits. These are the SingleSpinOperators
, or Pauli matrices, that are available for PauliProducts:
-
I: identity matrix \[ I = \begin{pmatrix} 1 & 0\\ 0 & 1 \end{pmatrix} \]
-
X: Pauli x matrix \[ X = \begin{pmatrix} 0 & 1\\ 1 & 0 \end{pmatrix} \]
-
Y: Pauli y matrix \[ Y = \begin{pmatrix} 0 & -i\\ i & 0 \end{pmatrix} \]
-
Z: Pauli z matrix \[ Z = \begin{pmatrix} 1 & 0\\ 0 & -1 \end{pmatrix} \]
DecoherenceProducts
DecoherenceProducts are products of a decoherence operators acting on single spins. These SingleDecoherenceOperators
are almost identical to the SinglePauliOperators
with the exception of an additional \(i\) factor and are well suited to represent decoherence properties
- I: identity matrix \[ \begin{pmatrix} 1 & 0\\ 0 & 1 \end{pmatrix} \]
- X: Pauli X matrix \[ \begin{pmatrix} 0 & 1\\ 1 & 0 \end{pmatrix} \]
- iY: Pauli Y matrix multiplied by i \[ \begin{pmatrix} 0 & 1 \\ -1 & 0 \end{pmatrix} \]
- Z: Pauli z matrix \[ \begin{pmatrix} 1 & 0\\ 0 & -1 \end{pmatrix} \]
Examples
In Python the separate operators can be set via functions. In the python interface a PauliProduct can often be replaced by its unique string representation.
from struqture_py.spins import PauliProduct, DecoherenceProduct
# A product of a X acting on spin 0, a Y acting on spin 3 and a Z acting on spin 20
pp = PauliProduct().x(0).y(3).z(20)
# Often equivalent the string representation
pp_string = str(pp)
# A product of a X acting on spin 0, a iY acting on spin 3 and a Z acting on spin 20
dp = DecoherenceProduct().x(0).iy(3).z(20)
# Often equivalent the string representation
dp_string = str(dp)
In Rust the user can also import enums for the operators acting on single spins. In Rust the equivalent string representation cannot be used in function and method arguments.
use struqture::prelude::*;
use struqture::spins::{
DecoherenceProduct, PauliProduct, SingleDecoherenceOperator, SingleSpinOperator,
};
// A product of a X acting on spin 0, a Y acting on spin 3 and a Z acting on spin 20
let pp = PauliProduct::new().x(0).y(3).z(20);
// Constructing with SingleSpinOperator
let pp_equivalent = PauliProduct::new()
.set_pauli(0, SingleSpinOperator::X)
.set_pauli(3, SingleSpinOperator::Y)
.set_pauli(20, SingleSpinOperator::Z);
// A product of a X acting on spin 0, a Y acting on spin 3 and a Z acting on spin 20
let dp = DecoherenceProduct::new().x(0).iy(3).z(20);
// Constructing with SingleSpinOperator
let dp_equivalent = DecoherenceProduct::new()
.set_pauli(0, SingleDecoherenceOperator::X)
.set_pauli(3, SingleDecoherenceOperator::IY)
.set_pauli(20, SingleDecoherenceOperator::Z);
Operators and Hamiltonians
A good example how complex objects are constructed from operator products are SpinOperators
and SpinHamiltonians
(for more information, see also).
These SpinOperators
and SpinHamiltonians
represent operators or Hamiltonians such as:
\[
\hat{O} = \sum_{j} \alpha_j \prod_{k=0}^N \sigma_{j, k} \\
\sigma_{j, k} \in \{ X_k, Y_k, Z_k, I_k \}
\]
where the \(\sigma_{j, k}\) are SinglePauliOperators
.
From a programming perspective the operators and Hamiltonians are HashMaps or Dictionaries with the PauliProducts
as keys and the coefficients \(\alpha_j\) as values.
In struqture we distinguish between spin operators and Hamiltonians to avoid introducing unphysical behaviour by accident.
While both are sums over PauliProducts, Hamiltonians are guaranteed to be hermitian. In a spin Hamiltonian , this means that the prefactor of each PauliProduct
has to be real.
Examples
Here is an example of how to build a SpinOperator
and a SpinHamiltonian
, in Rust:
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::spins::{PauliProduct, SpinOperator, SpinHamiltonian};
// Building the term sigma^x_0 * sigma^z_2: sigma_x acting on qubit 0
// and sigma_z acting on qubit 2
let pp = PauliProduct::new().x(0).z(2);
// O = (1 + 1.5 * i) * sigma^x_0 * sigma^z_2
let mut operator = SpinOperator::new();
operator.add_operator_product(pp.clone(), CalculatorComplex::new(1.0, 1.5)).unwrap();
assert_eq!(operator.get(&pp), &CalculatorComplex::new(1.0, 1.5));
println!("{}", operator);
// Or when overwriting the previous value
let mut operator = SpinOperator::new();
operator.set(pp.clone(), CalculatorComplex::new(1.0, 1.5)).unwrap();
println!("{}", operator);
// A complex entry is not valid for a SpinHamiltonian
let mut hamiltonian = SpinHamiltonian::new();
// This would fail
hamiltonian.add_operator_product(pp, CalculatorComplex::new(1.0, 1.5)).unwrap();
// This is possible
hamiltonian.add_operator_product(pp, 1.0.into()).unwrap();
println!("{}", hamiltonian);
In python, we need to use a SpinSystem
and SpinHamiltonianSystem
instead of anSpinOperator
and SpinHamiltonian
. See next section for more details.
Systems and HamiltonianSystems
Following the intention to avoid unphysical behaviour, SpinSystems and SpinHamiltonianSystems are wrappers around SpinOperators and SpinHamiltonians that allow to explicitly set the number of spins of the systems. When setting or adding a PauliProduct to the systems, it is guaranteed that the spin indices involved cannot exceed the number of spins in the system. Note that the user can decide to explicitly set the number of spins to be variable.
Examples
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::spins::{PauliProduct, SpinSystem};
let mut system = SpinSystem::new(Some(3));
// This will work
let pp = PauliProduct::new().x(0).z(2);
system
.add_operator_product(pp, CalculatorComplex::new(1.0, 1.5))
.unwrap();
println!("{}", system);
// This will not work, as the spin index of the PauliProduct is larger than
// the number of the spins in the system (the spin with the smallest index is 0).
let pp_error = PauliProduct::new().z(3);
let error = system.add_operator_product(pp_error, CalculatorComplex::new(1.0, 1.5));
println!("{:?}", error);
// This will work because we leave the number of spins dynamic
let mut system = SpinSystem::new(None);
system
.add_operator_product(PauliProduct::new().z(3), CalculatorComplex::new(1.0, 1.5))
.unwrap();
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import spins
system = spins.SpinSystem(3)
# This will work
pp = spins.PauliProduct().x(0).z(2)
system.add_operator_product(pp, CalculatorComplex.from_pair(1.0, 1.5))
print(system)
# This will not work, as the spin index of the PauliProduct is larger
# than the number of the spins in the system (the spin with the smallest index is 0).
pp_error = spins.PauliProduct().z(3)
value = CalculatorComplex.from_pair(1.0, 1.5)
# system.add_operator_product(pp_error, value) # Uncomment me!
# This will work because we leave the number of spins dynamic
system = spins.SpinSystem()
system.add_operator_product(spins.PauliProduct().z(3), 1.0)
Noise operators and systems
We describe decoherence by representing it with the Lindblad equation. The Lindblad equation is a master equation determining the time evolution of the density matrix. It is given by \[ \dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right) \] with the rate matrix \(\Gamma_{j,k}\) and the Lindblad operator \(L_{j}\).
To describe spin noise we use the Lindblad equation with \(\hat{H}=0\).
Therefore, to describe the pure noise part of the Lindblad equation one needs the rate matrix in a well defined basis of Lindblad operators.
We use DecoherenceProducts
as the operator basis.
The rate matrix and with it the Lindblad noise model is saved as a sum over pairs of DecoherenceProducts
, giving the operators acting from the left and right on the density matrix.
In programming terms the object SpinLindbladNoiseOperator
is given by a HashMap or Dictionary with the tuple (DecoherenceProduct
, DecoherenceProduct
) as keys and the entries in the rate matrix as values.
Similarly to SpinOperators, SpinLindbladNoiseOperators have a system equivalent: SpinLindbladNoiseSystem
, with a number of involved spins defined by the user. For more information on these, see the noise container chapter.
Examples
Here, we add the terms \( L_0 = \sigma_0^{x} \sigma_2^{z} \) and \( L_1 = \sigma_0^{x} \sigma_2^{z} \) with coefficient 1.0: \( 1.0 \left( L_0 \rho L_1^{\dagger} - \frac{1}{2} \{ L_1^{\dagger} L_0, \rho \} \right) \)
use struqture::prelude::*;
use struqture::spins::{DecoherenceProduct, SpinLindbladNoiseSystem};
// Constructing the system and product to be added to it
let mut system = SpinLindbladNoiseSystem::new(Some(3));
let dp = DecoherenceProduct::new().x(0).z(2);
// Adding in the 0X2Z term
system.add_operator_product((dp.clone(), dp.clone()), 1.0.into()).unwrap();
// Checking our system
assert_eq!(system.get(&(dp.clone(), dp)), &CalculatorComplex::new(1.0, 0.0));
println!("{}", system);
The equivalent code in python:
from struqture_py import spins
# Constructing the system and product to be added to it
system = spins.SpinLindbladNoiseSystem(3)
dp = spins.DecoherenceProduct().x(0).z(2)
# Adding in the 0X2Z term
system.add_operator_product((dp, dp), 1.0+1.5*1j)
print(system)
# In python we can also use the string representation
system = spins.SpinLindbladNoiseSystem(3)
system.add_operator_product(("0X2Z", "0X2Z"), 1.0+1.5*1j)
print(system)
Open systems
Physically open systems are quantum systems coupled to an environment that can often be described using Lindblad type of noise.
The Lindblad master equation is given by
\[
\dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right)
\]
In struqture they are composed of a hamiltonian (SpinHamiltonianSystem
) and noise (SpinLindbladNoiseSystem
). They have different ways to set terms in Rust and Python:
Examples
use qoqo_calculator::{CalculatorComplex, CalculatorFloat};
use struqture::prelude::*;
use struqture::spins::{DecoherenceProduct, PauliProduct, SpinLindbladOpenSystem};
let mut open_system = SpinLindbladOpenSystem::new(Some(3));
let pp = PauliProduct::new().z(1);
let dp = DecoherenceProduct::new().x(0).z(2);
// Add the Z_1 term into the system part of the open system
let system = open_system.system_mut();
system.add_operator_product(pp, CalculatorFloat::from(2.0)).unwrap();
// Add the X_0 Z_2 term into the noise part of the open system
let noise = open_system.noise_mut();
noise
.add_operator_product((dp.clone(), dp), CalculatorComplex::new(1.0, 0.0))
.unwrap();
println!("{}", open_system);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex, CalculatorFloat
from struqture_py import spins
open_system = spins.SpinLindbladOpenSystem(3)
pp = spins.PauliProduct().z(1)
dp = spins.DecoherenceProduct().x(0).z(2)
# Add the Z_1 term into the system part of the open system
open_system.system_add_operator_product(pp, CalculatorFloat(2.0))
# Add the X_0 Z_2 term into the noise part of the open system
open_system.noise_add_operator_product(
(dp, dp), CalculatorComplex.from_pair(0.0, 1.0))
print(open_system)
Matrix representation: spin objects only
All spin-objects can be converted into sparse matrices with the following convention.
If \(M_2\) corresponds to the matrix acting on spin 2 and \(M_1\) corresponds to the matrix acting on spin 1 the total matrix \(M\) acting on spins 0 to 2 is given by
\[
M = M_2 \otimes M_1 \otimes \mathbb{1}
\]
For an \(N\)-spin system an operator acts on the \(2^N\) dimensional space of state vectors.
A superoperator operates on the \(4^N\) dimensional space of flattened density-matrices.
struqture uses the convention that density matrices are flattened in row-major order
\[
\rho = \begin{pmatrix} a & b \\ c & d \end{pmatrix} => \vec{\rho} = \begin{pmatrix} a \\ b \\ c \\ d \end{pmatrix}
\]
For noiseless objects (SpinSystem
, SpinHamiltonianSystem
), sparse operators and sparse superoperators can be constructed, as we can represent the system as a wavefunction. For systems with noise (SpinLindbladNoiseSystem
, SpinLindbladOpenSystem
), however, we can only represent them as density matrices and can therefore only construct sparse superoperators.
Note that the matrix representation functionality exists only for spin objects, and can't be generated for bosonic, fermionic or mixed system objects.
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::spins::{DecoherenceProduct, SpinLindbladNoiseSystem};
let mut system = SpinLindbladNoiseSystem::new(Some(3));
let dp = DecoherenceProduct::new().x(0).z(2);
system
.add_operator_product((dp.clone(), dp), CalculatorComplex::new(1.0, 0.0))
.unwrap();
// Here we have a noise system, so we can only construct a superoperator
let matrix = system.sparse_matrix_superoperator(Some(3)).unwrap();
println!("{:?}", matrix);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import spins
from scipy.sparse import coo_matrix
system = spins.SpinLindbladNoiseSystem(3)
dp = spins.DecoherenceProduct().x(0).z(2)
system.add_operator_product((dp, dp), CalculatorComplex.from_pair(1.0, 1.5))
# Using the `sparse_matrix_superoperator_coo` function, you can also
# return the information in scipy coo_matrix form, which can be directly fed in:
python_coo = coo_matrix(system.sparse_matrix_superoperator_coo())
print(python_coo.todense())
Fermions
Building blocks
All fermionic objects in struqture are expressed based on products of fermionic creation and annihilation operators, which respect fermionic anti-commutation relations \[ \lbrace c_k^{\dagger}, c_j^{\dagger} \rbrace = 0, \\ \lbrace c_k, c_j \rbrace = 0, \\ \lbrace c_k, c_j^{\dagger} \rbrace = \delta_{k, j}. \]
FermionProducts
FermionProducts are simple combinations of fermionic creation and annihilation operators.
HermitianFermionProducts
HermitianFermionProducts are the hermitian equivalent of FermionProducts. This means that even though they are constructed the same (see the next section, Examples
), they internally store both that term and its hermitian conjugate. For instance, given the term \(c^{\dagger}_0 c_1 c_2\), a FermionProduct would represent \(c^{\dagger}_0 c_1 c_2\) while a HermitianFermionProduct would represent \(c^{\dagger}_0 c_1 c_2 + c^{\dagger}_2 c^{\dagger}_1 c_0\).
Examples
In both Python and Rust, the operator product is constructed by passing an array or a list of integers to represent the creation indices, and an array or a list of integers to represent the annihilation indices.
Note: (Hermitian)FermionProducts can only been created from the correct ordering of indices (the wrong sequence will return an error) but we have the create_valid_pair
function to create a valid Product from arbitrary sequences of operators which also transforms an index value according to the anti-commutation and hermitian conjugation rules.
from struqture_py.fermions import FermionProduct, HermitianFermionProduct
from qoqo_calculator_pyo3 import CalculatorComplex
# A product of a creation operator acting on fermionic mode 0 and an
# annihilation operator acting on fermionic mode 20
fp = FermionProduct([0], [20])
# Building the term c^{\dagger}_1 * c^{\dagger}_3 * c_0
fp = FermionProduct.create_valid_pair(
[3, 1], [0], CalculatorComplex.from_pair(1.0, 0.0))
# A product of a creation operator acting on fermionic mode 0 and an annihilation
# operator acting on fermionic mode 20, as well as a creation operator acting on
# fermionic mode 20 and an annihilation operator acting on fermionic mode 0
hfp = HermitianFermionProduct([0], [20])
# Building the term c^{\dagger}_0 * c^{\dagger}_3 * c_0 + c^{\dagger}_0 * c_3 * c_0
hfp = HermitianFermionProduct.create_valid_pair(
[3, 0], [0], CalculatorComplex.from_pair(1.0, 0.0))
In Rust the equivalent string representation cannot be used in function and method arguments.
use qoqo_calculator::CalculatorComplex;
use struqture::fermions::{FermionProduct, HermitianFermionProduct};
use struqture::prelude::*;
// Building the term c^{\dagger}_0 c_20
let fp_0 = FermionProduct::new([0], [20]).unwrap();
// Building the term c^{\dagger}_1 * c^{\dagger}_3 * c_0
let (fp_1, coeff) = FermionProduct::create_valid_pair(
[3, 1], [0], CalculatorComplex::from(1.0)).unwrap();
// A product of a creation operator acting on fermionic mode 0 and an annihilation
// operator acting on fermionic mode 20, as well as a creation operator acting on
// fermionic mode 20 and an annihilation operator acting on fermionic mode 0
let fp_0 = HermitianFermionProduct::new([0], [20]).unwrap();
// Building the term c^{\dagger}_0 * c^{\dagger}_3 * c_0 + c^{\dagger}_0 * c_3 * c_0
let (fp_1, coeff) = HermitianFermionProduct::create_valid_pair(
[3, 0], [0], CalculatorComplex::from(1.0)).unwrap();
Operators and Hamiltonians
Complex objects are constructed from operator products are FermionOperators
and FermionHamiltonians
(for more information, see also).
These FermionOperators
and FermionHamiltonians
represent operators or Hamiltonians such as:
\[ \hat{O} = \sum_{j=0}^N \alpha_j \left( \prod_{k=0}^N f(j, k) \right) \left( \prod_{l=0}^N g(j, l) \right) \]
with
\[ f(j, k) = \begin{cases} c_k^{\dagger} \\ \mathbb{1} \end{cases} , \]
\[ g(j, l) = \begin{cases} c_l \\ \mathbb{1} \end{cases} , \]
and
\(c^{\dagger}\) the fermionionic creation operator, \(c\) the fermionionic annihilation operator
\[ \lbrace c_k^{\dagger}, c_j^{\dagger} \rbrace = 0, \\
\lbrace c_k, c_j \rbrace = 0, \\
\lbrace c_k^{\dagger}, c_j \rbrace = \delta_{k, j}. \]
For instance, \(c^{\dagger}_0 c^{\dagger}_1 c_1\) is a term with a \(c^{\dagger}\) term acting on 0, and both a \(c^{\dagger}\) term and a \(c\) term acting on 1.
From a programming perspective the operators and Hamiltonians are HashMaps or Dictionaries with FermionProducts
or HermitianFermionProducts
(respectively) as keys and the coefficients \(\alpha_j\) as values.
In struqture we distinguish between fermionic operators and Hamiltonians to avoid introducing unphysical behaviour by accident.
While both are sums over normal ordered fermionic products (stored as HashMaps of products with a complex prefactor), Hamiltonians are guaranteed to be hermitian. In a fermionic Hamiltonian , this means that the sums of products are sums of hermitian fermionic products (we have not only the \(c^{\dagger}c\) terms but also their hermitian conjugate) and the on-diagonal terms are required to have real prefactors.
In the HermitianFermionProducts
, we only explicitly store one part of the hermitian fermionic product, and we have chosen to store the one which has the smallest index of the creators that is smaller than the smallest index of the annihilators.
Examples
Here is an example of how to build a product and using it to build an operator, in Rust:
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::fermions::{
FermionProduct, FermionOperator, HermitianFermionProduct, FermionHamiltonian
};
// Building the term c^{\dagger}_1 * c^{\dagger}_2 * c_0 * c_1
let fp = FermionProduct::new([1, 2], [0, 1]).unwrap();
// O = (1 + 1.5 * i) * c^{\dagger}_1 * c^{\dagger}_2 * c_0 * c_1
let mut operator = FermionOperator::new();
operator.add_operator_product(fp.clone(), CalculatorComplex::new(1.0, 1.5)).unwrap();
println!("{}", operator);
// Or when overwriting the previous value
let mut operator = FermionOperator::new();
operator.set(fp.clone(), CalculatorComplex::new(1.0, 1.5)).unwrap();
println!("{}", operator);
// A FermionProduct entry is not valid for a FermionHamiltonian
let mut hamiltonian = FermionHamiltonian::new();
// This would fail, as it uses HermitianFermionProducts, not FermionProducts
hamiltonian.add_operator_product(fp, CalculatorComplex::new(1.0, 1.5)).unwrap();
// This is possible
let hfp = HermitianFermionProduct::new([0, 2], [0, 1]).unwrap();
hamiltonian.add_operator_product(hfp, CalculatorComplex::new(1.5, 0.0)).unwrap();
println!("{}", hamiltonian);
In python, we need to use a FermionSystem
and FermionHamiltonianSystem
instead of a FermionOperator
and FermionHamiltonian
. See next section for more details.
Systems and HamiltonianSystems
Following the intention to avoid unphysical behaviour, FermionSystems and FermionHamiltonianSystems are wrappers around FermionOperators and FermionHamiltonians that allow to explicitly set the number of spins of the systems. When setting or adding a FermionProduct/HermitianFermionProduct to the systems, it is guaranteed that the fermionic indices involved cannot exceed the number of fermionic modes in the system. Note that the user can decide to explicitly set the number of fermionic modes to be variable.
Examples
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::fermions::{HermitianFermionProduct, FermionHamiltonianSystem};
let mut system = FermionHamiltonianSystem::new(Some(3));
// This will work
let hfp = HermitianFermionProduct::new([0, 1], [0, 2]).unwrap();
system.add_operator_product(hfp, CalculatorComplex::new(1.0, 1.5)).unwrap();
println!("{}", system);
// This will not work, as the fermionic index of the HermitianFermionProduct is larger
// than the number of the fermionic modes in the system (the fermionic mode with the
// smallest index is 0).
let hfp_error = HermitianFermionProduct::new([3], [3]).unwrap();
let error = system.add_operator_product(hfp_error, CalculatorComplex::new(1.0, 1.5));
println!("{:?}", error);
// This will work because we leave the number of spins dynamic
let hbf = HermitianFermionProduct::new([0, 1], [0, 2]).unwrap();
let mut system = FermionHamiltonianSystem::new(None);
system.add_operator_product(hbf, CalculatorComplex::new(1.0, 1.5)).unwrap();
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import fermions
system = fermions.FermionHamiltonianSystem(3)
# This will work
hfp = fermions.HermitianFermionProduct([0, 1], [0, 2])
system.add_operator_product(hfp, CalculatorComplex.from_pair(1.0, 1.5))
print(system)
# This will not work, as the fermioncic index of the HermitianFermionProduct is larger
# than the number of the fermionic modes in the system (the fermionic mode with the
# smallest index is 0).
hfp_error = fermions.HermitianFermionProduct([3], [3])
value = CalculatorComplex.from_pair(1.0, 1.5)
# system.add_operator_product(hfp_error, value) # Uncomment me!
# This will work because we leave the number of spins dynamic
system = fermions.FermionHamiltonianSystem()
hfp = fermions.HermitianFermionProduct([3], [3])
system.add_operator_product(hfp, CalculatorComplex.from_pair(1.0, 0.0))
Noise operators and systems
We describe decoherence by representing it with the Lindblad equation. The Lindblad equation is a master equation determining the time evolution of the density matrix. It is given by \[ \dot{\rho} = \mathcal{L}(\rho) = -i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right) \] with the rate matrix \(\Gamma_{j,k}\) and the Lindblad operator \(L_{j}\).
To describe fermionic noise we use the Lindblad equation with \(\hat{H}=0\).
Therefore, to describe the pure noise part of the Lindblad equation one needs the rate matrix in a well defined basis of Lindblad operators.
We use FermionProducts
as the operator basis.
The rate matrix and with it the Lindblad noise model is saved as a sum over pairs of FermionProducts
, giving the operators acting from the left and right on the density matrix.
In programming terms the object FermionLindbladNoiseOperator
is given by a HashMap or Dictionary with the tuple (FermionProduct
, FermionProduct
) as keys and the entries in the rate matrix as values.
Similarly to FermionOperators, FermionLindbladNoiseOperators have a system equivalent: FermionLindbladNoiseSystem
, with a number of involved fermionic modes defined by the user. For more information on these, see the noise container chapter.
Examples
Here, we add the terms \(L_0 = c^{\dagger}_0 c_0\) and \(L_1 = c^{\dagger}_0 c_0\) with coefficient 1.0: \( 1.0 \left( L_0 \rho L_1^{\dagger} - \frac{1}{2} \{ L_1^{\dagger} L_0, \rho \} \right) \)
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::fermions::{FermionProduct, FermionLindbladNoiseSystem};
// Setting up the system and the product we want to add to it
let mut system = FermionLindbladNoiseSystem::new(Some(3));
let fp = FermionProduct::new([0], [0]).unwrap();
// Adding the product to the system
system
.add_operator_product(
(fp.clone(), fp.clone()),
CalculatorComplex::new(1.0, 0.0)
).unwrap();
assert_eq!(system.get(&(fp.clone(), fp)), &CalculatorComplex::new(1.0, 0.0));
println!("{}", system);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import fermions
# Setting up the system and the product we want to add to it
system = fermions.FermionLindbladNoiseSystem(3)
fp = fermions.FermionProduct([0], [0])
# Adding the product to the system
system.add_operator_product((fp, fp), CalculatorComplex.from_pair(1.0, 1.5))
print(system)
# In python we can also use the string representation
system = fermions.FermionLindbladNoiseSystem(3)
system.add_operator_product((str(fp), str(fp)), 1.0+1.5*1j)
print(system)
Open systems
Physically open systems are quantum systems coupled to an environment that can often be described using Lindblad type of noise.
The Lindblad master equation is given by
\[
\dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right)
\]
In struqture they are composed of a Hamiltonian (FermionHamiltonianSystem
) and noise (FermionLindbladNoiseSystem
). They have different ways to set terms in Rust and Python:
Examples
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::fermions::{
FermionProduct, HermitianFermionProduct, FermionLindbladOpenSystem
};
let mut open_system = FermionLindbladOpenSystem::new(Some(3));
let hfp = HermitianFermionProduct::new([0, 1], [0, 2]).unwrap();
let fp = FermionProduct::new([0], [0]).unwrap();
// Adding the c^{\dagger}_0 c^{\dagger}_1 c_0 c_2 term to the system part of the open system
let system = open_system.system_mut();
system.add_operator_product(hfp, CalculatorComplex::new(2.0, 0.0)).unwrap();
// Adding the c^{\dagger}_0 c_0 part to the noise part of the open system
let noise = open_system.noise_mut();
noise
.add_operator_product(
(fp.clone(), fp), CalculatorComplex::new(1.0, 0.0)
).unwrap();
println!("{}", open_system);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import fermions
open_system = fermions.FermionLindbladOpenSystem(3)
hfp = fermions.HermitianFermionProduct([0, 1], [0, 2])
fp = fermions.FermionProduct([0], [0])
# Adding the c_0^dag c_1^dag c_0 c_2 term to the system part of the open system
open_system.system_add_operator_product(hfp, CalculatorComplex.from_pair(2.0, 0.0))
# Adding the c_0^dag c_0 part to the noise part of the open system
open_system.noise_add_operator_product(
(fp, fp), CalculatorComplex.from_pair(0.0, 1.0))
print(open_system)
Bosons
Building blocks
All bosonic objects in struqture are expressed based on products of bosonic creation and annihilation operators, which respect bosonic commutation relations \[ \lbrack c_k^{\dagger}, c_j^{\dagger} \rbrack = 0, \\ \lbrack c_k, c_j \rbrack = 0, \\ \lbrack c_k, c_j^{\dagger} \rbrack = \delta_{k, j}. \]
BosonProducts
BosonProducts are simple combinations of bosonic creation and annihilation operators.
HermitianBosonProducts
HermitianBosonProducts are the hermitian equivalent of BosonProducts. This means that even though they are constructed the same (see the next section, Examples
), they internally store both that term and its hermitian conjugate. For instance, given the term \(c^{\dagger}_0 c_1 c_2\), a BosonProduct would represent \(c^{\dagger}_0 c_1 c_2\) while a HermitianBosonProduct would represent \(c^{\dagger}_0 c_1 c_2 + c^{\dagger}_2 c^{\dagger}_1 c_0\).
Examples
In both Python and Rust, the operator product is constructed by passing an array or a list of integers to represent the creation indices, and an array or a list of integers to represent the annihilation indices.
Note: (Hermitian)BosonProducts can only been created from the correct ordering of indices (the wrong sequence will return an error) but we have the create_valid_pair
function to create a valid Product from arbitrary sequences of operators which also transforms an index value according to the commutation and hermitian conjugation rules.
from struqture_py.bosons import BosonProduct, HermitianBosonProduct
from qoqo_calculator_pyo3 import CalculatorComplex
# A product of a creation operator acting on bosonic mode 0 and an annihilation operator
# acting on bosonic mode 20
bp = BosonProduct([0], [20])
# Building the term c^{\dagger}_1 * c^{\dagger}_3 * c_0
bp = BosonProduct.create_valid_pair([3, 1], [0], CalculatorComplex.from_pair(1.0, 0.0))
# A product of a creation operator acting on bosonic mode 0 and an annihilation
# operator acting on bosonic mode 20, as well as a creation operator acting on
# bosonic mode 20 and an annihilation operator acting on bosonic mode 0
hbp = HermitianBosonProduct([0], [20])
# Building the term c^{\dagger}_0 * c^{\dagger}_3 * c_0 + c^{\dagger}_0 * c_3 * c_0
hbp = HermitianBosonProduct.create_valid_pair(
[3, 0], [0], CalculatorComplex.from_pair(1.0, 0.0))
In Rust the equivalent string representation cannot be used in function and method arguments.
use struqture::bosons::{BosonProduct, HermitianBosonProduct};
use struqture::ModeIndex;
use qoqo_calculator::CalculatorComplex;
// Building the term c^{\dagger}_0 c_20
let bp_0 = BosonProduct::new([0], [20]).unwrap();
// Building the term c^{\dagger}_1 * c^{\dagger}_3 * c_0
let (bp_1, coeff) = BosonProduct::create_valid_pair(
[3, 1], [0], CalculatorComplex::from(1.0)
).unwrap();
// A product of a creation operator acting on bosonic mode 0 and an annihilation operator
// acting on bosonic mode 20, as well as a creation operator acting on bosonic mode 20
// and an annihilation operator acting on bosonic mode 0
let bp_0 = HermitianBosonProduct::new([0], [20]).unwrap();
// Building the term c^{\dagger}_0 * c^{\dagger}_3 * c_0 + c^{\dagger}_0 * c_3 * c_0
let (bp_1, coeff) = HermitianBosonProduct::create_valid_pair(
[3, 0], [0], CalculatorComplex::from(1.0)
).unwrap();
Operators and Hamiltonians
Complex objects are constructed from operator products are BosonOperators
and BosonHamiltonians
(for more information, see also).
These BosonOperators
and BosonHamiltonians
represent operators or Hamiltonians such as:
\[ \hat{O} = \sum_{j=0}^N \alpha_j \left( \prod_{k=0}^N f(j, k) \right) \left( \prod_{l=0}^N g(j, l) \right) \]
with
\[ f(j, k) = \begin{cases} c_k^{\dagger} \\ \mathbb{1} \end{cases} , \]
\[ g(j, l) = \begin{cases} c_l \\ \mathbb{1} \end{cases} , \]
and
\(c^{\dagger}\) the bosonic creation operator, \(c\) the bosonic annihilation operator
\[ \lbrack c_k^{\dagger}, c_j^{\dagger} \rbrack = 0, \\
\lbrack c_k, c_j \rbrack = 0, \\
\lbrack c_k^{\dagger}, c_j \rbrack = \delta_{k, j}. \]
From a programming perspective the operators and Hamiltonians are HashMaps or Dictionaries with BosonProducts
or HermitianBosonProducts
(respectively) as keys and the coefficients \(\alpha_j\) as values.
In struqture we distinguish between bosonic operators and Hamiltonians to avoid introducing unphysical behaviour by accident.
While both are sums over normal ordered bosonic products (stored as HashMaps of products with a complex prefactor), Hamiltonians are guaranteed to be hermitian. In a bosonic Hamiltonian , this means that the sums of products are sums of hermitian bosonic products (we have not only the \(c^{\dagger}c\) terms but also their hermitian conjugate) and the on-diagonal terms are required to have real prefactors.
In the HermitianBosonProducts
, we only explicitly store one part of the hermitian bosonic product, and we have chosen to store the one which has the smallest index of the creators that is smaller than the smallest index of the annihilators.
Examples
Here is an example of how to build a BosonOperator
and a BosonHamiltonian
, in Rust:
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::bosons::{
BosonProduct, BosonOperator, HermitianBosonProduct, BosonHamiltonian
};
// Building the term c^{\dagger}_1 * c^{\dagger}_2 * c_0 * c_1
let bp = BosonProduct::new([1, 2], [0, 1]).unwrap();
// O = (1 + 1.5 * i) * c^{\dagger}_1 * c^{\dagger}_2 * c_0 * c_1
let mut operator = BosonOperator::new();
operator.add_operator_product(bp.clone(), CalculatorComplex::new(1.0, 1.5)).unwrap();
println!("{}", operator);
// Or when overwriting the previous value
let mut operator = BosonOperator::new();
operator.set(bp.clone(), CalculatorComplex::new(1.0, 1.5)).unwrap();
println!("{}", operator);
// A BosonProduct entry is not valid for a BosonHamiltonian
let mut hamiltonian = BosonHamiltonian::new();
// This would fail, as it uses HermitianBosonProducts, not BosonProducts
hamiltonian.add_operator_product(bp, CalculatorComplex::new(1.0, 1.5)).unwrap();
// This is possible
let hbp = HermitianBosonProduct::new([0, 2], [0, 1]).unwrap();
hamiltonian.add_operator_product(hbp, CalculatorComplex::new(1.5, 0.0)).unwrap();
println!("{}", hamiltonian);
In python, we need to use a BosonSystem
and BosonHamiltonianSystem
instead of a BosonOperator
and BosonHamiltonian
. See next section for more details.
Systems and HamiltonianSystems
Following the intention to avoid unphysical behaviour, BosonSystems and BosonHamiltonianSystems are wrappers around BosonOperators and BosonHamiltonians that allow to explicitly set the number of spins of the systems. When setting or adding a BosonProduct/HermitianBosonProduct to the systems, it is guaranteed that the bosonic indices involved cannot exceed the number of bosonic modes in the system. Note that the user can decide to explicitly set the number of bosonic modes to be variable.
Examples
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::bosons::{HermitianBosonProduct, BosonHamiltonianSystem};
let mut system = BosonHamiltonianSystem::new(Some(3));
// This will work
let hbp = HermitianBosonProduct::new([0, 1], [0, 2]).unwrap();
system.add_operator_product(hbp, CalculatorComplex::new(1.0, 1.5)).unwrap();
println!("{}", system);
// This will not work, as the bosonic index of the HermitianBosonProduct is larger than
// the number of the bosonic modes in the system (the bosonic mode with the
// smallest index is 0).
let hbp_error = HermitianBosonProduct::new([3], [3]).unwrap();
let error = system.add_operator_product(hbp_error, CalculatorComplex::new(1.0, 1.5));
println!("{:?}", error);
// This will work because we leave the number of spins dynamic
let hbp = HermitianBosonProduct::new([0, 1], [0, 2]).unwrap();
let mut system = BosonHamiltonianSystem::new(None);
system.add_operator_product(hbp, CalculatorComplex::new(1.0, 1.5)).unwrap();
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import bosons
system = bosons.BosonHamiltonianSystem(3)
# This will work
hbp = bosons.HermitianBosonProduct([0, 1], [0, 2])
system.add_operator_product(hbp, CalculatorComplex.from_pair(1.0, 1.5))
print(system)
# This will not work, as the bosoncic index of the HermitianBosonProduct
# is larger than the number of the bosonic modes in the system (the bosonic
# mode with the smallest index is 0).
hbp_error = bosons.HermitianBosonProduct([3], [3])
value = CalculatorComplex.from_pair(1.0, 1.5)
# system.add_operator_product(hbp_error, value) # Uncomment me!
# This will work because we leave the number of spins dynamic
system = bosons.BosonHamiltonianSystem()
hbp = bosons.HermitianBosonProduct([3], [3])
system.add_operator_product(hbp, CalculatorComplex.from_pair(1.0, 0.0))
Noise operators and systems
We describe decoherence by representing it with the Lindblad equation. The Lindblad equation is a master equation determining the time evolution of the density matrix. It is given by \[ \dot{\rho} = \mathcal{L}(\rho) = -i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right) \] with the rate matrix \(\Gamma_{j,k}\) and the Lindblad operator \(L_{j}\).
To describe bosonic noise we use the Lindblad equation with \(\hat{H}=0\).
Therefore, to describe the pure noise part of the Lindblad equation one needs the rate matrix in a well defined basis of Lindblad operators.
We use BosonProducts
as the operator basis.
The rate matrix and with it the Lindblad noise model is saved as a sum over pairs of BosonProducts
, giving the operators acting from the left and right on the density matrix.
In programming terms the object BosonLindbladNoiseOperator
is given by a HashMap or Dictionary with the tuple (BosonProduct
, BosonProduct
) as keys and the entries in the rate matrix as values.
Similarly to BosonOperators, BosonLindbladNoiseOperators have a system equivalent: BosonLindbladNoiseSystem
, with a number of involved bosonic modes defined by the user. For more information on these, see the noise container chapter.
Examples
Here, we add the terms \(L_0 = c^{\dagger}_0 c_0\) and \(L_1 = c^{\dagger}_0 c_0\) with coefficient 1.0: \( 1.0 \left( L_0 \rho L_1^{\dagger} - \frac{1}{2} \{ L_1^{\dagger} L_0, \rho \} \right) \)
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::bosons::{BosonProduct, BosonLindbladNoiseSystem};
// Setting up the system and the product we want to add to it
let mut system = BosonLindbladNoiseSystem::new(Some(3));
let bp = BosonProduct::new([0], [0]).unwrap();
// Adding the product to the system
system
.add_operator_product(
(bp.clone(), bp.clone()), CalculatorComplex::new(1.0, 0.0)
).unwrap();
assert_eq!(system.get(&(bp.clone(), bp)), &CalculatorComplex::new(1.0, 0.0));
println!("{}", system);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import bosons
# Setting up the system and the product we want to add to it
system = bosons.BosonLindbladNoiseSystem(3)
bp = bosons.BosonProduct([0], [0])
# Adding the product to the system
system.add_operator_product((bp, bp), CalculatorComplex.from_pair(1.0, 1.5))
print(system)
# In python we can also use the string representation
system = bosons.BosonLindbladNoiseSystem(3)
system.add_operator_product((str(bp), str(bp)), 1.0+1.5*1j)
print(system)
Open systems
Physically open systems are quantum systems coupled to an environment that can often be described using Lindblad type of noise.
The Lindblad master equation is given by
\[
\dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right)
\]
In struqture they are composed of a Hamiltonian (BosonHamiltonianSystem
) and noise (BosonLindbladNoiseSystem
). They have different ways to set terms in Rust and Python:
Examples
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::bosons::{BosonProduct, HermitianBosonProduct, BosonLindbladOpenSystem};
let mut open_system = BosonLindbladOpenSystem::new(Some(3));
let hbp = HermitianBosonProduct::new([0, 1], [0, 2]).unwrap();
let bp = BosonProduct::new([0], [0]).unwrap();
// Adding the c_0^dag c_1^dag c_0 c_2 term to the system part of the open system
let system = open_system.system_mut();
system.add_operator_product(hbp, CalculatorComplex::new(2.0, 0.0)).unwrap();
// Adding the c_0^dag c_0 part to the noise part of the open system
let noise = open_system.noise_mut();
noise
.add_operator_product((bp.clone(), bp), CalculatorComplex::new(1.0, 0.0))
.unwrap();
println!("{}", open_system);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import bosons
open_system = bosons.BosonLindbladOpenSystem(3)
hbp = bosons.HermitianBosonProduct([0, 1], [0, 2])
bp = bosons.BosonProduct([0], [0])
# Adding the c_0^dag c_1^dag c_0 c_2 term to the system part of the open system
open_system.system_add_operator_product(hbp, CalculatorComplex.from_pair(2.0, 0.0))
# Adding the c_0^dag c_0 part to the noise part of the open system
open_system.noise_add_operator_product(
(bp, bp), CalculatorComplex.from_pair(0.0, 1.0))
print(open_system)
Mixed Systems
Building blocks
All the mixed operators are expressed based on products of mixed indices which contain spin terms, bosonic terms and fermionic terms. The spin terms respect Pauli operator cyclicity, the bosonic terms respect bosonic commutation relations, and the fermionic terms respect fermionic anti-commutation relations.
These products respect the following relations: \[ -i \sigma^x \sigma^y \sigma^z = I \] \[ \lbrack c_{b, k}^{\dagger}, c_{b, j}^{\dagger} \rbrack = 0, \\ \lbrack c_{b, k}, c_{b, j} \rbrack = 0, \\ \lbrack c_{b, k}, c_{b, j}^{\dagger} \rbrack = \delta_{k, j}. \] \[ \lbrace c_{f, k}^{\dagger}, c_{f, j}^{\dagger} \rbrace = 0, \\ \lbrace c_{f, k}, c_{f, j} \rbrace = 0, \\ \lbrace c_{f, k}, c_{f, j}^{\dagger} \rbrace = \delta_{k, j}. \]
with \(c_b^{\dagger}\) the bosonic creation operator, \(c_b\) the bosonic annihilation operator, \(\lbrack ., . \rbrack\) the bosonic commutation relations, \(c_f^{\dagger}\) the fermionic creation operator, \(c_f\) the fermionic annihilation operator, and \(\lbrace ., . \rbrace\) the fermionic anti-commutation relations.
MixedProducts
MixedProducts are combinations of PauliProducts
, BosonProducts
and FermionProducts
.
HermitianMixedProducts
HermitianMixedProducts are the hermitian equivalent of MixedProducts. This means that even though they are constructed the same (see the Examples
section), they internally store both that term and its hermitian conjugate.
MixedDecoherenceProducts
MixedDecoherenceProducts are combinations of DecoherenceProducts
, BosonProducts
and FermionProducts
.
Examples
In both Python and Rust, the operator product is constructed by passing an array/a list of spin terms, an array/a list of bosonic terms and an array/a list of fermionic terms.
from struqture_py import mixed_systems, bosons, spins, fermions
# Building the spin term sigma^x_0 sigma^z_1
pp = spins.PauliProduct().x(0).z(1)
# Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0 * c_b_1
bp = bosons.BosonProduct([1, 2], [1])
# Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
fp = fermions.FermionProduct([0, 1], [0, 1])
# Building the term sigma^x_0 sigma^z_1 c_b^{\dagger}_1 * c_b^{\dagger}_2
# * c_b_0 * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
hmp = mixed_systems.MixedProduct([pp], [bp], [fp])
# Building the term sigma^x_0 sigma^z_1 c_b^{\dagger}_1 * c_b^{\dagger}_2 *
# c_b_0 * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1 + h.c.
hmp = mixed_systems.HermitianMixedProduct([pp], [bp], [fp])
# Building the spin term sigma^x_0 sigma^z_1
dp = spins.DecoherenceProduct().x(0).z(1)
# Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0 * c_b_1
bp = bosons.BosonProduct([1, 2], [0, 1])
# Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
fp = fermions.FermionProduct([0, 1], [0, 1])
# This will work
mdp = mixed_systems.MixedDecoherenceProduct([dp], [bp], [fp])
In Rust the equivalent string representation cannot be used in function and method arguments.
use struqture::prelude::*;
use struqture::spins::{DecoherenceProduct, PauliProduct};
use struqture::bosons::BosonProduct;
use struqture::fermions::FermionProduct;
use struqture::mixed_systems::{
MixedProduct, HermitianMixedProduct, MixedDecoherenceProduct
};
// Building the spin term sigma^x_0 sigma^z_1
let pp = PauliProduct::new().x(0).z(1);
// Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0 * c_b_1
let bp = BosonProduct::new([1, 2], [0, 1]).unwrap();
// Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let fp = FermionProduct::new([0, 1], [0, 1]).unwrap();
// Building the term sigma^x_0 sigma^z_1 c_b^{\dagger}_1 * c_b^{\dagger}_2
// * c_b_0 * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let mp = MixedProduct::new([pp.clone()], [bp.clone()], [fp.clone()]).unwrap();
// Building the term sigma^x_0 sigma^z_1 c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0
// * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1 + h.c.
let hmp = HermitianMixedProduct::new([pp], [bp], [fp]).unwrap();
// Building the spin term sigma^x_0 sigma^z_1
let dp = DecoherenceProduct::new().x(0).z(1);
// Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0 * c_b_1
let bp = BosonProduct::new([1, 2], [0, 1]).unwrap();
// Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let fp = FermionProduct::new([0, 1], [0, 1]).unwrap();
// Building the term sigma^x_0 sigma^z_1 c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0
// * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let mdp = MixedDecoherenceProduct::new([dp], [bp], [fp]).unwrap();
Operators and Hamiltonians
Complex objects are constructed from operator products are MixedOperators
and MixedHamiltonians
(for more information, see also).
These MixedOperators
and MixedHamiltonians
represent operators or Hamiltonians such as:
\[ \hat{H} = \sum_j \alpha_j \prod_k \sigma_{j, k} \prod_{l, m} c_{b, l, j}^{\dagger} c_{b, m, j} \prod_{r, s} c_{f, r, j}^{\dagger} c_{f, s, j} \]
with commutation relations and cyclicity respected.
From a programming perspective the operators and Hamiltonians are HashMaps or Dictionaries with MixedProducts
or HermitianMixedProducts
(respectively) as keys and the coefficients \(\alpha_j\) as values.
In struqture we distinguish between mixed operators and Hamiltonians to avoid introducing unphysical behaviour by accident. While both are sums over normal ordered mixed products (stored as HashMaps of products with a complex prefactor), Hamiltonians are guaranteed to be hermitian to avoid introducing unphysical behaviour by accident. In a mixed Hamiltonian , this means that the sums of products are sums of hermitian mixed products (we have not only the \(c^{\dagger}c\) terms but also their hermitian conjugate) and the on-diagonal terms are required to have real prefactors. We also require the smallest index of the creators to be smaller than the smallest index of the annihilators.
For MixedOperators
and MixedHamiltonians
, we need to specify the number of spin subsystems, bosonic subsystems and fermionic subsystems exist in the operator/Hamiltonian . See the example for more information.
Examples
Here is an example of how to build a MixedOperator
and a MixedHamiltonian
, in Rust:
use qoqo_calculator::CalculatorComplex;
use struqture::prelude::*;
use struqture::bosons::BosonProduct;
use struqture::fermions::FermionProduct;
use struqture::spins::PauliProduct;
use struqture::mixed_systems::{
MixedOperator, MixedProduct, HermitianMixedProduct, MixedHamiltonian
};
// Building the spin term sigma^x_0 sigma^z_1
let pp_0 = PauliProduct::new().x(0).z(1);
// Building the spin term sigma^y_0
let pp_1 = PauliProduct::new().y(0);
// Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0 * c_b_1
let bp = BosonProduct::new([1, 2], [1]).unwrap();
// Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let fp = FermionProduct::new([0, 1], [0, 1]).unwrap();
// Building the term sigma^x_0 sigma^z_1 c_b^{\dagger}_1 * c_b^{\dagger}_2
// * c_b_0 * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let mp_0 = MixedProduct::new([pp_0.clone()], [bp.clone()], [fp.clone()]).unwrap();
// Building the term sigma^y_0
let mp_1 = MixedProduct::new(
[pp_1.clone()],
[BosonProduct::new([], []).unwrap()],
[FermionProduct::new([], []).unwrap()]
).unwrap();
// Building the operator
let mut operator = MixedOperator::new(1, 1, 1);
operator.add_operator_product(mp_0.clone(), CalculatorComplex::new(1.0, 1.5)).unwrap();
operator.add_operator_product(mp_1.clone(), CalculatorComplex::from(2.0)).unwrap();
println!("{}", operator);
// Or when overwriting the previous value
let mut operator = MixedOperator::new(1, 1, 1);
operator.set(mp_0, CalculatorComplex::new(1.0, 1.5)).unwrap();
operator.set(mp_1, CalculatorComplex::from(2.0)).unwrap();
println!("{}", operator);
// Building the term sigma^x_0 sigma^z_1 c_b^{\dagger}_1 * c_b^{\dagger}_2
// * c_b_0 * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let hmp_0 = HermitianMixedProduct::new([pp_0], [bp], [fp]).unwrap();
// Building the term sigma^y_0
let hmp_1 = HermitianMixedProduct::new(
[pp_1],
[BosonProduct::new([], []).unwrap()],
[FermionProduct::new([], []).unwrap()]
).unwrap();
// Building the operator
let mut hamiltonian = MixedHamiltonian::new(1, 1, 1);
hamiltonian.add_operator_product(hmp_0, CalculatorComplex::new(1.0, 1.5)).unwrap();
hamiltonian.add_operator_product(hmp_1, CalculatorComplex::from(2.0)).unwrap();
println!("{}", hamiltonian);
In python, we need to use a MixedSystem
and MixedHamiltonianSystem
instead of a MixedOperator
and MixedHamiltonian
. See next section for more details.
Systems and HamiltonianSystems
Following the intention to avoid unphysical behaviour, MixedSystems and MixedHamiltonianSystems are wrappers around MixedOperators and MixedHamiltonians that allow to explicitly set the number of spins, the number of bosonic modes and the number of fermionic modes in each subsystem of the systems. When setting or adding a MixedProduct/HermitianMixedProduct to the systems, it is guaranteed that the mixed indices involved cannot exceed the number of spins, bosonic modes, and fermionic modes (for each subsystem) in the system. Additionally, the correct number of subsystems needs to have been specified for each type. Note that the user can decide to explicitly set the number of of each particle type to be variable.
Examples
use qoqo_calculator::CalculatorComplex;
use struqture::bosons::BosonProduct;
use struqture::fermions::FermionProduct;
use struqture::mixed_systems::{HermitianMixedProduct, MixedHamiltonianSystem};
use struqture::prelude::*;
use struqture::spins::PauliProduct;
let mut system = MixedHamiltonianSystem::new([Some(2), Some(1)], [Some(3)], [Some(3)]);
// Building the spin term sigma^x_0 sigma^z_1
let pp_0 = PauliProduct::new().x(0).z(1);
// Building the spin term sigma^y_0
let pp_1 = PauliProduct::new().y(0);
// Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0 * c_b_1
let bp = BosonProduct::new([1, 2], [1]).unwrap();
// Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let fp = FermionProduct::new([0, 1], [0, 1]).unwrap();
// This will work
let hmp = HermitianMixedProduct::new(
[pp_0.clone(), pp_1.clone()],
[bp],
[fp.clone()]
).unwrap();
system
.add_operator_product(hmp, CalculatorComplex::new(1.0, 1.5))
.unwrap();
println!("{}", system);
// This will not work, as the number of subsystems of the system and
// product do not match.
let hmp_error = HermitianMixedProduct::new(
[pp_0, pp_1.clone()],
[],
[fp.clone()]
).unwrap();
let error = system.add_operator_product(hmp_error, CalculatorComplex::new(1.0, 1.5));
println!("{:?}", error);
// This will not work, as the number of spins in the first
// subsystem of the system and product do not match.
let pp_error = PauliProduct::new().x(10);
let hmp_error = HermitianMixedProduct::new(
[pp_error.clone(), pp_1.clone()],
[BosonProduct::new([], []).unwrap()],
[fp.clone()],
)
.unwrap();
let error = system.add_operator_product(hmp_error, CalculatorComplex::new(1.0, 1.5));
println!("{:?}", error);
// This will work because we leave the number of spins dynamic
let hmp = HermitianMixedProduct::new(
[pp_error, pp_1],
[BosonProduct::new([], []).unwrap()],
[fp]
).unwrap();
let mut system = MixedHamiltonianSystem::new([None, None], [None], [None]);
system
.add_operator_product(hmp, CalculatorComplex::new(1.0, 0.0))
.unwrap();
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import bosons, fermions, spins, mixed_systems
system = mixed_systems.MixedHamiltonianSystem([2, 1], [3], [3])
# Building the spin term sigma^x_0 sigma^z_1
pp_0 = spins.PauliProduct().x(0).z(1)
# Building the spin term sigma^y_0
pp_1 = spins.PauliProduct().y(0)
# Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_0 * c_b_1
bp = bosons.BosonProduct([1, 2], [1])
# Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
fp = fermions.FermionProduct([0, 1], [0, 1])
# This will work
hmp = mixed_systems.HermitianMixedProduct([pp_0, pp_1], [bp], [fp])
system.add_operator_product(hmp, CalculatorComplex.from_pair(1.0, 1.5))
print(system)
# This will not work, as the number of subsystems of the
# system and product do not match.
hmp_error = mixed_systems.HermitianMixedProduct([pp_0, pp_1], [], [fp])
value = CalculatorComplex.from_pair(1.0, 1.5)
# system.add_operator_product(hmp_error, value) # Uncomment me!
# This will not work, as the number of spins in the first subsystem
# of the system and product do not match.
pp_error = spins.PauliProduct().x(10)
hmp_error = mixed_systems.HermitianMixedProduct(
[pp_error, pp_1], [bosons.BosonProduct([], [])], [fp])
# system.add_operator_product(hmp_error, value) # Uncomment me!
# This will work because we leave the number of spins dynamic.
hmp = mixed_systems.HermitianMixedProduct(
[pp_error, pp_1], [bosons.BosonProduct([], [])], [fp])
system = mixed_systems.MixedHamiltonianSystem([None, None], [None], [None])
system.add_operator_product(hmp_error, CalculatorComplex.from_pair(1.0, 0.0))
Noise operators and systems
We describe decoherence by representing it with the Lindblad equation.
The Lindblad equation is a master equation determining the time evolution of the density matrix.
It is given by
\[
\dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right)
\]
with the rate matrix \(\Gamma_{j,k}\) and the Lindblad operator \(L_{j}\).
To describe the pure noise part of the Lindblad equation one needs the rate matrix in a well defined basis of Lindblad operators.
We use MixedDecoherenceProducts
as the operator basis. To describe mixed noise we use the Lindblad equation with \(\hat{H}=0\).
The rate matrix and with it the Lindblad noise model is saved as a sum over pairs of MixedDecoherenceProducts
, giving the operators acting from the left and right on the density matrix.
In programming terms the object MixedLindbladNoiseOperators
is given by a HashMap or Dictionary with the tuple (MixedDecoherenceProduct
, MixedDecoherenceProduct
) as keys and the entries in the rate matrix as values.
Similarly to MixedOperators, MixedLindbladNoiseOperators have a system equivalent: MixedLindbladNoiseOperators
, with a number of involved spins, bosonic modes and fermionic modes (for each subsystem) defined by the user. For more information on these, see the noise container chapter.
Examples
Here, we add the terms \(L_0 = \left( \sigma_0^x \sigma_1^z \right) \left( c_{b, 1}^{\dagger} c_{b, 1} \right) \left( c_{f, 0}^{\dagger} c_{f, 1}^{\dagger} c_{f, 0} c_{f, 1} \right)\) and \(L_1 = \left( \sigma_0^x \sigma_1^z \right) \left( c_{b, 1}^{\dagger} c_{b, 1} \right) \left( c_{f, 0}^{\dagger} c_{f, 1}^{\dagger} c_{f, 0} c_{f, 1} \right)\) with coefficient 1.0: \( 1.0 \left( L_0 \rho L_1^{\dagger} - \frac{1}{2} \{ L_1^{\dagger} L_0, \rho \} \right) \)
use qoqo_calculator::CalculatorComplex;
use struqture::bosons::BosonProduct;
use struqture::fermions::FermionProduct;
use struqture::mixed_systems::{MixedDecoherenceProduct, MixedLindbladNoiseSystem};
use struqture::prelude::*;
use struqture::spins::DecoherenceProduct;
let mut system = MixedLindbladNoiseSystem::new([Some(3)], [Some(3)], [Some(3)]);
// Building the spin term sigma^x_0 sigma^z_1
let pp = DecoherenceProduct::new().x(0).z(1);
// Building the bosonic term c_b^{\dagger}_1 * c_b_1
let bp = BosonProduct::new([1], [1]).unwrap();
// Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let fp = FermionProduct::new([0, 1], [0, 1]).unwrap();
// Building the term sigma^x_0 sigma^z_1 * c_b^{\dagger}_1
// * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let mdp = MixedDecoherenceProduct::new(
[pp],
[bp],
[fp]
).unwrap();
// Adding in the mixed decoherence product
system
.add_operator_product(
(mdp.clone(), mdp.clone()),
CalculatorComplex::new(1.0, 1.5),
)
.unwrap();
// Checking the system
assert_eq!(
system.get(&(mdp.clone(), mdp)),
&CalculatorComplex::new(1.0, 1.5)
);
println!("{}", system);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import bosons, fermions, spins, mixed_systems
system = mixed_systems.MixedLindbladNoiseSystem([3], [3], [3])
# Building the spin term sigma^x_0 sigma^z_1
pp_0 = spins.DecoherenceProduct().x(0).z(1)
# Building the bosonic term c_b^{\dagger}_0 * c_b_0
bp = bosons.BosonProduct([0], [0])
# Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
fp = fermions.FermionProduct([0, 1], [0, 1])
# Building the term sigma^x_0 sigma^z_1
# * c_b^{\dagger}_0 * c_b_0 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
mdp = mixed_systems.MixedDecoherenceProduct([pp], [bp], [fp])
# Adding in the mixed decoherence product
system.add_operator_product(
(mdp, mdp), CalculatorComplex.from_pair(1.0, 1.5))
print(system)
# In python we can also use the string representation
system = mixed_systems.MixedLindbladNoiseSystem([3], [3], [3])
system.add_operator_product((str(mdp), str(mdp)), 1.0+1.5*1j)
print(system)
Open systems
Physically open systems are quantum systems coupled to an environment that can often be described using Lindblad type of noise. The Lindblad master equation is given by \[ \dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right) \]
In struqture they are composed of a Hamiltonian (MixedHamiltonianSystem) and noise (MixedLindbladNoiseSystem). They have different ways to set terms in Rust and Python:
Examples
use qoqo_calculator::CalculatorComplex;
use struqture::bosons::BosonProduct;
use struqture::fermions::FermionProduct;
use struqture::mixed_systems::{
HermitianMixedProduct, MixedDecoherenceProduct, MixedLindbladOpenSystem,
};
use struqture::prelude::*;
use struqture::spins::{DecoherenceProduct, PauliProduct};
let mut open_system = MixedLindbladOpenSystem::new([Some(3)], [Some(3)], [Some(3)]);
// Building the spin term sigma^x_0 sigma^z_1
let pp = PauliProduct::new().x(0).z(1);
// Building the bosonic term c_b^{\dagger}_0 * c_b_0
let bp = BosonProduct::new([0], [0]).unwrap();
// Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let fp = FermionProduct::new([0, 1], [0, 1]).unwrap();
// Building the term (sigma^x_0 sigma^z_1 * c_b^{\dagger}_0
// * c_b_0 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1)
// + h.c.
let hmp = HermitianMixedProduct::new([pp], [bp], [fp]).unwrap();
// Building the spin term sigma^x_0 sigma^z_1
let pp = DecoherenceProduct::new().x(0).z(1);
// Building the bosonic term c_b^{\dagger}_1 * c_b_1
let bp = BosonProduct::new([1], [1]).unwrap();
// Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let fp = FermionProduct::new([0, 1], [0, 1]).unwrap();
// Building the term sigma^x_0 sigma^z_1 * c_b^{\dagger}_1
// * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
let mdp = MixedDecoherenceProduct::new([pp], [bp], [fp]).unwrap();
// Adding in the system term
let system = open_system.system_mut();
system
.add_operator_product(hmp, CalculatorComplex::new(2.0, 0.0))
.unwrap();
// Adding in the noise term
let noise = open_system.noise_mut();
noise
.add_operator_product((mdp.clone(), mdp), CalculatorComplex::new(1.0, 0.0))
.unwrap();
println!("{}", open_system);
The equivalent code in python:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py import bosons, fermions, spins, mixed_systems
open_system = mixed_systems.MixedLindbladOpenSystem([3], [3], [3])
# Building the spin term sigma^x_0 sigma^z_1
pp = spins.PauliProduct().x(0).z(1)
# Building the bosonic term c_b^{\dagger}_1 * c_b^{\dagger}_2 * c_b_1
bp = bosons.BosonProduct([1, 2], [1])
# Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
fp = fermions.FermionProduct([0, 1], [0, 1])
# Building the term sigma^x_0 sigma^z_1 * c_b^{\dagger}_1
# * c_b^{\dagger}_2 * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
# + h.c.
hmp = mixed_systems.HermitianMixedProduct([pp], [bp], [fp])
# Building the spin term sigma^x_0 sigma^z_1
dp = spins.DecoherenceProduct().x(0).z(1)
# Building the bosonic term c_b^{\dagger}_1 * c_b_1
bp = bosons.BosonProduct([1], [1])
# Building the fermionic term c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
fp = fermions.FermionProduct([0, 1], [0, 1])
# Building the term sigma^x_0 sigma^z_1 * c_b^{\dagger}_1
# * c_b_1 * c_f^{\dagger}_0 * c_f^{\dagger}_1 * c_f_0 * c_f_1
mdp = mixed_systems.MixedDecoherenceProduct([dp], [bp], [fp])
# Adding in the system term
open_system.system_add_operator_product(hmp, CalculatorComplex.from_pair(2.0, 0.0))
# Adding in the noise term
open_system.noise_add_operator_product(
(mdp, mdp), CalculatorComplex.from_pair(0.0, 1.0))
print(open_system)
Container Types
This part of the user documentation focuses on the shared patterns between all physical types: spins, fermions, bosons and mixed systems. All container types for operators, Hamiltonians and open systems behave like hash maps or dictionaries with products of fundamental quantum operators as keys.
The following container types are available, regardless of physical type:
- indices
- operators
- Hamiltonians
- systems
- Hamiltonian systems
- noise operators
- noise systems
- open systems
Products and Indices
The fundamental design of struqture uses products of quantum operators acting on single spins or modes to build up all represented objects. For spins those are SinglePauliOperator
and SingleDecoherenceOperator
and for Fermions and Bosons those are simply fermionic creation and annihilation operators.
Since these operators on single modes or spins form a complete basis of the operator space, each physical object that is represented in struqture can be built up from sum over products of these operators, be it an operator, a Hamiltonian or a noise description.
These sum objects can then be represented in a sparse fashion by saving the sum as a HashMap or Dictionary where the values are the prefactors of the operator products in the sum. The keys of the HashMap are the operator products or for noise objects tuples of operator products.
One of the goals of struqture is to avoid introducing unphysical behaviour by encoding guarantees into the types of operators. For operator products that are not always Hermitian, struqture provides a Hermitian variant of the operator product. This variant picks by design one of the two hermitian conjugated versions of the operator product. It can be used to uniquely represent the coefficient in sum objects that are themselves Hermitian (Hamiltonians) where the coefficients of Hermitian conjugated operator products in the sum also need to be Hermitian conjugated.
The operator products in struqture are
PauliProduct
DecoherenceProduct
FermionProduct
HermitianFermionProduct
BosonProduct
HermitianBosonProdcut
MixedProduct
HermitianMixedProduct
MixedDecoherenceProduct
Operators and Systems
Operators and systems act on a state space using HashMaps (Dictionaries) of operator products and values. The distinction between operators and systems is that systems are given a fixed system size by the user when creating the object.
For spins, the operators and systems represent
\[
\hat{O} = \sum_{j} \alpha_j \prod_{k=0}^N \sigma_{j, k} \\
\sigma_{j, k} \in \{ X_k, Y_k, Z_k, I_k \}
\]
where the \(\sigma_{j, k}\) are SinglePauliOperators
.
For bosons, the operators and systems represent \[ \hat{O} = \sum_{j=0}^N \alpha_j \prod_{k, l} c_{k, j}^{\dagger} c_{l, j} \] with \(c^{\dagger}\) the bosonic creation operator, \(c\) the bosonic annihilation operator \[ \lbrack c_k^{\dagger}, c_j^{\dagger} \rbrack = 0, \\ \lbrack c_k, c_j \rbrack = 0, \\ \lbrack c_k^{\dagger}, c_j \rbrack = \delta_{k, j}. \]
For fermions, the operators and systems represent \[ \hat{O} = \sum_{j=0}^N \alpha_j \prod_{k, l} c_{k, j}^{\dagger} c_{l,j} \] with \(c^{\dagger}\) the fermionionic creation operator, \(c\) the fermionionic annihilation operator \[ \lbrace c_k^{\dagger}, c_j^{\dagger} \rbrace = 0, \\ \lbrace c_k, c_j \rbrace = 0, \\ \lbrace c_k^{\dagger}, c_j \rbrace = \delta_{k, j}. \]
The operators and systems in struqture are
SpinOperator
SpinSystem
DecoherenceOperator
FermionOperator
FermionSystem
BosonOperator
BosonSystem
MixedOperator
MixedSystem
Hamiltonians and HamiltonianSystems
Hamiltonians are hermitian equivalents to Operators, and HamiltonionSystems are the hermitian equivalents to Systems. The operator products for Hamiltonian and Hamiltonian systems are hermitian, meaning that the term is stored, as well as its hermitian conjugate. Also, in order for the Hamiltonian to be hermitian, any operator product on the diagonal of the matrix of interactions must be real.
The Hamiltonians and Hamiltonian systems in struqture are
SpinHamiltonian
SpinHamiltonianSystem
FermionHamiltonian
FermionHamiltonianSystem
BosonHamiltonian
BosonHamiltonianSystem
MixedHamiltonian
MixedHamiltonianSystem
Operators and Systems
Operators and systems act on a state space using HashMaps (Dictionaries) of operator products and values. The distinction between operators and systems is that systems are given a fixed system size by the user when creating the object.
For spins, the operators and systems represent
\[
\hat{O} = \sum_{j} \alpha_j \prod_{k=0}^N \sigma_{j, k} \\
\sigma_{j, k} \in \{ X_k, Y_k, Z_k, I_k \}
\]
where the \(\sigma_{j, k}\) are SinglePauliOperators
.
For bosons, the operators and systems represent \[ \hat{O} = \sum_{j=0}^N \alpha_j \prod_{k, l} c_{k, j}^{\dagger} c_{l, j} \] with \(c^{\dagger}\) the bosonic creation operator, \(c\) the bosonic annihilation operator \[ \lbrack c_k^{\dagger}, c_j^{\dagger} \rbrack = 0, \\ \lbrack c_k, c_j \rbrack = 0, \\ \lbrack c_k^{\dagger}, c_j \rbrack = \delta_{k, j}. \]
For fermions, the operators and systems represent \[ \hat{O} = \sum_{j=0}^N \alpha_j \prod_{k, l} c_{k, j}^{\dagger} c_{l,j} \] with \(c^{\dagger}\) the fermionionic creation operator, \(c\) the fermionionic annihilation operator \[ \lbrace c_k^{\dagger}, c_j^{\dagger} \rbrace = 0, \\ \lbrace c_k, c_j \rbrace = 0, \\ \lbrace c_k^{\dagger}, c_j \rbrace = \delta_{k, j}. \]
The operators and systems in struqture are
SpinOperator
SpinSystem
DecoherenceOperator
FermionOperator
FermionSystem
BosonOperator
BosonSystem
MixedOperator
MixedSystem
Hamiltonians and HamiltonianSystems
Hamiltonians are hermitian equivalents to Operators, and HamiltonionSystems are the hermitian equivalents to Systems. The operator products for Hamiltonian and Hamiltonian systems are hermitian, meaning that the term is stored, as well as its hermitian conjugate. Also, in order for the Hamiltonian to be hermitian, any operator product on the diagonal of the matrix of interactions must be real.
The Hamiltonians and Hamiltonian systems in struqture are
SpinHamiltonian
SpinHamiltonianSystem
FermionHamiltonian
FermionHamiltonianSystem
BosonHamiltonian
BosonHamiltonianSystem
MixedHamiltonian
MixedHamiltonianSystem
Noise Operators
We describe decoherence by representing it with the Lindblad equation. The Lindblad equation is a master equation determining the time evolution of the density matrix. For pure noise terms it is given by \[ \dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right) \] with the rate matrix \(\Gamma_{j,k}\) and the Lindblad operator \(L_{j}\).
Each Lindblad operator is an operator product (in the spins case, a decoherence operator product - for more information see spins container chapter). LindbladNoiseOperators are built as HashMaps (Dictionaries) of Lindblad operators and values, in order to build the non-coherent part of the Lindblad master equation: \[ \sum_{j,k} \Gamma_{j,k} \left( L_{j} \rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho\} \right) \].
The noise operators in struqture are
SpinLindbladNoiseOperator
BosonLindbladNoiseOperator
FermionLindbladNoiseOperator
MixedLindbladNoiseOperator
Noise Systems
Noise systems are noise operators with a fixed size, which is given by the user in the new
function.
The noise systems in struqture are
SpinLindbladNoiseSystem
BosonLindbladNoiseSystem
FermionLindbladNoiseSystem
MixedLindbladNoiseSystem
Open Systems
Open systems represent a full system and environment. Mathematically, this means that a LindbladOpenSystem represents the entire Lindblad equation. The Lindblad equation is a master equation determining the time evolution of the density matrix: \[ \dot{\rho} = \mathcal{L}(\rho) =-i [\hat{H}, \rho] + \sum_{j,k} \Gamma_{j,k} \left( L_{j}\rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho \} \right) \] with the Hamiltonian of the system \(\hat{H}\), the rate matrix \(\Gamma_{j,k}\), and the Lindblad operator \(L_{j}\).
Each LindbladOpenSystem is therefore composed of a HamiltonianSystem: \[ -i [\hat{H}, \rho] \]
and a LindbladNoiseSystem: \[ \sum_{j,k} \Gamma_{j,k} \left( L_{j} \rho L_{k}^{\dagger} - \frac{1}{2} \{ L_k^{\dagger} L_j, \rho\} \right) \]
The open systems in struqture are
SpinLindbladOpenSystem
BosonLindbladOpenSystem
FermionLindbladOpenSystem
MixedLindbladOpenSystem
Applied example
In this example, we will create the spin-boson Hamiltonian we have used for open-system research in our paper, for 1 spin and 3 bosonic modes.
The Hamiltonian is as follows: \[ \hat{H} = \hat{H}_S + \hat{H}_B + \hat{H}_C \]
with the spin system Hamiltonian \(\hat{H}_S\) :
\[ \hat{H} = \frac {\hbar \Delta} {2} \sigma^z_0, \]
the bosonic bath Hamiltonian \(\hat{H}_B\) :
\[ \hat{H} = \sum_{k=0}^2 \hbar \omega_k c_k^{\dagger} c_k, \]
and the coupling between system and bath \(\hat{H}_C\) :
\[ \hat{H} = \sigma_0^x \sum_{k=0}^2 \frac {v_k} {2} \left( c_k + c_k^{\dagger} \right) \]
For simplicity, we will set \(\hbar\) to 1.0 for this example.
Rust implementation:
use qoqo_calculator::CalculatorComplex;
use struqture::bosons::BosonProduct;
use struqture::mixed_systems::{
HermitianMixedProduct, MixedHamiltonianSystem,
};
use struqture::prelude::*;
use struqture::spins::PauliProduct;
let mut system = MixedHamiltonianSystem::new([Some(1)], [Some(3)], []);
// Setting up constants:
let delta = 1.0;
let omega_k = [2.0, 3.0, 4.0];
let v_k = [5.0, 6.0, 7.0];
// First, H_S:
let pp = PauliProduct::new().z(1);
let hmp = HermitianMixedProduct::new(
[pp], [BosonProduct::new([], []).unwrap()], []
).unwrap();
system
.add_operator_product(hmp, CalculatorComplex::new(delta / 2.0, 0.0))
.unwrap();
// Second, H_B:
for k in 0..3 {
let bp = BosonProduct::new([k], [k]).unwrap();
let hmp = HermitianMixedProduct::new(
[PauliProduct::new()], [bp], []
).unwrap();
system
.add_operator_product(
hmp, CalculatorComplex::new(v_k[k] / 2.0, 0.0)
).unwrap();
}
// Third, H_C: the hermitian conjugate is implicitly stored,
// we don't need to add it manually
let pp = PauliProduct::new().x(0);
for k in 0..3 {
let bp = BosonProduct::new([], [k]).unwrap();
let hmp = HermitianMixedProduct::new([pp.clone()], [bp], []).unwrap();
system
.add_operator_product(
hmp, CalculatorComplex::new(omega_k[k], 0.0)
).unwrap();
}
// Our resulting H:
println!("{}", system);
Python implementation:
from qoqo_calculator_pyo3 import CalculatorComplex
from struqture_py.bosons import BosonProduct
from struqture_py.mixed_systems import (
HermitianMixedProduct, HermitianMixedProduct, MixedHamiltonianSystem,
)
from struqture_py.spins import (PauliProduct, PauliProduct)
system = MixedHamiltonianSystem([1], [3], [])
# Setting up constants:
delta = 1.0
omega_k = [2.0, 3.0, 4.0]
v_k = [5.0, 6.0, 7.0]
# First, H_S:
pp = PauliProduct().z(1)
hmp = HermitianMixedProduct([pp], [BosonProduct([], [])], [])
system.add_operator_product(
hmp, CalculatorComplex.from_pair(delta / 2.0, 0.0)
)
# Second, H_B:
for k in range(3):
bp = BosonProduct([k], [k])
hmp = HermitianMixedProduct([PauliProduct()], [bp], [])
system.add_operator_product(
hmp, CalculatorComplex.from_pair(v_k[k] / 2.0, 0.0)
)
# Third, H_C: the hermitian conjugate is implicitly stored,
# we don't need to add it manually
pp = PauliProduct().x(0)
for k in range(3):
bp = BosonProduct([], [k])
hmp = HermitianMixedProduct([pp], [bp], [])
system.add_operator_product(
hmp, CalculatorComplex.from_pair(omega_k[k], 0.0)
)
# Our resulting H:
print(system)