Utilities

Collection of some utilities for simulations.

ITensor Utilities

QuantumDynamics.Utilities.MPO_to_MPSFunction
MPO_to_MPS(ρ::MPO, fbcombiner)

Convert a given MPO to an MPS by combining the two site indices on every tensor into a single one using the vector of combiner tensors.

source
ITensorMPS.expectFunction
ITensorMPS.expect(ρ::MPO, ops; kwargs...)

Extends the ITensorMPS expect function to handle density matrices in the form of MPOs.

source
QuantumDynamics.Utilities.operator_times_identityFunction
operator_times_identity(osum::Opsum, space, iden_space, combiners)

Calculates as an MPO the the operator corresponding to $osum \otimes 1$ where the first space is given by the indices in space and the second space is given by the indices in iden_space. Finally the MPO is put together in a single forward-backward MPO using the combiners .

source

ITensor Utilities for Density Matrix Simulations

QuantumDynamics.Utilities.FBSitesType
FBSites

Stores site information for density matrix simulations.

  • sites_forward: Forward path sites tagged with +
  • sites_backward: Backward path sites tagged with -
  • fb_combiners: ITensor combiners for combining the forward and backward sites
  • sites_fb: Forward-backward sites tagged FBSite
source
QuantumDynamics.Utilities.fb_siteindsFunction
fb_siteinds(sitetype, Nsites::Int; kwargs...)

Generates sites necessary for simulations in the forward-backward Liouville space. This function returns variables of the FBSites type, with forward sites being tagged with + and backward sites being tagged with -. The forward-backward sites are tagged as FBSite, and a combiner is also provided.

Arguments:

  • sitetype: Type of the individual sites. Supports all standard ITensor site types and some custom ones defined here
  • Nsites: Number of sites
  • kwargs...: Other keyword arguments, most important among which are whether quantum numbers should be used
source
QuantumDynamics.Utilities.density_matrix_mpsFunction
density_matrix_mps(ψ::MPS, fbsites::FBSites)

Convert a given MPS wave function defined on fbsites.sites_forward into a density matrix represented by an MPS over the forward-backward sites.

source

HDF5 Utilities

Generic Utilities

QuantumDynamics.Utilities.get_BLAS_implementationFunction
get_BLAS_implementation()

Reports the BLAS implementation under use. The default implementation used by Julia is OpenBlas. MKL is supported through the external package MKL.jl, which needs to be installed and loaded before the loading of QuantumDynamics.jl

source
QuantumDynamics.Utilities.trapezoidFunction
trapezoid(x, y; discrete::Bool=false, exec=ThreadedEx())

Returns the trapezoidal integration of y with respect to x. If discrete is set to true, then returns sum of y. exec encodes the execution paradigm and is one of seq or par.

source
QuantumDynamics.Utilities.fourier_transformFunction
fourier_transform(x::AbstractArray{<:Real}, f::AbstractArray; neg::Bool=true, unitary::Bool=true)

Returns the Fourier transform of f:

$F(k) = \mathcal{N}\,\int_{x_\text{min}}^{x_\text{max}} f(x)\,e^{\pm i k x}\,dx$

Arguments:

  • x: range over which to do the Riemann sum
  • f: function whose Fourier transform is required
  • neg [Default: true]: pick the negative sign in the exponent if neg==true
  • unitary [Default: true]: choose $\mathcal{N}=\frac{1}{\sqrt{2\pi}}$ if unitary==true, or $\mathcal{N}=1$ otherwise
source
QuantumDynamics.Utilities.inverse_fourier_transformFunction
inverse_fourier_transform(k::AbstractArray{<:Real}, F::AbstractArray; x0::Float64=0.0, neg::Bool=true, unitary::Bool=true)

Returns the Fourier transform of f:

$f(x) = \mathcal{N}\,\int_{k_\text{min}}^{k_\text{max}} F(k)\,e^{\mp i k x}\,dk$

Arguments:

  • k: range over which to do the Riemann sum
  • F: function whose inverse Fourier transform is required
  • x0 [Default: 0.0]: the first value of the output coordinate
  • neg [Default: true]: pick the positive sign in the exponent if neg==true. Notice the reversed convention because this is the inverse Fourier transform.
  • unitary [Default: true]: choose $\mathcal{N}=\frac{1}{\sqrt{2\pi}}$ if unitary==true, or $\mathcal{N}=\frac{1}{2\pi}$ otherwise
source
QuantumDynamics.Utilities.calculate_LiouvillianFunction
calculate_Liouvillian(H::OpSum, forward_sites, backward_sites, combiners)

Returns the Liouvillian MPO corresponding to the Hamiltonian provided as an ITensor OpSum, to be built on the forward_sites and backward_sites, and then combined using the combiners provided.

source
calculate_Liouvillian(Hamiltonian::AbstractMatrix{Complex})

Returns the Liouvillian corresponding to the given Hamiltonian.

source
QuantumDynamics.Utilities.unhash_pathFunction
unhash_path(path_num::Int, ntimes::Int, sdim::Int)

Construct a path for a system with sdim dimensions, corresponding to the number path_num, with ntimes time steps.

source

There are a few utilities for creating specific kinds of Hamiltonians.

QuantumDynamics.Utilities.create_nn_hamiltonianFunction
create_nn_hamiltonian(; site_energies::AbstractVector{AbstractFloat}, couplings::AbstractVector{AbstractFloat}, periodic::Bool)

Creates a nearest neighbour Hamiltonian with the given site_energies and couplings. Periodic boundary conditions can also be used by passing true into the periodic argument.

source

Extra arguments

Many of the algorithms require extra, method-specific arguments. These are implemented as subtypes of Utilities.ExtraArgs.

QuantumDynamics.Utilities.DiffEqArgsType

Extra parameters for solving differential equations. Currently has a threshold for magnitude-based filtering. The default values are:

  • reltol = 1e-10
  • abstol = 1e-10
  • solver = Tsit5()
source
QuantumDynamics.Utilities.TensorNetworkArgsType

Extra parameters for tensor network algorithms. Currently has an SVD cutoff, a maximum bond dimension maxdim, and the contraction algorithm. The default values are as follows:

  • cutoff = 1e-8
  • maxdim = 500
  • algorithm = "naive"

Other options for algorithm are "densitymatrix" and "fit".

source