Utilities
Collection of some utilities for simulations.
ITensor Utilities
QuantumDynamics.Utilities.convert_ITensor_to_matrix — Function
convert_ITensor_to_matrix(tens, sinit, sterm)Converts an ITensor with two indices to a matrix. The index sinit is mapped to the column and sterm is mapped to the row.
QuantumDynamics.Utilities.identity_MPO — Function
identity_MPO(sites)Returns the identity MPO based on the given sites.
QuantumDynamics.Utilities.MPO_to_MPS — Function
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.
QuantumDynamics.Utilities.MPS_to_MPO — Function
MPS_to_MPO(ρ::MPS, fbcombiner)Split a given MPS to an MPO by using the vector of combiner tensors.
ITensorMPS.expect — Function
ITensorMPS.expect(ρ::MPO, ops; kwargs...)Extends the ITensorMPS expect function to handle density matrices in the form of MPOs.
QuantumDynamics.Utilities.operator_times_identity — Function
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 .
ITensor Utilities for Density Matrix Simulations
QuantumDynamics.Utilities.forward_backward_combiner — Function
forward_backward_combiner(sites1::Vector, sites2::Vector)Generates combiners for generating forward-backward paths in the Liouville space from separate forward and backward path indices.
QuantumDynamics.Utilities.FBSites — Type
FBSitesStores 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 sitessites_fb: Forward-backward sites tagged FBSite
QuantumDynamics.Utilities.fb_siteinds — Function
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 hereNsites: Number of siteskwargs...: Other keyword arguments, most important among which are whether quantum numbers should be used
QuantumDynamics.Utilities.density_matrix_mps — Function
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.
HDF5 Utilities
QuantumDynamics.Utilities.create_and_select_group — Function
create_and_select_group(base, new_group)Checks if new_group exists in base. Selects it if it does, else creates it.
QuantumDynamics.Utilities.check_or_insert_value — Function
check_or_insert_value(base, variable, value)Inserts value into base[variable]. If it already exists, checks if the value is correct, and throws if not.
QuantumDynamics.Utilities.merge_into — Function
merge_into(source::String, destination::String)Merge data from the HDF5 file at source to the one at destination.
Generic Utilities
QuantumDynamics.Utilities.get_BLAS_implementation — Function
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
QuantumDynamics.Utilities.trapezoid — Function
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.
QuantumDynamics.Utilities.fourier_transform — Function
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 sumf: function whose Fourier transform is requiredneg[Default:true]: pick the negative sign in the exponent ifneg==trueunitary[Default:true]: choose $\mathcal{N}=\frac{1}{\sqrt{2\pi}}$ ifunitary==true, or $\mathcal{N}=1$ otherwise
QuantumDynamics.Utilities.inverse_fourier_transform — Function
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 sumF: function whose inverse Fourier transform is requiredx0[Default:0.0]: the first value of the output coordinateneg[Default:true]: pick the positive sign in the exponent ifneg==true. Notice the reversed convention because this is the inverse Fourier transform.unitary[Default:true]: choose $\mathcal{N}=\frac{1}{\sqrt{2\pi}}$ ifunitary==true, or $\mathcal{N}=\frac{1}{2\pi}$ otherwise
QuantumDynamics.Utilities.commutator — Function
commutator(A, B)Returns the commutator A and B: AB - BA.
QuantumDynamics.Utilities.nh_commutator — Function
nh_commutator(A, B)Returns the commutator A and B: AB - BA'.
QuantumDynamics.Utilities.calculate_Liouvillian — Function
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.
calculate_Liouvillian(Hamiltonian::AbstractMatrix{Complex})Returns the Liouvillian corresponding to the given Hamiltonian.
QuantumDynamics.Utilities.ExternalField — Type
ExternalField provides an abstract interface for encoding an external field, V(t), interacting with the system through the operator, coupling_op.
QuantumDynamics.Utilities.density_matrix_to_vector — Function
density_matrix_to_vector(ρ::AbstractMatrix{<:Complex})Returns the vector representation of the density matrix ρ compatible with the forward-backward propagators.
QuantumDynamics.Utilities.density_matrix_vector_to_matrix — Function
density_matrix_vector_to_matrix(ρvec::AbstractVector{<:Complex})Returns the matrix form of the vector ρvec.
QuantumDynamics.Utilities.hash_path — Function
hash_path(states, sdim)Returns the hashed location of a path for a system with sdim dimensions.
QuantumDynamics.Utilities.unhash_path — Function
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.
QuantumDynamics.Utilities.unhash_path_blips — Function
unhash_path_blips(ntimes::Int, sdim::Int, nblips::Int)Construct all the paths for a system with sdim dimensions with ntimes time steps and nblips blips.
QuantumDynamics.Utilities.apply_propagator — Function
apply_propagator(; propagators, ρ0, ntimes, dt)Apply a series of ntimes propagators to an initial reduced density matrix ρ0 and return the result as a tuple of (time, ρs).
There are a few utilities for creating specific kinds of Hamiltonians.
QuantumDynamics.Utilities.create_tls_hamiltonian — Function
create_tls_hamiltonian(; ϵ::AbstractFloat, Δ::AbstractFloat)Creates a two-level system Hamiltonian:
$H = \frac{ϵ}{2}σ_z - \frac{Δ}{2}σ_x$
QuantumDynamics.Utilities.create_nn_hamiltonian — Function
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.
Extra arguments
Many of the algorithms require extra, method-specific arguments. These are implemented as subtypes of Utilities.ExtraArgs.
QuantumDynamics.Utilities.ExtraArgs — Type
Abstract type for encoding all the method specific numerical parameters.
QuantumDynamics.Utilities.DiffEqArgs — Type
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()
QuantumDynamics.Utilities.TensorNetworkArgs — Type
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".